1#ifndef DUNE_FEM_SOLVERPARAMETER_HH 
    2#define DUNE_FEM_SOLVERPARAMETER_HH 
    4#include <dune/fem/io/parameter.hh> 
    5#include <dune/fem/common/staticlistofint.hh> 
   12    namespace LinearSolver
 
   14      struct ToleranceCriteria {
 
   16        static const int relative = 1;
 
   17        static const int residualReduction = 2;
 
   21    struct SolverParameter
 
   23    : 
public LocalParameter< SolverParameter, SolverParameter >
 
   28      const std::string keyPrefix_;
 
   30      ParameterReader parameter_;
 
   34      LIST_OF_INT_FORWARDED(Solvers,
 
   45      static const std::string solverMethodTable(
int id)
 
   47        return Solvers::to_string( 
id );
 
   53      LIST_OF_INT_FORWARDED(Preconditioners,
 
   65      static const std::string preconditionMethodTable(
int id)
 
   67        return Preconditioners::to_string( 
id );
 
   70      SolverParameter ( 
const ParameterReader ¶meter = Parameter::container() )
 
   71        : keyPrefix_( 
"fem.solver." ), parameter_( parameter )
 
   74      explicit SolverParameter ( 
const std::string keyPrefix, 
const ParameterReader ¶meter = Parameter::container() )
 
   75        : keyPrefix_( keyPrefix ), parameter_( parameter )
 
   78      const std::string& keyPrefix()
 const { 
return keyPrefix_; }
 
   80      const ParameterReader& parameter()
 const { 
return parameter_; }
 
   90      virtual bool verbose()
 const 
   94          verbose_ = parameter_.getValue< 
bool >( keyPrefix_ + 
"verbose", false ) ? 1 : 0;
 
   96        return bool(verbose_);
 
   99      virtual void setVerbose( 
const bool verb )
 
  101        verbose_ = verb ? 1 : 0;
 
  106      int defaultErrorMeasure = 0;
 
  107      virtual void setDefaultErrorMeasure(
int def) { defaultErrorMeasure = def; }
 
  108      virtual int errorMeasure()
 const 
  110        const std::string errorTypeTable[] =
 
  111        { 
"absolute", 
"relative", 
"residualreduction" };
 
  112        const int errorType = parameter_.getEnum( keyPrefix_ + 
"errormeasure", errorTypeTable,
 
  113                              defaultErrorMeasure );
 
  117      virtual double tolerance (  )
 const 
  121          double defaultTol = 1e-8;
 
  122          if(parameter_.exists(keyPrefix_ + 
"tolerance"))
 
  123            tolerance_ = parameter_.getValue< 
double >( keyPrefix_ + 
"tolerance", defaultTol );
 
  126            if( parameter_.exists(keyPrefix_ + 
"absolutetol") ||
 
  127                parameter_.exists(keyPrefix_ + 
"reductiontol") )
 
  130              int measure = errorMeasure();
 
  132                tolerance_ = absoluteTol__();
 
  134                tolerance_ = reductionTol__();
 
  136            else tolerance_ = defaultTol; 
 
  142      virtual void setTolerance ( 
const double eps )
 
  144        assert( eps >= 0.0 );
 
  148      virtual int maxIterations ()
 const 
  150        if( maxIterations_ < 0 )
 
  152          if(parameter_.exists(keyPrefix_ + 
"maxlineariterations"))
 
  154            std::cout << 
"WARNING: Parameter " + keyPrefix_ +
"maxlineariterations is deprecated. Please use " + keyPrefix_ + 
"maxiterations instead." << std::endl;
 
  155            maxIterations_ =  parameter_.getValue< 
double >(keyPrefix_ + 
"maxlineariterations");
 
  160        return maxIterations_;
 
  163      virtual void  setMaxIterations ( 
const int maxIter )
 
  165        assert( maxIter >= 0 );
 
  166        maxIterations_ = maxIter;
 
  169      virtual int solverMethod(
 
  170            const std::vector<int> standardMethods,
 
  171            const std::vector<std::string> &additionalMethods = {},
 
  172            int defaultMethod = 0   
 
  175        std::vector<std::string> methodTable(standardMethods.size()+additionalMethods.size());
 
  176        for (std::size_t i=0;i<standardMethods.size();++i)
 
  177          methodTable[i] = solverMethodTable(standardMethods[i]);
 
  178        for (std::size_t i=0;i<additionalMethods.size();++i)
 
  179          methodTable[standardMethods.size()+i] = additionalMethods[i];
 
  181        if( parameter_.exists( keyPrefix_ + 
"method" ) ||
 
  182           !parameter_.exists( 
"method" ) )
 
  183          method = parameter_.getEnum( keyPrefix_ + 
"method", methodTable, defaultMethod );
 
  186          method = parameter_.getEnum( 
"krylovmethod", methodTable, defaultMethod );
 
  187          std::cout << 
"WARNING: using old parameter name 'krylovmethod' " 
  188                    << 
"please switch to '" << keyPrefix_ << 
"method'\n";
 
  190        if (method < standardMethods.size())
 
  191          return standardMethods[method];
 
  193          return -(method-standardMethods.size()); 
 
  196      virtual int gmresRestart()
 const 
  198        int defaultRestart = 20;
 
  199        return parameter_.getValue< 
int >( keyPrefix_ + 
"gmres.restart", defaultRestart );
 
  202      virtual int preconditionMethod(
 
  203            const std::vector<int> standardMethods,
 
  204            const std::vector<std::string> &additionalMethods = {},
 
  205            int defaultMethod = 0   
 
  208        std::vector<std::string> methodTable(standardMethods.size()+additionalMethods.size());
 
  209        for (std::size_t i=0;i<standardMethods.size();++i)
 
  210          methodTable[i] = preconditionMethodTable(standardMethods[i]);
 
  211        for (std::size_t i=0;i<additionalMethods.size();++i)
 
  212          methodTable[standardMethods.size()+i] = additionalMethods[i];
 
  213        std::size_t method = parameter_.getEnum( keyPrefix_ + 
"preconditioning.method", methodTable, defaultMethod );
 
  214        if (method < standardMethods.size())
 
  215          return standardMethods[method];
 
  217          return -(method-standardMethods.size()); 
 
  220      virtual double relaxation ()
 const 
  222        return parameter_.getValue< 
double >( keyPrefix_ + 
"preconditioning.relaxation", 1.1 );
 
  225      virtual int preconditionerIteration ()
 const 
  227        return parameter_.getValue< 
int >( keyPrefix_ + 
"preconditioning.iterations", 1 );
 
  230      virtual int preconditionerLevel ()
 const 
  232        return parameter_.getValue< 
int >( keyPrefix_ + 
"preconditioning.level", 0 );
 
  235      virtual bool threading ()
 const 
  237        return parameter_.getValue< 
bool >( keyPrefix_ + 
"threading", true );
 
  240      virtual bool knollTrick()
 const 
  242        return parameter_.getValue< 
bool >( keyPrefix_ + 
"knolltrick", false );
 
  246      virtual double absoluteTol__ ( )
  const 
  248        if( absoluteTol_ < 0 )
 
  250          if(parameter_.exists(keyPrefix_ + 
"linabstol"))
 
  252            std::cout << 
"WARNING: Parameter " + keyPrefix_ + 
"linabstol is deprecated. Please use " + keyPrefix_ + 
"absolutetol instead." << std::endl;
 
  253            absoluteTol_ =  parameter_.getValue< 
double >(keyPrefix_ + 
"linabstol");
 
  257            std::cout << 
"WARNING: Parameter " + keyPrefix_ + 
"absolutetol is deprecated. Please use " + keyPrefix_ + 
"tolerance instead." << std::endl;
 
  258            absoluteTol_ =  parameter_.getValue< 
double >(keyPrefix_ +  
"absolutetol", 1e-8 );
 
  264      virtual double reductionTol__ (  )
 const 
  266        if( reductionTol_ < 0 )
 
  268          if(parameter_.exists(keyPrefix_ + 
"linreduction"))
 
  270            std::cout << 
"WARNING: Parameter " + keyPrefix_ +
"linreduction is deprecated. Please use " + keyPrefix_ + 
"reductiontol instead." << std::endl;
 
  271            reductionTol_ =  parameter_.getValue< 
double >(keyPrefix_ + 
"linreduction");
 
  275            std::cout << 
"WARNING: Parameter " + keyPrefix_ +
"reductiontol is deprecated. Please use " + keyPrefix_ + 
"tolerance instead." << std::endl;
 
  276            reductionTol_ = parameter_.getValue< 
double >( keyPrefix_ + 
"reductiontol", 1e-8 );
 
  279        return reductionTol_;
 
  282      mutable int    verbose_       = -1;
 
  283      mutable int    maxIterations_ = -1;
 
  284      mutable double absoluteTol_   = -1.;
 
  285      mutable double reductionTol_  = -1.;
 
  286      mutable double tolerance_     = -1.;
 
@ absolute
|a-b| <= epsilon
Definition: float_cmp.hh:110
 
constexpr GeometryType none(unsigned int dim)
Returns a GeometryType representing a singular of dimension dim.
Definition: type.hh:471
 
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:485
 
Dune namespace.
Definition: alignedallocator.hh:13