DUNE-FEM (unstable)

implicit.hh
1#ifndef DUNE_FEM_SOLVER_RUNGEKUTTA_IMPLICIT_HH
2#define DUNE_FEM_SOLVER_RUNGEKUTTA_IMPLICIT_HH
3
4//- system includes
5#include <sstream>
6#include <vector>
7
8//- dune-common includes
10
11//- dune-fem includes
12#include <dune/fem/solver/rungekutta/basicimplicit.hh>
13#include <dune/fem/solver/rungekutta/butchertable.hh>
14#include <dune/fem/solver/rungekutta/timestepcontrol.hh>
15
16namespace DuneODE
17{
18
19 // ImplicitRungeKuttaSolver
20 // ------------------------
21
23 template< class HelmholtzOperator, class NonlinearSolver, class TimeStepControl = ImplicitRungeKuttaTimeStepControl >
25 : public BasicImplicitRungeKuttaSolver< HelmholtzOperator, NonlinearSolver, TimeStepControl >
26 {
29
30 public:
31 typedef HelmholtzOperator HelmholtzOperatorType;
32 typedef typename BaseType::TimeStepControlType TimeStepControlType;
33
34 typedef typename TimeStepControlType::TimeProviderType TimeProviderType;
35 typedef typename BaseType::ParameterType ParameterType;
36 typedef typename BaseType::NonlinearSolverParameterType NonlinearSolverParameterType;
37
45 ImplicitRungeKuttaSolver ( HelmholtzOperatorType &helmholtzOp,
46 TimeProviderType &timeProvider,
47 const SimpleButcherTable< double >& butcherTable,
48 const Dune::Fem::ParameterReader &parameter = Dune::Fem::Parameter::container() )
49 : BaseType( helmholtzOp,
50 butcherTable,
51 TimeStepControlType( timeProvider, ParameterType( parameter ) ),
52 NonlinearSolverParameterType( parameter )
53 )
54 {}
55
64 ImplicitRungeKuttaSolver ( HelmholtzOperatorType &helmholtzOp,
65 TimeProviderType &timeProvider,
66 int order,
67 const ParameterType& tscParam,
68 const NonlinearSolverParameterType& nlsParam )
69 : BaseType( helmholtzOp,
70 defaultButcherTables( tscParam.selectedSolver( order ) ),
71 TimeStepControlType( timeProvider, tscParam ),
72 nlsParam )
73 {}
74
83 ImplicitRungeKuttaSolver ( HelmholtzOperatorType &helmholtzOp,
84 TimeProviderType &timeProvider,
85 int order,
86 const Dune::Fem::ParameterReader &parameter = Dune::Fem::Parameter::container() )
87 : BaseType( helmholtzOp,
88 defaultButcherTables( ParameterType( parameter ).selectedSolver( order ) ),
89 TimeStepControlType( timeProvider, ParameterType( parameter ) ),
90 NonlinearSolverParameterType( parameter ) )
91 {}
92
100 ImplicitRungeKuttaSolver ( HelmholtzOperatorType &helmholtzOp,
101 TimeProviderType &timeProvider,
102 const ParameterType& tscParam,
103 const NonlinearSolverParameterType& nlsParam )
104 : BaseType( helmholtzOp,
105 defaultButcherTables( tscParam.selectedSolver( 1 ) ),
106 TimeStepControlType( timeProvider, tscParam ),
107 nlsParam )
108 {}
109
116 ImplicitRungeKuttaSolver ( HelmholtzOperatorType &helmholtzOp,
117 TimeProviderType &timeProvider,
118 const Dune::Fem::ParameterReader &parameter = Dune::Fem::Parameter::container() )
119 : BaseType( helmholtzOp,
120 defaultButcherTables( ParameterType( parameter ).selectedSolver( 1 ) ),
121 TimeStepControlType( timeProvider, ParameterType( parameter ) ),
122 NonlinearSolverParameterType( parameter ) )
123 {}
124
125 protected:
126 static SimpleButcherTable< double > defaultButcherTables ( const int solverId )
127 {
128 switch( solverId )
129 {
130 case 1:
131 return implicitEulerButcherTable();
132 case 2:
133 return gauss2ButcherTable();
134 case 3:
135 return implicit3ButcherTable();
136 case 4:
137 return implicit34ButcherTable();
138 default:
139 DUNE_THROW( NotImplemented, "Implicit Runge-Kutta method with id " << solverId << " not implemented." );
140 }
141 }
142 };
143
144} // namespace DuneODE
145
146#endif // #ifndef DUNE_FEM_SOLVER_RUNGEKUTTA_IMPLICIT_HH
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 &parameter=Dune::Fem::Parameter::container())
constructor
Definition: implicit.hh:116
ImplicitRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, TimeProviderType &timeProvider, const SimpleButcherTable< double > &butcherTable, const Dune::Fem::ParameterReader &parameter=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 &parameter=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
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 27, 22:29, 2024)