LifeV
PreconditionerAztecOO.cpp
Go to the documentation of this file.
1 //@HEADER
2 /*
3 *******************************************************************************
4 
5  Copyright (C) 2004, 2005, 2007 EPFL, Politecnico di Milano, INRIA
6  Copyright (C) 2010 EPFL, Politecnico di Milano, Emory University
7 
8  This file is part of LifeV.
9 
10  LifeV is free software; you can redistribute it and/or modify
11  it under the terms of the GNU Lesser General Public License as published by
12  the Free Software Foundation, either version 3 of the License, or
13  (at your option) any later version.
14 
15  LifeV is distributed in the hope that it will be useful,
16  but WITHOUT ANY WARRANTY; without even the implied warranty of
17  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18  Lesser General Public License for more details.
19 
20  You should have received a copy of the GNU Lesser General Public License
21  along with LifeV. If not, see <http://www.gnu.org/licenses/>.
22 
23 *******************************************************************************
24 */
25 //@HEADER
26 
27 /*!
28  @file
29  @brief AztecOO preconditioner
30 
31  @author Cristiano Malossi <cristiano.malossi@epfl.ch>
32  @contributor Gwenol Grandperrin <gwenol.grandperrin@epfl.ch>
33  @maintainer Gwenol Grandperrin <gwenol.grandperrin@epfl.ch>
34 
35  @date 17-11-2009
36  */
37 
39 #include <lifev/core/LifeV.hpp>
40 
41 namespace LifeV
42 {
43 
44 // ===================================================
45 // Constructors & Destructor
46 // ===================================================
48  super (),
49  M_solver()
50 {
51 
52 #ifdef HAVE_LIFEV_DEBUG
53  debugStream ( 7100 ) << "PreconditionerAztecOO::PreconditionerAztecOO() \n";
54 #endif
55 
56 }
57 
58 
59 // ===================================================
60 // Methods
61 // ===================================================
62 Int
63 PreconditionerAztecOO::buildPreconditioner ( operator_type& matrix )
64 {
65 
66 #ifdef HAVE_LIFEV_DEBUG
67  debugStream ( 7100 ) << "PreconditionerAztecOO::buildPreconditioner( Operator ) \n";
68 #endif
69 
70  if ( this->M_preconditionerCreated )
71  {
73  }
74 
75  M_solver->solver().SetPrecMatrix ( matrix->matrixPtr().get() );
76 
77  M_solver->solver().SetAztecOption ( AZ_pre_calc, AZ_calc );
78  M_solver->solver().SetAztecOption ( AZ_keep_info, 1 );
79 
80  Real estimateConditionNumber;
81  M_solver->solver().ConstructPreconditioner ( estimateConditionNumber );
82 
83  this->M_preconditionerCreated = true;
84 
85  return ( EXIT_SUCCESS );
86 }
87 
88 void
90 {
91 
92 #ifdef HAVE_LIFEV_DEBUG
93  debugStream ( 7100 ) << "PreconditionerAztecOO::precReset() \n";
94 #endif
95 
96  M_solver->solver().SetAztecOption ( AZ_keep_info, 0 );
97  M_solver->solver().SetAztecOption ( AZ_pre_calc, AZ_calc );
98 
99  // Perform one "fake" iteration to delete the preconditioner
100  Int AZoutputOption = M_solver->solver().GetAztecOption ( AZ_output );
101  M_solver->solver().SetAztecOption ( AZ_output, AZ_none );
102  M_solver->solver().Iterate ( 0, 1.e14 );
103  M_solver->solver().SetAztecOption ( AZ_output, AZoutputOption );
104 
105  this->M_preconditionerCreated = false;
106 }
107 
108 void
109 PreconditionerAztecOO::createParametersList ( list_Type& /*list*/, const GetPot& dataFile, const std::string& section, const std::string& subSection )
110 {
111  // Preconditioner
112  M_solver->getParametersList().set ( "precond", dataFile ( ( section + "/" + subSection + "/precond" ).data(), "dom_decomp" ) );
113 
114  // Compute the preconditioner
115  M_solver->getParametersList().set ( "pre_calc", dataFile ( ( section + "/" + subSection + "/pre_calc" ).data(), "calc" ) );
116 
117  // Reordering
118  M_solver->getParametersList().set ( "reorder", dataFile ( ( section + "/" + subSection + "/reorder" ).data(), 1 ) );
119 
120  // Keep the information
121  M_solver->getParametersList().set ( "keep_info", dataFile ( ( section + "/" + subSection + "/keep_info" ).data(), 1 ) );
122 
123  // Polynomial order when using Jacobi and GS preconditioners
124  //M_solver->getParametersList().set( "poly_ord", dataFile( ( section + "/" + subSection + "/poly_ord" ).data(), 3 ) );
125 
126  // Overlap level
127  M_solver->getParametersList().set ( "overlap", dataFile ( ( section + "/" + subSection + "/overlap" ).data(), AZ_none ) );
128 
129  // Overlap type
130  M_solver->getParametersList().set ( "type_overlap", dataFile ( ( section + "/" + subSection + "/type_overlap" ).data(), AZ_standard ) );
131 
132 
133  // SUBDOMAIN SOLVER
134 
135  M_solver->getParametersList().set ( "subdomain_solve", dataFile ( ( section + "/" + subSection + "/subdomain_solve" ).data(), "ILUT" ) );
136 
137  M_solver->getParametersList().set ( "drop", dataFile ( ( section + "/" + subSection + "/drop" ).data(), 1.e-5 ) );
138 
139  //M_solver->getParametersList().set( "graph_fill", dataFile( ( section + "/" + subSection + "/graph_fill" ).data(), 6. ) );
140 
141  M_solver->getParametersList().set ( "ilut_fill", dataFile ( ( section + "/" + subSection + "/ilut_fill" ).data(), 4. ) );
142 
143  //M_solver->getParametersList().set( "init_guess", dataFile( ( section + "/" + subSection + "/init_guess" ).data(), AZ_NOT_ZERO ) );
144 
145  //M_solver->getParametersList().set( "keep_kvecs", dataFile( ( section + "/" + subSection + "/keep_kvecs" ).data(), 0 ) );
146 
147  //M_solver->getParametersList().set( "apply_kvecs", dataFile( ( section + "/" + subSection + "/apply_kvecs" ).data(), AZ_FALSE ) );
148 
149  //M_solver->getParametersList().set( "orth_kvecs", dataFile( ( section + "/" + subSection + "/orth_kvecs" ).data(), AZ_FALSE ) );
150 
151  //M_solver->getParametersList().set( "ignore_scaling", dataFile( ( section + "/" + subSection + "/ignore_scaling" ).data(), AZ_FALSE ) );
152 
153  //M_solver->getParametersList().set( "check_update_size", dataFile( ( section + "/" + subSection + "/check_update_size" ).data(), AZ_FALSE ) );
154 
155  //M_solver->getParametersList().set( "omega", dataFile( ( section + "/" + subSection + "/omega" ).data(), 1. ) );
156 
157  //M_solver->getParametersList().set( "update_reduction", dataFile( ( section + "/" + subSection + "/update_reduction" ).data(), 10e10 ) );
158 
159  // SET PARAMETERS
161 
162  // DISPLAY LIST
163  if ( dataFile ( (section + "/displayList").data(), false ) )
164  {
165  M_solver->getParametersList().print ( std::cout );
166  }
167 }
168 
169 void
170 PreconditionerAztecOO::showMe ( std::ostream& output ) const
171 {
172  output << "showMe must be implemented for the PreconditionerAztecOO class" << std::endl;
173 }
174 
175 Real
177 {
178 
179 #ifdef HAVE_LIFEV_DEBUG
180  debugStream ( 7100 ) << "PreconditionerAztecOO::Condest() \n";
181 #endif
182 
183  return M_solver->solver().Condest();
184 }
185 
186 // ===================================================
187 // Set Methods
188 // ===================================================
189 void
191  const std::string& section )
192 {
193 
194 #ifdef HAVE_LIFEV_DEBUG
195  debugStream ( 7100 ) << "PreconditionerAztecOO::setDataFromGetPot(dataFile, section) \n";
196 #endif
197 
198  createParametersList ( M_list, dataFile, section, "AztecOO" );
199 }
200 
201 
202 // ===================================================
203 // Get Methods
204 // ===================================================
207 {
208 
209 #ifdef HAVE_LIFEV_DEBUG
210  debugStream ( 7100 ) << "PreconditionerAztecOO::getPrec() \n";
211 #endif
212 
213  if ( this->M_preconditionerCreated )
214  {
215  return M_solver->solver().GetPrecMatrix();
216  }
217 
218  return 0;
219 }
220 
223 {
224 
225 #ifdef HAVE_LIFEV_DEBUG
226  debugStream ( 7100 ) << "PreconditionerAztecOO::getPrec() \n";
227 #endif
228 
229  std::shared_ptr<Epetra_RowMatrix> prec;
230  if ( this->M_preconditionerCreated )
231  {
232  prec.reset (M_solver->solver().GetPrecMatrix() );
233  return prec;
234  }
235  return Preconditioner::prec_type();
236 }
237 
238 } // namespace LifeV
Int buildPreconditioner(operator_type &matrix)
Build the preconditioner.
void createParametersList(list_Type &list, const GetPot &dataFile, const std::string &section, const std::string &subSection)
Create the list of parameters of the preconditioner.
void setDataFromGetPot(const GetPot &dataFile, const std::string &section)
Set the data of the preconditioner using a GetPot object.
void resetPreconditioner()
Reset the preconditioner.
Real condest()
Compute the condition number of the preconditioner.
int32_type Int
Generic integer data.
Definition: LifeV.hpp:188
void updateInverseJacobian(const UInt &iQuadPt)
PreconditionerAztecOO - The implementation of Preconditioner for AztecOO preconditioners.
Epetra_Operator prec_raw_type
Preconditioner(const commPtr_Type &comm=commPtr_Type())
Constructor.
double Real
Generic real data.
Definition: LifeV.hpp:175
Teuchos::ParameterList list_Type
Preconditioner - Abstract preconditioner class.
std::shared_ptr< prec_raw_type > prec_type
void setParameters(bool cerrWarningIfUnused=false)
Set the current parameters with the internal parameters list.
super::prec_raw_type * preconditioner()
Return the pointer to the preconditioner.
super::prec_type preconditionerPtr()
Return the shared pointer to the preconditioner.
virtual void showMe(std::ostream &output=std::cout) const
Show informations about the preconditioner.