diff --git a/opm/autodiff/CPRPreconditioner.hpp b/opm/autodiff/CPRPreconditioner.hpp index 84afac7f5..27e81592e 100644 --- a/opm/autodiff/CPRPreconditioner.hpp +++ b/opm/autodiff/CPRPreconditioner.hpp @@ -52,7 +52,8 @@ #include namespace Opm { -namespace + +namespace ISTLUtility { /// /// \brief A traits class for selecting the types of the preconditioner. @@ -85,6 +86,7 @@ struct CPRSelector { return new Operator(m); } + }; #if HAVE_MPI @@ -266,8 +268,49 @@ createEllipticPreconditionerPointer(const M& Ae, double relax, ::EllipticPreconditionerPointer EllipticPreconditionerPointer; return EllipticPreconditionerPointer(new ParallelPreconditioner(Ae, comm, relax)); } -#endif -} // end namespace + +#endif // #if HAVE_MPI + +/// \brief Creates the elliptic preconditioner (ILU0) +/// \param opA The operator representing the matrix of the system. +/// \param relax The relaxation parameter for ILU0. +/// \param comm The object describing the parallelization information and communication. +// \param amgPtr The unique_ptr to be filled (return) +template +inline void +createAMGPreconditionerPointer( Op& opA, const double relax, const P& comm, std::unique_ptr< AMG >& amgPtr ) +{ + // type of matrix + typedef typename Op::matrix_type M; + + // The coupling metric used in the AMG + typedef Dune::Amg::FirstDiagonal CouplingMetric; + + // The coupling criterion used in the AMG + typedef Dune::Amg::SymmetricCriterion CritBase; + + // The coarsening criterion used in the AMG + typedef Dune::Amg::CoarsenCriterion Criterion; + + // TODO: revise choice of parameters + int coarsenTarget=1200; + Criterion criterion(15,coarsenTarget); + criterion.setDebugLevel( 0 ); // no debug information, 1 for printing hierarchy information + criterion.setDefaultValuesIsotropic(2); + criterion.setNoPostSmoothSteps( 1 ); + criterion.setNoPreSmoothSteps( 1 ); + + // for DUNE 2.2 we also need to pass the smoother args + typedef typename AMG::Smoother Smoother; + typedef typename Dune::Amg::SmootherTraits::Arguments SmootherArgs; + SmootherArgs smootherArgs; + smootherArgs.iterations = 1; + smootherArgs.relaxationFactor = relax; + + amgPtr.reset( new AMG(opA, criterion, smootherArgs, comm ) ); +} + +} // end namespace ISTLUtility struct CPRParameter { @@ -348,19 +391,21 @@ createEllipticPreconditionerPointer(const M& Ae, double relax, Dune::SolverCategory::sequential:Dune::SolverCategory::overlapping }; + typedef ISTLUtility::CPRSelector CPRSelectorType ; + //! \brief Elliptic Operator - typedef typename CPRSelector::Operator Operator; + typedef typename CPRSelectorType::Operator Operator; //! \brief preconditioner for the whole system (here either ILU(0) or ILU(n) typedef Dune::Preconditioner WholeSystemPreconditioner; //! \brief type of the unique pointer to the ilu-0 preconditioner //! used the for the elliptic system - typedef typename CPRSelector::EllipticPreconditionerPointer + typedef typename CPRSelectorType::EllipticPreconditionerPointer EllipticPreconditionerPointer; //! \brief amg preconditioner for the elliptic system - typedef typename CPRSelector::AMG AMG; + typedef typename CPRSelectorType::AMG AMG; /*! \brief Constructor. @@ -382,7 +427,7 @@ createEllipticPreconditionerPointer(const M& Ae, double relax, de_( Ae_.N() ), ve_( Ae_.M() ), dmodified_( A_.N() ), - opAe_(CPRSelector::makeOperator(Ae_, commAe)), + opAe_(CPRSelectorType::makeOperator(Ae_, commAe)), precond_(), // ilu0 preconditioner for elliptic system amg_(), // amg preconditioner for elliptic system pre_(), // copy A will be made be the preconditioner @@ -395,10 +440,10 @@ createEllipticPreconditionerPointer(const M& Ae, double relax, // create the preconditioner for the whole system. if( param_.cpr_ilu_n_ == 0 ) { - pre_ = createILU0Ptr( A_, param_.cpr_relax_, comm_ ); + pre_ = ISTLUtility::createILU0Ptr( A_, param_.cpr_relax_, comm_ ); } else { - pre_ = createILUnPtr( A_, param_.cpr_ilu_n_, param_.cpr_relax_, comm_ ); + pre_ = ISTLUtility::createILUnPtr( A_, param_.cpr_ilu_n_, param_.cpr_relax_, comm_ ); } } @@ -552,35 +597,11 @@ createEllipticPreconditionerPointer(const M& Ae, double relax, { if( amg ) { - // The coupling metric used in the AMG - typedef Dune::Amg::FirstDiagonal CouplingMetric; - - // The coupling criterion used in the AMG - typedef Dune::Amg::SymmetricCriterion CritBase; - - // The coarsening criterion used in the AMG - typedef Dune::Amg::CoarsenCriterion Criterion; - - // TODO: revise choice of parameters - int coarsenTarget=1200; - Criterion criterion(15,coarsenTarget); - criterion.setDebugLevel( 0 ); // no debug information, 1 for printing hierarchy information - criterion.setDefaultValuesIsotropic(2); - criterion.setNoPostSmoothSteps( 1 ); - criterion.setNoPreSmoothSteps( 1 ); - - // for DUNE 2.2 we also need to pass the smoother args - typedef typename AMG::Smoother Smoother; - typedef typename Dune::Amg::SmootherTraits::Arguments SmootherArgs; - SmootherArgs smootherArgs; - smootherArgs.iterations = 1; - smootherArgs.relaxationFactor = param_.cpr_relax_; - - amg_ = std::unique_ptr< AMG > (new AMG(*opAe_, criterion, smootherArgs, comm )); + ISTLUtility::createAMGPreconditionerPointer( *opAe_ , param_.cpr_relax_, comm, amg_ ); } else { - precond_ = createEllipticPreconditionerPointer( Ae_, param_.cpr_relax_, comm); + precond_ = ISTLUtility::createEllipticPreconditionerPointer( Ae_, param_.cpr_relax_, comm); } } };