1 #ifndef DUNE_FEM_PETSCSOLVER_HH 2 #define DUNE_FEM_PETSCSOLVER_HH 31 template<
class DF,
class Op >
32 class PetscInverseOperator
33 :
public Operator< DF, DF >
38 monitor (KSP ksp,
int it, PetscReal rnorm,
void *mctx)
42 std::cout <<
"PETSc::KSP: it = " << it <<
" res = " << rnorm << std::endl;
44 return PetscErrorCode(0);
48 typedef DF DiscreteFunctionType;
49 typedef DiscreteFunctionType DestinationType;
50 typedef Op OperatorType;
63 PetscInverseOperator (
const OperatorType &op,
70 reduction_( reduction ),
71 absLimit_( absLimit ),
78 initialize( parameter );
88 PetscInverseOperator (
const OperatorType &op,
94 reduction_( reduction ),
95 absLimit_ ( absLimit ),
97 verbose_( parameter.getValue< bool >(
"fem.solver.verbose", false ) ),
102 initialize( parameter );
105 PetscInverseOperator (
const OperatorType &op,
110 reduction_( reduction ),
111 absLimit_ ( absLimit ),
112 maxIter_(
std::numeric_limits< int >::
max()),
113 verbose_( parameter.getValue< bool >(
"fem.solver.verbose", false ) ),
118 initialize( parameter );
122 ~PetscInverseOperator()
125 ::Dune::Petsc::PCDestroy( &pc_ );
127 ::Dune::Petsc::KSPDestroy( &ksp_ );
133 ::Dune::Petsc::KSPCreate( &ksp_ );
135 enum PetscSolver { petsc_cg = 0,
142 const std::string solverNames [] = {
"cg" ,
"bicg",
"bicgstab",
"gmres" };
143 PetscSolver solverType = (PetscSolver) parameter.getEnum(
"petsc.kspsolver.method", solverNames, 0 );
146 if ( solverType == petsc_cg )
147 ::Dune::Petsc::KSPSetType( ksp_, KSPCG );
148 else if ( solverType == petsc_bicg )
149 ::Dune::Petsc::KSPSetType( ksp_, KSPBICG );
150 else if ( solverType == petsc_bicgstab )
151 ::Dune::Petsc::KSPSetType( ksp_, KSPBCGS );
152 else if ( solverType == petsc_gmres )
154 PetscInt restart = parameter.getValue<
int>(
"petsc.gmresrestart", 10 );
155 ::Dune::Petsc::KSPSetType( ksp_, KSPGMRES );
156 ::Dune::Petsc::KSPGMRESSetRestart( ksp_, restart );
160 DUNE_THROW(InvalidStateException,
"PetscInverseOperator: wrong solver choosen" );
164 solverName_ = solverNames[ solverType ];
170 enum PetscPcType { petsc_none = 0,
187 const PetscPcType serialStart = petsc_ilu;
188 const PetscPcType serialEnd = petsc_icc;
190 PCType type = PCNONE ;
193 const std::string preconNames [end] = {
"none",
"asm",
"sor",
"jacobi",
195 "ilu-n",
"lu",
"icc",
196 "mumps",
"superlu" };
197 PetscPcType pcType = (PetscPcType) parameter.getEnum(
"petsc.preconditioning.method", preconNames, 0 );
199 if( pcType == petsc_none )
201 else if( pcType == petsc_asm )
203 else if( pcType == petsc_sor )
205 else if ( pcType == petsc_jacobi )
207 else if ( pcType == petsc_hypre )
211 else if ( pcType == petsc_ml )
213 else if ( pcType == petsc_ilu )
215 else if( pcType == petsc_lu )
217 else if( pcType == petsc_icc )
219 else if( pcType == petsc_mumps )
226 else if( pcType == petsc_superlu )
234 DUNE_THROW(InvalidStateException,
"PetscInverseOperator: wrong preconditiong choosen" );
241 std::cerr <<
"WARNING: PetscInverseOperator: " << preconNames[ pcType ] <<
" preconditioning does not work in parallel, falling back to Additive Schwarz!!" << std::endl;
247 PetscInt pcLevel = 1 + parameter.getValue<
int>(
"petsc.preconditioning.iterations", 0 );
250 ::Dune::Petsc::PCCreate( &pc_ );
251 ::Dune::Petsc::PCSetType( pc_, type );
252 ::Dune::Petsc::PCFactorSetLevels( pc_, pcLevel );
254 if ( pcType == petsc_hypre )
258 ::Dune::Petsc::PCHYPRESetType( pc_,
"boomeramg" );
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);
269 Mat& A =
const_cast< Mat &
> (op_.petscMatrix());
273 #if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 5 274 ::Dune::Petsc::KSPSetOperators( ksp_, A, A, SAME_PRECONDITIONER);
276 ::Dune::Petsc::KSPSetOperators( ksp_, A, A );
279 PetscInt maxits = maxIter_ ;
280 PetscReal reduc = reduction_;
281 ::Dune::Petsc::KSPSetTolerances(ksp_, reduc, 1.e-50, PETSC_DEFAULT, maxits);
286 ::Dune::Petsc::KSPView( ksp_ );
287 ::Dune::Petsc::KSPMonitorSet( ksp_, &monitor, PETSC_NULL, PETSC_NULL);
291 void prepare (
const DiscreteFunctionType& Arg, DiscreteFunctionType& Dest)
const 295 void finalize ()
const 298 void printTexInfo(std::ostream& out)
const 300 out <<
"Solver: " << solverName_ <<
" eps = " << reduction_ ;
308 void apply(
const DiscreteFunctionType& arg, DiscreteFunctionType& dest )
const 311 ::Dune::Petsc::KSPSolve(ksp_, *arg.petscVec() , *dest.petscVec() );
314 if( dest.space().continuous() )
321 ::Dune::Petsc::KSPGetIterationNumber( ksp_, &its );
326 int iterations()
const 332 double averageCommTime()
const 341 void operator() (
const DiscreteFunctionType& arg, DiscreteFunctionType& dest )
const 348 PetscInverseOperator ();
349 PetscInverseOperator(
const PetscInverseOperator& ) ;
350 PetscInverseOperator& operator= (
const PetscInverseOperator& );
353 const OperatorType &op_;
360 mutable int iterations_;
361 std::string solverName_;
362 std::string precondName_;
371 #endif // #if HAVE_PETSC 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
static bool verbose()
obtain the cached value for fem.verbose
Definition: io/parameter.hh:444