1#ifndef DUNE_FEM_SOLVER_RUNGEKUTTA_BASICROW_HH
2#define DUNE_FEM_SOLVER_RUNGEKUTTA_BASICROW_HH
18#include <dune/fem/solver/odesolver.hh>
19#include <dune/fem/solver/parameter.hh>
27 struct NoROWRungeKuttaSourceTerm
30 bool operator() (
double time,
double timeStepSize,
int stage,
const T &u,
const std::vector< T * > &update, T &source )
36 void limit( T& update,
const double time ) {}
39 double initialTimeStepEstimate (
double time,
const T &u )
const
45 double timeStepEstimate ()
const
53 template<
class HelmholtzOperator,
class NonlinearSolver,
class TimeStepControl,
class SourceTerm = NoROWRungeKuttaSourceTerm >
64 typedef HelmholtzOperator HelmholtzOperatorType;
65 typedef NonlinearSolver NonlinearSolverType;
66 typedef TimeStepControl TimeStepControlType;
67 typedef SourceTerm SourceTermType;
69 typedef typename NonlinearSolver::ParameterType NonlinearSolverParameterType;
70 typedef typename NonlinearSolverType::LinearInverseOperatorType LinearInverseOperatorType;
72 typedef typename HelmholtzOperator::SpaceOperatorType::PreconditionOperatorType PreconditionOperatorType;
83 template<
class ButcherTable >
86 const ButcherTable &butcherTable,
87 const TimeStepControlType &timeStepControl,
88 const SourceTermType &sourceTerm,
89 const NonlinearSolverParameterType& parameter )
90 : helmholtzOp_( helmholtzOp ),
91 linearSolver_( parameter.linear() ),
92 timeStepControl_( timeStepControl ),
93 sourceTerm_( sourceTerm ),
94 stages_( butcherTable.
stages() ),
95 ord_( butcherTable.
order() ),
96 alpha_( butcherTable.A() ),
97 alpha2_( butcherTable.B() ),
100 c_( butcherTable.c() ),
101 rhs_(
"RK rhs", helmholtzOp_.space() ),
102 temp_(
"RK temp", helmholtzOp_.space() ),
103 update_(
stages(), nullptr ),
104 maxLinearIterations_( parameter.linear().maxIterations() ),
105 preconditioner_(helmholtzOp.spaceOperator().preconditioner())
107 setup( butcherTable );
110 template<
class ButcherTable >
112 TimeProviderType& timeProvider,
113 const ButcherTable &butcherTable,
114 const TimeStepControlType &timeStepControl,
115 const SourceTermType &sourceTerm,
116 const Dune::Fem::ParameterReader ¶meter = Dune::Fem::Parameter::container() )
117 :
BasicROWRungeKuttaSolver( helmholtzOp, timeProvider, butcherTable, timeStepControl, sourceTerm, NonlinearSolverParameterType( parameter ) )
126 template<
class ButcherTable >
129 const ButcherTable &butcherTable,
130 const TimeStepControlType &timeStepControl,
131 const NonlinearSolverParameterType& parameter )
135 template<
class ButcherTable >
137 TimeProviderType& timeProvider,
138 const ButcherTable &butcherTable,
139 const TimeStepControlType &timeStepControl,
140 const Dune::Fem::ParameterReader ¶meter = Dune::Fem::Parameter::container() )
141 :
BasicROWRungeKuttaSolver( helmholtzOp, timeProvider, butcherTable, timeStepControl, SourceTermType(), NonlinearSolverParameterType( parameter ) )
149 template<
class ButcherTable >
152 const ButcherTable &butcherTable,
153 const NonlinearSolverParameterType& parameter )
154 :
BasicROWRungeKuttaSolver( helmholtzOp, timeProvider, butcherTable, TimeStepControlType(), SourceTermType(), parameter )
157 template<
class ButcherTable >
159 TimeProviderType& timeProvider,
160 const ButcherTable &butcherTable,
161 const Dune::Fem::ParameterReader ¶meter = Dune::Fem::Parameter::container() )
162 :
BasicROWRungeKuttaSolver( helmholtzOp, timeProvider, butcherTable, TimeStepControlType(), SourceTermType(), NonlinearSolverParameterType( parameter ) )
165 template<
class ButcherTable >
166 void setup(
const ButcherTable& butcherTable )
170 std::cout <<
"ROW method of order=" << butcherTable.
order() <<
" with " << stages_ <<
" stages" << std::endl;
174 for(
int i = 0; i <
stages(); ++i )
176 std::ostringstream name;
177 name <<
"RK stage " << i;
178 update_[ i ] =
new DestinationType( name.str(), helmholtzOp_.space() );
182 for(
int i = 0; i <
stages(); ++i )
183 gamma_[ i ] = alpha_[ i ][ i ];
186 alpha_.
mtv( butcherTable.b(), beta_ );
193 for(
int i = 0; i <
stages(); ++i )
200 const double time = timeStepControl_.time();
202 helmholtzOp_.setTime( time );
203 helmholtzOp_.initializeTimeStepSize( U0 );
204 const double helmholtzEstimate = helmholtzOp_.timeStepEstimate();
206 double sourceTermEstimate = sourceTerm_.initialTimeStepEstimate( time, U0 );
208 if( sourceTermEstimate < 0.0 ) sourceTermEstimate = helmholtzEstimate ;
210 timeStepControl_.initialTimeStepSize( helmholtzEstimate, sourceTermEstimate );
216 void solve ( DestinationType &U, MonitorType &monitor )
220 typename HelmholtzOperatorType::JacobianOperatorType jOp(
"jacobianOperator", U.space(), U.space() );
222 const double time = timeStepControl_.time();
223 const double timeStepSize = timeStepControl_.timeStepSize();
224 assert( timeStepSize > 0.0 );
226 for(
int s = 0; s <
stages(); ++s )
229 DestinationType& updateStage = *update_[ s ];
232 for(
int k = 0; k < s; ++k )
233 rhs_.axpy( alpha2_[ s ][ k ], *update_[ k ] );
234 helmholtzOp_.spaceOperator()(rhs_,updateStage);
235 updateStage *= (gamma_[s]*timeStepSize);
236 for(
int k = 0; k < s; ++k )
237 updateStage.axpy( -gamma_[s]*alpha_[ s ][ k ], *update_[ k ] );
239 rhs_.assign( updateStage );
242 const double stageTime = time + c_[ s ]*timeStepSize;
243 helmholtzOp_.setTime( stageTime );
248 helmholtzOp_.setLambda( gamma_[ s ]*timeStepSize );
249 helmholtzOp_.jacobian( U, jOp );
251 const int remLinearIts = maxLinearIterations_;
253 linearSolver_.parameter().setMaxIterations( remLinearIts );
256 linearSolver_.bind( jOp, *preconditioner_ );
258 linearSolver_.bind( jOp );
260 linearSolver_( rhs_, updateStage );
261 monitor.linearSolverIterations_ += linearSolver_.iterations();
263 linearSolver_.unbind();
267 if(0 && timeStepControl_.computeError() )
270 DestinationType Uerr( U );
274 for(
int s = 0; s <
stages(); ++s )
275 U.axpy( beta_[ s ], *update_[ s ] );
278 Uerr.axpy( -1.0, U );
279 const double errorU = Uerr.scalarProductDofs( Uerr );
280 const double normU = U.scalarProductDofs( U );
282 if( normU > 0 && errorU > 0 )
284 error = std::sqrt( errorU / normU );
286 std::cout << std::scientific <<
"Error in RK = " << error <<
" norm " << errorU <<
" " << normU << std::endl;
292 for(
int s = 0; s <
stages(); ++s )
293 U.axpy( beta_[ s ], *update_[ s ] );
296 monitor.error_ = error;
299 timeStepControl_.timeStepEstimate( helmholtzOp_.timeStepEstimate(), sourceTerm_.timeStepEstimate(), monitor );
308 out <<
"Generic " <<
stages() <<
"-stage implicit Runge-Kutta solver.\\\\" << std::endl;
313 double infNorm(
const DestinationType& U,
const DestinationType& Uerr )
const
315 typedef typename DestinationType :: ConstDofIteratorType ConstDofIteratorType ;
316 const ConstDofIteratorType uend = U.dend();
318 for( ConstDofIteratorType u = U.dbegin(), uerr = Uerr.dbegin(); u != uend; ++u, ++uerr )
321 double uerrval = *uerr ;
322 double div = std::abs( std::max( uval, uerrval ) );
324 double norm = std::abs( uval - uerrval );
325 if( std::abs(div) > 1e-12 )
327 res = std::max( res, norm );
332 HelmholtzOperatorType& helmholtzOp_;
333 LinearInverseOperatorType linearSolver_;
335 TimeStepControl timeStepControl_;
336 SourceTerm sourceTerm_;
345 DestinationType rhs_,temp_;
346 std::vector< DestinationType * > update_;
348 const int maxLinearIterations_;
350 const PreconditionOperatorType *preconditioner_;
ROW RungeKutta ODE solver.
Definition: basicrow.hh:56
void solve(DestinationType &U, MonitorType &monitor)
solve the system
Definition: basicrow.hh:216
BasicROWRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, TimeProviderType &timeProvider, const ButcherTable &butcherTable, const TimeStepControlType &timeStepControl, const SourceTermType &sourceTerm, const NonlinearSolverParameterType ¶meter)
constructor
Definition: basicrow.hh:84
void description(std::ostream &out) const
print description of ODE solver to out stream
Definition: basicrow.hh:306
int stages() const
return stages of RK solver (-1 if not implemented)
Definition: basicrow.hh:302
int order() const
return order of RK solver (-1 if not implemented)
Definition: basicrow.hh:304
void initialize(const DestinationType &U0)
apply operator once to get dt estimate
Definition: basicrow.hh:198
BasicROWRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, TimeProviderType &timeProvider, const ButcherTable &butcherTable, const NonlinearSolverParameterType ¶meter)
constructor
Definition: basicrow.hh:150
BasicROWRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, TimeProviderType &timeProvider, const ButcherTable &butcherTable, const TimeStepControlType &timeStepControl, const NonlinearSolverParameterType ¶meter)
constructor
Definition: basicrow.hh:127
~BasicROWRungeKuttaSolver()
destructor
Definition: basicrow.hh:191
Interface class for ODE Solver.
Definition: odesolverinterface.hh:21
Monitor MonitorType
monitor type
Definition: odesolverinterface.hh:59
virtual void solve(DestinationType &u)
solve where is the internal operator.
Definition: odesolverinterface.hh:75
DestinationImp DestinationType
type of destination
Definition: odesolverinterface.hh:62
void invert(bool doPivoting=true)
Compute inverse.
MAT & rightmultiply(const DenseMatrix< M2 > &M)
Multiplies M from the right to this matrix.
Definition: densematrix.hh:650
constexpr void mtv(const X &x, Y &y) const
y = A^T x
Definition: densematrix.hh:392
static bool verbose()
obtain the cached value for fem.verbose with default verbosity level 2
Definition: parameter.hh:466
general base for time providers
Definition: timeprovider.hh:36
This file implements a dense matrix with dynamic numbers of rows and columns.
This file implements a dense vector with a dynamic size.
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:489