1#ifndef DUNE_FEM_SOLVER_RUNGEKUTTA_IMPLICIT_HH
2#define DUNE_FEM_SOLVER_RUNGEKUTTA_IMPLICIT_HH
12#include <dune/fem/solver/rungekutta/basicimplicit.hh>
13#include <dune/fem/solver/rungekutta/butchertable.hh>
14#include <dune/fem/solver/rungekutta/timestepcontrol.hh>
23 template<
class HelmholtzOperator,
class NonlinearSolver,
class TimeStepControl = ImplicitRungeKuttaTimeStepControl >
31 typedef HelmholtzOperator HelmholtzOperatorType;
32 typedef typename BaseType::TimeStepControlType TimeStepControlType;
34 typedef typename TimeStepControlType::TimeProviderType TimeProviderType;
35 typedef typename BaseType::ParameterType ParameterType;
36 typedef typename BaseType::NonlinearSolverParameterType NonlinearSolverParameterType;
46 TimeProviderType &timeProvider,
47 const SimpleButcherTable< double >& butcherTable,
48 const Dune::Fem::ParameterReader ¶meter = Dune::Fem::Parameter::container() )
51 TimeStepControlType( timeProvider, ParameterType( parameter ) ),
52 NonlinearSolverParameterType( parameter )
65 TimeProviderType &timeProvider,
67 const ParameterType& tscParam,
68 const NonlinearSolverParameterType& nlsParam )
70 defaultButcherTables( tscParam.selectedSolver( order ) ),
71 TimeStepControlType( timeProvider, tscParam ),
84 TimeProviderType &timeProvider,
86 const Dune::Fem::ParameterReader ¶meter = Dune::Fem::Parameter::container() )
88 defaultButcherTables( ParameterType( parameter ).selectedSolver( order ) ),
89 TimeStepControlType( timeProvider, ParameterType( parameter ) ),
90 NonlinearSolverParameterType( parameter ) )
101 TimeProviderType &timeProvider,
102 const ParameterType& tscParam,
103 const NonlinearSolverParameterType& nlsParam )
105 defaultButcherTables( tscParam.selectedSolver( 1 ) ),
106 TimeStepControlType( timeProvider, tscParam ),
117 TimeProviderType &timeProvider,
118 const Dune::Fem::ParameterReader ¶meter = Dune::Fem::Parameter::container() )
120 defaultButcherTables( ParameterType( parameter ).selectedSolver( 1 ) ),
121 TimeStepControlType( timeProvider, ParameterType( parameter ) ),
122 NonlinearSolverParameterType( parameter ) )
126 static SimpleButcherTable< double > defaultButcherTables (
const int solverId )
131 return implicitEulerButcherTable();
133 return gauss2ButcherTable();
135 return implicit3ButcherTable();
137 return implicit34ButcherTable();
Implicit RungeKutta ODE solver.
Definition: basicimplicit.hh:54
Implicit RungeKutta ODE solver.
Definition: implicit.hh:26
ImplicitRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, TimeProviderType &timeProvider, const Dune::Fem::ParameterReader ¶meter=Dune::Fem::Parameter::container())
constructor
Definition: implicit.hh:116
ImplicitRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, TimeProviderType &timeProvider, const SimpleButcherTable< double > &butcherTable, const Dune::Fem::ParameterReader ¶meter=Dune::Fem::Parameter::container())
constructor
Definition: implicit.hh:45
ImplicitRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, TimeProviderType &timeProvider, const ParameterType &tscParam, const NonlinearSolverParameterType &nlsParam)
constructor
Definition: implicit.hh:100
ImplicitRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, TimeProviderType &timeProvider, int order, const Dune::Fem::ParameterReader ¶meter=Dune::Fem::Parameter::container())
constructor
Definition: implicit.hh:83
ImplicitRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, TimeProviderType &timeProvider, int order, const ParameterType &tscParam, const NonlinearSolverParameterType &nlsParam)
constructor
Definition: implicit.hh:64
Default exception for dummy implementations.
Definition: exceptions.hh:263
A few common exception classes.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218