dune-fem  2.4.1-rc
semiimplicit.hh
Go to the documentation of this file.
1 #ifndef DUNE_FEM_SOLVER_RUNGEKUTTA_SEMIIMPLICIT_HH
2 #define DUNE_FEM_SOLVER_RUNGEKUTTA_SEMIIMPLICIT_HH
3 
4 //- system includes
5 #include <sstream>
6 #include <vector>
7 
8 //- dune-common includes
9 #include <dune/common/exceptions.hh>
10 
11 //- dune-fem includes
15 
16 namespace DuneODE
17 {
18 
19  // SemiImplicitRungeKuttaSourceTerm
20  // --------------------------------
21 
22  template< class ExplicitOperator >
24  {
26 
27  public:
28  typedef ExplicitOperator ExplicitOperatorType;
29 
30  typedef typename ExplicitOperatorType::DestinationType DestinationType;
31 
32  template< class ButcherTable >
33  SemiImplicitRungeKuttaSourceTerm ( ExplicitOperatorType &explicitOp,
34  const ButcherTable &butcherTable,
35  const Dune::DynamicMatrix< double > &implicitA )
36  : explicitOp_( explicitOp ),
37  alpha_( butcherTable.A() ),
38  gamma_( butcherTable.stages() ),
39  c_( butcherTable.c() ),
40  uex_( "SIRK u-explicit", explicitOp_.space() ),
41  limiter_( explicitOp_.hasLimiter() )
42  {
43  Dune::DynamicMatrix< double > Ainv( implicitA );
44  Ainv.invert();
45  alpha_.rightmultiply( Ainv );
46 
47  for( int i = 0; i < butcherTable.stages(); ++i )
48  {
49  gamma_[ i ] = 1.0;
50  for( int j = 0; j < i; ++j )
51  gamma_[ i ] -= alpha_[ i ][ j ];
52  }
53  }
54 
55  bool operator() ( double time, double timeStepSize, int stage,
56  const DestinationType &u, const std::vector< DestinationType * > &update,
57  DestinationType &source )
58  {
59  uex_.assign( u );
60  uex_ *= gamma_[ stage ];
61  for( int k = 0; k < stage; ++k )
62  uex_.axpy( alpha_[ stage ][ k ], *update[ k ] );
63  explicitOp_.setTime( time + c_[ stage ]*timeStepSize );
64  explicitOp_( uex_, source );
65 
66  return true;
67  }
68 
69  void limit( DestinationType& update, const double time )
70  {
71  if( limiter_ )
72  {
73  // set correct time
74  explicitOp_.setTime( time );
75  // copy given function
76  uex_.assign( update );
77  // apply limiter
78  explicitOp_.limit( uex_, update );
79  }
80  }
81 
82  double initialTimeStepEstimate ( double time, const DestinationType &u ) const
83  {
84  explicitOp_.setTime( time );
85  explicitOp_.initializeTimeStepSize( u );
86  return explicitOp_.timeStepEstimate();
87  }
88 
89  double timeStepEstimate () const
90  {
91  return explicitOp_.timeStepEstimate();
92  }
93 
94  private:
95  ExplicitOperatorType &explicitOp_;
96  Dune::DynamicMatrix< double > alpha_;
97  Dune::DynamicVector< double > gamma_, c_;
98  DestinationType uex_;
99  const bool limiter_ ;
100  };
101 
102 
103 
104  // SemiImplicitRungeKuttaSolver
105  // ----------------------------
106 
108  template< class ExplicitOperator, class HelmholtzOperator, class NonlinearSolver >
110  : public BasicImplicitRungeKuttaSolver< HelmholtzOperator, NonlinearSolver, ImplicitRungeKuttaTimeStepControl, SemiImplicitRungeKuttaSourceTerm< ExplicitOperator > >
111  {
114 
115  public:
116  typedef ExplicitOperator ExplicitOperatorType;
117  typedef HelmholtzOperator HelmholtzOperatorType;
120 
124 
134  SemiImplicitRungeKuttaSolver ( ExplicitOperatorType &explicitOp,
135  HelmholtzOperatorType &helmholtzOp,
136  TimeProviderType &timeProvider, int order,
137  const ParametersType& tscParams,
138  const NonlinearSolverParametersType& nlsParams )
139  : BaseType( helmholtzOp,
140  butcherTable( order, false ),
141  TimeStepControlType( timeProvider, tscParams ),
142  SourceTermType( explicitOp, butcherTable( order, true ), butcherTable( order, false ).A() ),
143  nlsParams )
144  {}
145 
146  SemiImplicitRungeKuttaSolver ( ExplicitOperatorType &explicitOp,
147  HelmholtzOperatorType &helmholtzOp,
148  TimeProviderType &timeProvider,
149  int order,
151  : BaseType( helmholtzOp,
152  butcherTable( order, false ),
153  TimeStepControlType( timeProvider, parameter ),
154  SourceTermType( explicitOp, butcherTable( order, true ), butcherTable( order, false ).A() ),
155  NonlinearSolverParametersType( parameter ) )
156  {}
157 
166  SemiImplicitRungeKuttaSolver ( ExplicitOperatorType &explicitOp,
167  HelmholtzOperatorType &helmholtzOp,
168  TimeProviderType &timeProvider,
169  const ParametersType& tscParams,
170  const NonlinearSolverParametersType& nlsParams )
171  : BaseType( helmholtzOp,
172  butcherTable( 1, false ),
173  TimeStepControlType( timeProvider, tscParams ),
174  SourceTermType( explicitOp, butcherTable( 1, true ), butcherTable( 1, false ).A() ),
175  nlsParams )
176  {}
177 
178  SemiImplicitRungeKuttaSolver ( ExplicitOperatorType &explicitOp,
179  HelmholtzOperatorType &helmholtzOp,
180  TimeProviderType &timeProvider,
182  : BaseType( helmholtzOp,
183  butcherTable( 1, false ),
184  TimeStepControlType( timeProvider, parameter ),
185  SourceTermType( explicitOp, butcherTable( 1, true ), butcherTable( 1, false ).A() ),
186  NonlinearSolverParametersType( parameter ) )
187  {}
188  protected:
189  static SimpleButcherTable< double > butcherTable ( int order, bool expl )
190  {
191  switch( order )
192  {
193  case 1:
194  return semiImplicitEulerButcherTable( expl );
195  case 2:
196  return semiImplicitSSP222ButcherTable( expl );
197  case 3:
198  return semiImplicit33ButcherTable( expl );
199  case 4:
200  return semiImplicitIERK45ButcherTable( expl );
201  default:
202  DUNE_THROW( NotImplemented, "Semi-implicit Runge-Kutta method of order " << order << " not implemented." );
203  }
204  }
205  };
206 
207 } // namespace DuneODE
208 
209 #endif // #ifndef DUNE_FEM_SOLVER_RUNGEKUTTA_IMPLICIT_HH
static SimpleButcherTable< double > butcherTable(int order, bool expl)
Definition: semiimplicit.hh:189
ExplicitOperator ExplicitOperatorType
Definition: semiimplicit.hh:116
BaseType::SourceTermType SourceTermType
Definition: semiimplicit.hh:119
Definition: butchertable.hh:15
SemiImplicitRungeKuttaSolver(ExplicitOperatorType &explicitOp, HelmholtzOperatorType &helmholtzOp, TimeProviderType &timeProvider, const Dune::Fem::ParameterReader &parameter=Dune::Fem::Parameter::container())
Definition: semiimplicit.hh:178
BaseType::ParametersType ParametersType
Definition: semiimplicit.hh:122
SemiImplicitRungeKuttaSolver(ExplicitOperatorType &explicitOp, HelmholtzOperatorType &helmholtzOp, TimeProviderType &timeProvider, const ParametersType &tscParams, const NonlinearSolverParametersType &nlsParams)
constructor
Definition: semiimplicit.hh:166
Definition: semiimplicit.hh:23
SimpleButcherTable< double > semiImplicitEulerButcherTable(bool expl)
Definition: butchertable.cc:93
Implicit RungeKutta ODE solver.
Definition: basicimplicit.hh:52
double timeStepEstimate() const
Definition: semiimplicit.hh:89
Definition: multistep.hh:16
SimpleButcherTable< double > semiImplicitSSP222ButcherTable(bool expl)
Definition: butchertable.cc:189
SemiImplicitRungeKuttaSourceTerm(ExplicitOperatorType &explicitOp, const ButcherTable &butcherTable, const Dune::DynamicMatrix< double > &implicitA)
Definition: semiimplicit.hh:33
SimpleButcherTable< double > semiImplicit33ButcherTable(bool expl)
Definition: butchertable.cc:159
TimeStepControlType::TimeProviderType TimeProviderType
Definition: semiimplicit.hh:121
general base for time providers
Definition: timeprovider.hh:35
Definition: timestepcontrol.hh:165
Implicit RungeKutta ODE solver.
Definition: semiimplicit.hh:109
ExplicitOperator ExplicitOperatorType
Definition: semiimplicit.hh:28
static ParameterContainer & container()
Definition: io/parameter.hh:190
void limit(DestinationType &update, const double time)
Definition: semiimplicit.hh:69
SemiImplicitRungeKuttaSolver(ExplicitOperatorType &explicitOp, HelmholtzOperatorType &helmholtzOp, TimeProviderType &timeProvider, int order, const Dune::Fem::ParameterReader &parameter=Dune::Fem::Parameter::container())
Definition: semiimplicit.hh:146
ExplicitOperatorType::DestinationType DestinationType
Definition: semiimplicit.hh:30
BaseType::NonlinearSolverParametersType NonlinearSolverParametersType
Definition: semiimplicit.hh:123
double initialTimeStepEstimate(double time, const DestinationType &u) const
Definition: semiimplicit.hh:82
HelmholtzOperator HelmholtzOperatorType
Definition: semiimplicit.hh:117
SemiImplicitRungeKuttaSolver(ExplicitOperatorType &explicitOp, HelmholtzOperatorType &helmholtzOp, TimeProviderType &timeProvider, int order, const ParametersType &tscParams, const NonlinearSolverParametersType &nlsParams)
constructor
Definition: semiimplicit.hh:134
Definition: timestepcontrol.hh:21
bool operator()(double time, double timeStepSize, int stage, const DestinationType &u, const std::vector< DestinationType * > &update, DestinationType &source)
Definition: semiimplicit.hh:55
BaseType::TimeStepControlType TimeStepControlType
Definition: semiimplicit.hh:118
SimpleButcherTable< double > semiImplicitIERK45ButcherTable(bool expl)
Definition: butchertable.cc:326