1 #ifndef DUNE_FEM_SOLVER_RUNGEKUTTA_BASICIMPLICIT_HH 2 #define DUNE_FEM_SOLVER_RUNGEKUTTA_BASICIMPLICIT_HH 11 #include <dune/common/dynmatrix.hh> 12 #include <dune/common/dynvector.hh> 13 #include <dune/common/nullptr.hh> 27 bool operator() (
double time,
double timeStepSize,
int stage,
const T &u,
const std::vector< T * > &update, T &source )
33 void limit( T& update,
const double time ) {}
51 template<
class HelmholtzOperator,
class NonlinearSolver,
class TimeStepControl,
class SourceTerm = NoImplicitRungeKuttaSourceTerm >
79 template<
class ButcherTable >
81 const ButcherTable &butcherTable,
82 const TimeStepControlType &timeStepControl,
83 const SourceTermType &sourceTerm,
84 const NonlinearSolverParametersType& parameters )
85 : helmholtzOp_( helmholtzOp ),
86 nonlinearSolver_( helmholtzOp_, parameters ),
87 timeStepControl_( timeStepControl ),
88 sourceTerm_( sourceTerm ),
89 stages_( butcherTable.stages() ),
90 alpha_( butcherTable.A() ),
93 c_( butcherTable.c() ),
94 rhs_(
"RK rhs", helmholtzOp_.space() ),
95 update_( stages(), nullptr )
97 setup( butcherTable );
106 template<
class ButcherTable >
108 const ButcherTable &butcherTable,
109 const TimeStepControlType &timeStepControl,
110 const NonlinearSolverParametersType& parameters )
111 : helmholtzOp_( helmholtzOp ),
112 nonlinearSolver_( helmholtzOp_, parameters ),
113 timeStepControl_( timeStepControl ),
114 stages_( butcherTable.stages() ),
115 alpha_( butcherTable.A() ),
118 c_( butcherTable.c() ),
119 rhs_(
"RK rhs", helmholtzOp_.space() ),
120 update_( stages(), nullptr )
122 setup( butcherTable );
125 template<
class ButcherTable >
127 const ButcherTable &butcherTable,
128 const TimeStepControlType &timeStepControl,
130 : helmholtzOp_( helmholtzOp ),
131 nonlinearSolver_( helmholtzOp_, parameter ),
132 timeStepControl_( timeStepControl ),
133 stages_( butcherTable.stages() ),
134 alpha_( butcherTable.A() ),
137 c_( butcherTable.c() ),
138 rhs_(
"RK rhs", helmholtzOp_.space() ),
139 update_( stages(), nullptr )
141 setup( butcherTable );
144 template<
class ButcherTable >
146 const ButcherTable &butcherTable,
148 : helmholtzOp_( helmholtzOp ),
149 nonlinearSolver_( helmholtzOp_, parameter ),
150 timeStepControl_( parameter ),
151 stages_( butcherTable.stages() ),
152 alpha_( butcherTable.A() ),
155 c_( butcherTable.c() ),
156 rhs_(
"RK rhs", helmholtzOp_.space() ),
157 update_( stages(), nullptr )
159 setup( butcherTable );
163 template<
class ButcherTable >
164 void setup(
const ButcherTable& butcherTable )
167 for(
int i = 0; i < stages(); ++i )
169 std::ostringstream name;
170 name <<
"RK stage " << i;
171 update_[ i ] =
new DestinationType( name.str(), helmholtzOp_.space() );
175 Dune::DynamicMatrix< double > AL( alpha_ );
176 for(
int i = 0; i < stages(); ++i )
178 gamma_[ i ] = AL[ i ][ i ];
183 alpha_.mtv( butcherTable.b(), beta_ );
185 alpha_.leftmultiply( AL );
186 for(
int i = 0; i < stages(); ++i )
187 alpha_[ i ][ i ] = gamma_[ i ];
189 for(
int i = 0; i < stages(); ++i )
192 for(
int j = 0; j < i; ++j )
193 gamma_[ i ] -= alpha_[ i ][ j ];
197 for(
int i = 0; i < stages(); ++i )
198 delta_ -= beta_[ i ];
207 for(
int i = 0; i < stages(); ++i )
214 const double time = timeStepControl_.time();
216 helmholtzOp_.setTime( time );
217 helmholtzOp_.initializeTimeStepSize( U0 );
218 const double helmholtzEstimate = helmholtzOp_.timeStepEstimate();
220 double sourceTermEstimate = sourceTerm_.initialTimeStepEstimate( time, U0 );
222 if( sourceTermEstimate < 0.0 ) sourceTermEstimate = helmholtzEstimate ;
224 timeStepControl_.initialTimeStepSize( helmholtzEstimate, sourceTermEstimate );
227 using BaseType::solve;
230 void solve ( DestinationType &U, MonitorType &monitor )
234 const double time = timeStepControl_.time();
235 const double timeStepSize = timeStepControl_.timeStepSize();
236 assert( timeStepSize > 0.0 );
237 for(
int s = 0; s < stages(); ++s )
240 DestinationType& updateStage = *update_[ s ];
243 updateStage.assign( U );
244 updateStage *= gamma_[ s ];
245 for(
int k = 0; k < s; ++k )
246 updateStage.axpy( alpha_[ s ][ k ], *update_[ k ] );
248 const double stageTime = time + c_[ s ]*timeStepSize;
249 if( sourceTerm_( time, timeStepSize, s, U, update_, rhs_ ) )
251 updateStage.axpy( alpha_[ s ][ s ]*timeStepSize, rhs_ );
252 sourceTerm_.limit( updateStage, stageTime );
256 helmholtzOp_.setTime( stageTime );
257 helmholtzOp_.setLambda( 0 );
258 helmholtzOp_( updateStage, rhs_ );
261 helmholtzOp_.setLambda( alpha_[ s ][ s ]*timeStepSize );
262 nonlinearSolver_( rhs_, updateStage );
265 monitor.newtonIterations_ += nonlinearSolver_.iterations();
266 monitor.linearSolverIterations_ += nonlinearSolver_.linearIterations();
269 if( !nonlinearSolver_.converged() )
270 return timeStepControl_.reduceTimeStep( helmholtzOp_.timeStepEstimate(), sourceTerm_.timeStepEstimate(), monitor );
274 if( timeStepControl_.computeError() )
277 DestinationType Uerr( U );
281 for(
int s = 0; s < stages(); ++s )
282 U.axpy( beta_[ s ], *update_[ s ] );
285 Uerr.axpy( -1.0, U );
286 const double errorU = Uerr.scalarProductDofs( Uerr );
287 const double normU = U.scalarProductDofs( U );
289 if( normU > 0 && errorU > 0 )
293 std::cout << std::scientific <<
"Error in RK = " << error <<
" norm " << errorU <<
" " << normU << std::endl;
300 for(
int s = 0; s < stages(); ++s )
301 U.axpy( beta_[ s ], *update_[ s ] );
304 monitor.error_ = error;
307 timeStepControl_.timeStepEstimate( helmholtzOp_.timeStepEstimate(), sourceTerm_.timeStepEstimate(), monitor );
314 out <<
"Generic " << stages() <<
"-stage implicit Runge-Kutta solver.\\\\" << std::endl;
318 double infNorm(
const DestinationType& U,
const DestinationType& Uerr )
const 320 typedef typename DestinationType :: ConstDofIteratorType ConstDofIteratorType ;
321 const ConstDofIteratorType uend = U.dend();
323 for( ConstDofIteratorType u = U.dbegin(), uerr = Uerr.dbegin(); u != uend; ++u, ++uerr )
326 double uerrval = *uerr ;
329 double norm =
std::abs( uval - uerrval );
345 Dune::DynamicVector< double >
gamma_, beta_, c_;
353 #endif // #ifndef DUNE_FEM_SOLVER_RUNGEKUTTA_BASICIMPLICIT_HH void description(std::ostream &out) const
print description of ODE solver to out stream
Definition: basicimplicit.hh:312
static double sqrt(const Double &v)
Definition: double.hh:870
BaseType::DestinationType DestinationType
Definition: basicimplicit.hh:60
Interface class for ODE Solver.
Definition: odesolverinterface.hh:11
SourceTerm sourceTerm_
Definition: basicimplicit.hh:340
DestinationType rhs_
Definition: basicimplicit.hh:347
NonlinearSolver::ParametersType NonlinearSolverParametersType
Definition: basicimplicit.hh:70
int stages_
Definition: basicimplicit.hh:342
BaseType::MonitorType MonitorType
Definition: basicimplicit.hh:59
static constexpr T max(T a)
Definition: utility.hh:65
void limit(T &update, const double time)
Definition: basicimplicit.hh:33
SourceTerm SourceTermType
Definition: basicimplicit.hh:65
NonlinearSolverType nonlinearSolver_
Definition: basicimplicit.hh:338
TimeStepControl TimeStepControlType
Definition: basicimplicit.hh:64
Implicit RungeKutta ODE solver.
Definition: basicimplicit.hh:52
Definition: multistep.hh:16
void initialize(const DestinationType &U0)
apply operator once to get dt estimate
Definition: basicimplicit.hh:212
double delta_
Definition: basicimplicit.hh:343
NonlinearSolver NonlinearSolverType
Definition: basicimplicit.hh:63
void setup(const ButcherTable &butcherTable)
Definition: basicimplicit.hh:164
Dune::DynamicVector< double > gamma_
Definition: basicimplicit.hh:345
Definition: basicimplicit.hh:24
~BasicImplicitRungeKuttaSolver()
destructor
Definition: basicimplicit.hh:205
general base for time providers
Definition: timeprovider.hh:35
Double abs(const Double &a)
Definition: double.hh:860
BasicImplicitRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, const ButcherTable &butcherTable, const TimeStepControlType &timeStepControl, const SourceTermType &sourceTerm, const NonlinearSolverParametersType ¶meters)
constructor
Definition: basicimplicit.hh:80
BasicImplicitRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, const ButcherTable &butcherTable, const TimeStepControlType &timeStepControl, const Dune::Fem::ParameterReader ¶meter=Dune::Fem::Parameter::container())
Definition: basicimplicit.hh:126
double initialTimeStepEstimate(double time, const T &u) const
Definition: basicimplicit.hh:36
Dune::Fem::TimeProviderBase TimeProviderType
Definition: basicimplicit.hh:67
bool operator()(double time, double timeStepSize, int stage, const T &u, const std::vector< T * > &update, T &source)
Definition: basicimplicit.hh:27
static ParameterContainer & container()
Definition: io/parameter.hh:190
TimeStepControl timeStepControl_
Definition: basicimplicit.hh:339
double infNorm(const DestinationType &U, const DestinationType &Uerr) const
Definition: basicimplicit.hh:318
HelmholtzOperatorType & helmholtzOp_
Definition: basicimplicit.hh:337
Definition: odesolverinterface.hh:17
int stages() const
Definition: basicimplicit.hh:310
std::vector< DestinationType * > update_
Definition: basicimplicit.hh:348
HelmholtzOperator HelmholtzOperatorType
Definition: basicimplicit.hh:62
Definition: timestepcontrol.hh:21
Dune::DynamicMatrix< double > alpha_
Definition: basicimplicit.hh:344
void solve(DestinationType &U, MonitorType &monitor)
solve the system
Definition: basicimplicit.hh:230
BasicImplicitRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, const ButcherTable &butcherTable, const Dune::Fem::ParameterReader ¶meter=Dune::Fem::Parameter::container())
Definition: basicimplicit.hh:145
BasicImplicitRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, const ButcherTable &butcherTable, const TimeStepControlType &timeStepControl, const NonlinearSolverParametersType ¶meters)
constructor
Definition: basicimplicit.hh:107
double timeStepEstimate() const
Definition: basicimplicit.hh:42
TimeStepControlType::ParametersType ParametersType
Definition: basicimplicit.hh:69
DestinationImp DestinationType
type of destination
Definition: odesolverinterface.hh:53