dune-fem  2.4.1-rc
petscsolver.hh
Go to the documentation of this file.
1 #ifndef DUNE_FEM_PETSCSOLVER_HH
2 #define DUNE_FEM_PETSCSOLVER_HH
3 
4 #include <limits>
5 
9 
10 #if HAVE_PETSC
12 
13 namespace Dune
14 {
15 
16  namespace Fem
17  {
18 
19  //=====================================================================
20  // Implementation for PETSc matrix based Krylov solvers
21  //=====================================================================
22 
27  // PETScSolver
28  // --------------
29 
31  template< class DF, class Op >
32  class PetscInverseOperator
33  : public Operator< DF, DF >
34  {
35  protected:
36  // monitor function for PETSc solvers
37  static PetscErrorCode
38  monitor (KSP ksp, int it, PetscReal rnorm, void *mctx)
39  {
40  if( Parameter :: verbose () )
41  {
42  std::cout << "PETSc::KSP: it = " << it << " res = " << rnorm << std::endl;
43  }
44  return PetscErrorCode(0);
45  }
46 
47  public:
48  typedef DF DiscreteFunctionType;
49  typedef DiscreteFunctionType DestinationType;
50  typedef Op OperatorType;
51 
52  public:
63  PetscInverseOperator ( const OperatorType &op,
64  double reduction,
65  double absLimit,
66  int maxIter,
67  bool verbose,
68  const ParameterReader &parameter = Parameter::container() )
69  : op_( op ),
70  reduction_( reduction ),
71  absLimit_( absLimit ),
72  maxIter_( maxIter ),
73  verbose_( verbose ),
74  iterations_( 0 ),
75  solverName_("none"),
76  precondName_("none")
77  {
78  initialize( parameter );
79  }
80 
88  PetscInverseOperator ( const OperatorType &op,
89  double reduction,
90  double absLimit,
91  int maxIter,
92  const ParameterReader &parameter = Parameter::container() )
93  : op_( op ),
94  reduction_( reduction ),
95  absLimit_ ( absLimit ),
96  maxIter_( maxIter ),
97  verbose_( parameter.getValue< bool >( "fem.solver.verbose", false ) ),
98  iterations_( 0 ),
99  solverName_("none"),
100  precondName_("none")
101  {
102  initialize( parameter );
103  }
104 
105  PetscInverseOperator ( const OperatorType &op,
106  double reduction,
107  double absLimit,
108  const ParameterReader &parameter = Parameter::container() )
109  : op_( op ),
110  reduction_( reduction ),
111  absLimit_ ( absLimit ),
112  maxIter_( std::numeric_limits< int >::max()),
113  verbose_( parameter.getValue< bool >( "fem.solver.verbose", false ) ),
114  iterations_( 0 ),
115  solverName_("none"),
116  precondName_("none")
117  {
118  initialize( parameter );
119  }
120 
122  ~PetscInverseOperator()
123  {
124  // destroy PC context
125  ::Dune::Petsc::PCDestroy( &pc_ );
126  // destroy solver context
127  ::Dune::Petsc::KSPDestroy( &ksp_ );
128  }
129 
130  void initialize ( const ParameterReader &parameter )
131  {
132  // Create linear solver context
133  ::Dune::Petsc::KSPCreate( &ksp_ );
134 
135  enum PetscSolver { petsc_cg = 0,
136  petsc_bicg = 1,
137  petsc_bicgstab = 2,
138  petsc_gmres = 3
139  };
140 
141  // see PETSc docu for more types
142  const std::string solverNames [] = { "cg" , "bicg", "bicgstab", "gmres" };
143  PetscSolver solverType = (PetscSolver) parameter.getEnum("petsc.kspsolver.method", solverNames, 0 );
144 
145  // select linear solver
146  if ( solverType == petsc_cg )
147  ::Dune::Petsc::KSPSetType( ksp_, KSPCG );
148  else if ( solverType == petsc_bicg )
149  ::Dune::Petsc::KSPSetType( ksp_, KSPBICG ); // this is the BiCG version of PETSc.
150  else if ( solverType == petsc_bicgstab )
151  ::Dune::Petsc::KSPSetType( ksp_, KSPBCGS ); // this is the BiCG-stab version of PETSc.
152  else if ( solverType == petsc_gmres )
153  {
154  PetscInt restart = parameter.getValue<int>("petsc.gmresrestart", 10 );
155  ::Dune::Petsc::KSPSetType( ksp_, KSPGMRES );
156  ::Dune::Petsc::KSPGMRESSetRestart( ksp_, restart );
157  }
158  else
159  {
160  DUNE_THROW(InvalidStateException,"PetscInverseOperator: wrong solver choosen" );
161  }
162 
163  // store solver name
164  solverName_ = solverNames[ solverType ];
165 
167  // preconditioning
169 
170  enum PetscPcType { petsc_none = 0, // no preconditioning
171  // parallel preconditioners
172  petsc_asm = 1, // Additive Schwarz
173  petsc_sor = 2, // SOR and SSOR
174  petsc_jacobi = 3, // Jacobi preconditioning
175  // requiring additional packages
176  petsc_hypre = 4, // Hypre preconditioning
177  petsc_ml = 5, // ML preconditioner (from Trilinos)
178  // serial preconditioners
179  petsc_ilu = 6, // ILU preconditioning
180  petsc_lu = 7, // LU factorization
181  petsc_icc = 8, // Incomplete Cholesky factorization
182  // direct solvers
183  petsc_mumps = 9, // use mumps
184  petsc_superlu = 10, // use superlu-dist
185  end = 11
186  };
187  const PetscPcType serialStart = petsc_ilu;
188  const PetscPcType serialEnd = petsc_icc;
189 
190  PCType type = PCNONE ;
191 
192  // see PETSc docu for more types
193  const std::string preconNames [end] = { "none", "asm", "sor", "jacobi",
194  "hypre", "ml",
195  "ilu-n", "lu", "icc",
196  "mumps", "superlu" };
197  PetscPcType pcType = (PetscPcType) parameter.getEnum("petsc.preconditioning.method", preconNames, 0 );
198 
199  if( pcType == petsc_none )
200  type = PCNONE ;
201  else if( pcType == petsc_asm )
202  type = PCASM;
203  else if( pcType == petsc_sor )
204  type = PCSOR;
205  else if ( pcType == petsc_jacobi )
206  type = PCJACOBI;
207  else if ( pcType == petsc_hypre )
208  {
209  type = PCHYPRE;
210  }
211  else if ( pcType == petsc_ml )
212  type = PCML;
213  else if ( pcType == petsc_ilu )
214  type = PCILU;
215  else if( pcType == petsc_lu )
216  type = PCLU;
217  else if( pcType == petsc_icc )
218  type = PCICC;
219  else if( pcType == petsc_mumps )
220  {
221  /* could use
222  KSPSetType(ksp,KSPPREONLY);
223  */
224  type = PCLU;
225  }
226  else if( pcType == petsc_superlu )
227  {
228  /* could use
229  KSPSetType(ksp,KSPPREONLY);
230  */
231  type = PCLU;
232  }
233  else
234  DUNE_THROW(InvalidStateException,"PetscInverseOperator: wrong preconditiong choosen" );
235 
236  // check whether the preconditioner can be used in parallel
237  if( MPIManager :: size() > 1 && pcType >= serialStart && pcType <= serialEnd )
238  {
239  if( MPIManager :: rank() == 0 )
240  {
241  std::cerr << "WARNING: PetscInverseOperator: " << preconNames[ pcType ] << " preconditioning does not work in parallel, falling back to Additive Schwarz!!" << std::endl;
242  }
243  type = PCASM;
244  }
245 
246  // get level of perconditioner iterations
247  PetscInt pcLevel = 1 + parameter.getValue<int>("petsc.preconditioning.iterations", 0 );
248 
249  // create preconditioning context
250  ::Dune::Petsc::PCCreate( &pc_ );
251  ::Dune::Petsc::PCSetType( pc_, type );
252  ::Dune::Petsc::PCFactorSetLevels( pc_, pcLevel );
253 
254  if ( pcType == petsc_hypre )
255  {
256  // set type of HYPRE preconditioner to boomer-amg
257  // there are also other preconditioners in this package
258  ::Dune::Petsc::PCHYPRESetType( pc_, "boomeramg" );
259  }
260 
261  // set preconditioning context
262  ::Dune::Petsc::KSPSetPC( ksp_, pc_ );
263  if( pcType == petsc_mumps )
264  ::Dune::Petsc::PCFactorSetMatSolverPackage(pc_,MATSOLVERMUMPS);
265  if( pcType == petsc_superlu )
266  ::Dune::Petsc::PCFactorSetMatSolverPackage(pc_,MATSOLVERSUPERLU_DIST);
267 
268  // get matrix from linear operator
269  Mat& A = const_cast< Mat & > (op_.petscMatrix());
270 
271  // set operator to PETSc solver context
272  // ::Dune::Petsc::KSPSetOperators( ksp_, A, A, DIFFERENT_NONZERO_PATTERN);
273 #if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 5
274  ::Dune::Petsc::KSPSetOperators( ksp_, A, A, SAME_PRECONDITIONER);
275 #else
276  ::Dune::Petsc::KSPSetOperators( ksp_, A, A );
277 #endif
278  // set prescribed tolerances
279  PetscInt maxits = maxIter_ ;
280  PetscReal reduc = reduction_;
281  ::Dune::Petsc::KSPSetTolerances(ksp_, reduc, 1.e-50, PETSC_DEFAULT, maxits);
282 
283  // set monitor in verbose mode
284  if( verbose_ )
285  {
286  ::Dune::Petsc::KSPView( ksp_ );
287  ::Dune::Petsc::KSPMonitorSet( ksp_, &monitor, PETSC_NULL, PETSC_NULL);
288  }
289  }
290 
291  void prepare (const DiscreteFunctionType& Arg, DiscreteFunctionType& Dest) const
292  {
293  }
294 
295  void finalize () const
296  {}
297 
298  void printTexInfo(std::ostream& out) const
299  {
300  out << "Solver: " << solverName_ << " eps = " << reduction_ ;
301  out << "\\\\ \n";
302  }
303 
308  void apply( const DiscreteFunctionType& arg, DiscreteFunctionType& dest ) const
309  {
310  // call PETSc solvers
311  ::Dune::Petsc::KSPSolve(ksp_, *arg.petscVec() , *dest.petscVec() );
312 
313  // for continuous solution we need a communication here
314  if( dest.space().continuous() )
315  {
316  dest.communicate();
317  }
318 
319  // get number of iterations
320  PetscInt its ;
321  ::Dune::Petsc::KSPGetIterationNumber( ksp_, &its );
322  iterations_ = its;
323  }
324 
325  // return number of iterations
326  int iterations() const
327  {
328  return iterations_;
329  }
330 
332  double averageCommTime() const
333  {
334  return -1;
335  }
336 
341  void operator() ( const DiscreteFunctionType& arg, DiscreteFunctionType& dest ) const
342  {
343  apply(arg,dest);
344  }
345 
346  private:
347  // no copying
348  PetscInverseOperator ();
349  PetscInverseOperator( const PetscInverseOperator& ) ;
350  PetscInverseOperator& operator= ( const PetscInverseOperator& );
351 
352  protected:
353  const OperatorType &op_; // linear operator
354  KSP ksp_; // PETSc Krylov Space solver context
355  PC pc_; // PETSc perconditioning context
356  double reduction_;
357  double absLimit_;
358  int maxIter_;
359  bool verbose_ ;
360  mutable int iterations_;
361  std::string solverName_;
362  std::string precondName_;
363  };
364 
366 
367  } // namespace Fem
368 
369 } // namespace Dune
370 
371 #endif // #if HAVE_PETSC
372 
373 #endif // #ifndef DUNE_FEM_PETSCSOLVER_HH
static int rank()
Definition: mpimanager.hh:116
BasicParameterReader< std::function< const std::string *(const std::string &, const std::string *) > > ParameterReader
Definition: reader.hh:264
static int size()
Definition: mpimanager.hh:121
static double max(const Double &v, const double p)
Definition: double.hh:387
Definition: coordinate.hh:4
static ParameterContainer & container()
Definition: io/parameter.hh:190
STL namespace.
static bool verbose()
obtain the cached value for fem.verbose
Definition: io/parameter.hh:444