DUNE-FEM (unstable)

basicrow.hh
1#ifndef DUNE_FEM_SOLVER_RUNGEKUTTA_BASICROW_HH
2#define DUNE_FEM_SOLVER_RUNGEKUTTA_BASICROW_HH
3
4//- system includes
5#include <algorithm>
6#include <cmath>
7#include <cassert>
8#include <limits>
9#include <sstream>
10#include <string>
11#include <vector>
12
13//- dune-common includes
16
17//- dune-fem includes
18#include <dune/fem/solver/odesolver.hh>
19#include <dune/fem/solver/parameter.hh>
20
21namespace DuneODE
22{
23
24 // NoROWRungeKuttaSourceTerm
25 // ------------------------------
26
27 struct NoROWRungeKuttaSourceTerm
28 {
29 template< class T >
30 bool operator() ( double time, double timeStepSize, int stage, const T &u, const std::vector< T * > &update, T &source )
31 {
32 return false;
33 }
34
35 template< class T >
36 void limit( T& update, const double time ) {}
37
38 template< class T >
39 double initialTimeStepEstimate ( double time, const T &u ) const
40 {
41 // return negative value to indicate that implicit time step should be used
42 return -1.0;
43 }
44
45 double timeStepEstimate () const
46 {
48 }
49
50 };
51
53 template< class HelmholtzOperator, class NonlinearSolver, class TimeStepControl, class SourceTerm = NoROWRungeKuttaSourceTerm >
55 : public OdeSolverInterface< typename HelmholtzOperator::DomainFunctionType >
56 {
59
60 public:
61 typedef typename BaseType::MonitorType MonitorType;
62 typedef typename BaseType::DestinationType DestinationType;
63
64 typedef HelmholtzOperator HelmholtzOperatorType;
65 typedef NonlinearSolver NonlinearSolverType;
66 typedef TimeStepControl TimeStepControlType;
67 typedef SourceTerm SourceTermType;
68
69 typedef typename NonlinearSolver::ParameterType NonlinearSolverParameterType;
70 typedef typename NonlinearSolverType::LinearInverseOperatorType LinearInverseOperatorType;
71
72 typedef typename HelmholtzOperator::SpaceOperatorType::PreconditionOperatorType PreconditionOperatorType;
73
75
83 template< class ButcherTable >
84 BasicROWRungeKuttaSolver ( HelmholtzOperatorType &helmholtzOp,
85 TimeProviderType& timeProvider,
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() ),
98 gamma_( stages() ),
99 beta_( stages() ),
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())
106 {
107 setup( butcherTable );
108 }
109
110 template< class ButcherTable >
111 BasicROWRungeKuttaSolver ( HelmholtzOperatorType &helmholtzOp,
112 TimeProviderType& timeProvider,
113 const ButcherTable &butcherTable,
114 const TimeStepControlType &timeStepControl,
115 const SourceTermType &sourceTerm,
116 const Dune::Fem::ParameterReader &parameter = Dune::Fem::Parameter::container() )
117 : BasicROWRungeKuttaSolver( helmholtzOp, timeProvider, butcherTable, timeStepControl, sourceTerm, NonlinearSolverParameterType( parameter ) )
118 {}
119
126 template< class ButcherTable >
127 BasicROWRungeKuttaSolver ( HelmholtzOperatorType &helmholtzOp,
128 TimeProviderType& timeProvider,
129 const ButcherTable &butcherTable,
130 const TimeStepControlType &timeStepControl,
131 const NonlinearSolverParameterType& parameter )
132 : BasicROWRungeKuttaSolver( helmholtzOp, timeProvider, butcherTable, timeStepControl, SourceTermType(), parameter )
133 {}
134
135 template< class ButcherTable >
136 BasicROWRungeKuttaSolver ( HelmholtzOperatorType &helmholtzOp,
137 TimeProviderType& timeProvider,
138 const ButcherTable &butcherTable,
139 const TimeStepControlType &timeStepControl,
140 const Dune::Fem::ParameterReader &parameter = Dune::Fem::Parameter::container() )
141 : BasicROWRungeKuttaSolver( helmholtzOp, timeProvider, butcherTable, timeStepControl, SourceTermType(), NonlinearSolverParameterType( parameter ) )
142 {}
143
149 template< class ButcherTable >
150 BasicROWRungeKuttaSolver ( HelmholtzOperatorType &helmholtzOp,
151 TimeProviderType& timeProvider,
152 const ButcherTable &butcherTable,
153 const NonlinearSolverParameterType& parameter )
154 : BasicROWRungeKuttaSolver( helmholtzOp, timeProvider, butcherTable, TimeStepControlType(), SourceTermType(), parameter )
155 {}
156
157 template< class ButcherTable >
158 BasicROWRungeKuttaSolver ( HelmholtzOperatorType &helmholtzOp,
159 TimeProviderType& timeProvider,
160 const ButcherTable &butcherTable,
161 const Dune::Fem::ParameterReader &parameter = Dune::Fem::Parameter::container() )
162 : BasicROWRungeKuttaSolver( helmholtzOp, timeProvider, butcherTable, TimeStepControlType(), SourceTermType(), NonlinearSolverParameterType( parameter ) )
163 {}
164
165 template< class ButcherTable >
166 void setup( const ButcherTable& butcherTable )
167 {
169 {
170 std::cout << "ROW method of order=" << butcherTable.order() << " with " << stages_ << " stages" << std::endl;
171 }
172
173 // create intermediate functions
174 for( int i = 0; i < stages(); ++i )
175 {
176 std::ostringstream name;
177 name << "RK stage " << i;
178 update_[ i ] = new DestinationType( name.str(), helmholtzOp_.space() );
179 }
180
181 // compute coefficients
182 for( int i = 0; i < stages(); ++i )
183 gamma_[ i ] = alpha_[ i ][ i ];
184
185 alpha_.invert();
186 alpha_.mtv( butcherTable.b(), beta_ );
187 alpha2_.rightmultiply( alpha_ );
188 }
189
192 {
193 for( int i = 0; i < stages(); ++i )
194 delete update_[ i ];
195 }
196
198 void initialize ( const DestinationType &U0 )
199 {
200 const double time = timeStepControl_.time();
201
202 helmholtzOp_.setTime( time );
203 helmholtzOp_.initializeTimeStepSize( U0 );
204 const double helmholtzEstimate = helmholtzOp_.timeStepEstimate();
205
206 double sourceTermEstimate = sourceTerm_.initialTimeStepEstimate( time, U0 );
207 // negative time step is given by the empty source term
208 if( sourceTermEstimate < 0.0 ) sourceTermEstimate = helmholtzEstimate ;
209
210 timeStepControl_.initialTimeStepSize( helmholtzEstimate, sourceTermEstimate );
211 }
212
213 using BaseType::solve;
214
216 void solve ( DestinationType &U, MonitorType &monitor )
217 {
218 monitor.reset();
219
220 typename HelmholtzOperatorType::JacobianOperatorType jOp( "jacobianOperator", U.space(), U.space() );
221
222 const double time = timeStepControl_.time();
223 const double timeStepSize = timeStepControl_.timeStepSize();
224 assert( timeStepSize > 0.0 );
225 // std::cout << "solving... " << time << " : " << U << std::endl;
226 for( int s = 0; s < stages(); ++s )
227 {
228 // update for stage s
229 DestinationType& updateStage = *update_[ s ];
230
231 rhs_.assign( U );
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 ] );
238
239 rhs_.assign( updateStage );
240
241 // solve the system
242 const double stageTime = time + c_[ s ]*timeStepSize;
243 helmholtzOp_.setTime( stageTime );
244 // compute jacobian if the diagonal entry in the butcher tableau changes
245 // if ( s==0 || (gamma_[s-1] != gamma_[s]) )
246 {
247 // std::cout << "lambda=" << gamma_[ s ]*timeStepSize << std::endl;
248 helmholtzOp_.setLambda( gamma_[ s ]*timeStepSize );
249 helmholtzOp_.jacobian( U, jOp );
250 }
251 const int remLinearIts = maxLinearIterations_;
252
253 linearSolver_.parameter().setMaxIterations( remLinearIts );
254
255 if (preconditioner_)
256 linearSolver_.bind( jOp, *preconditioner_ );
257 else
258 linearSolver_.bind( jOp );
259
260 linearSolver_( rhs_, updateStage );
261 monitor.linearSolverIterations_ += linearSolver_.iterations();
262
263 linearSolver_.unbind();
264 }
265
266 double error = 0.0;
267 if(0 && timeStepControl_.computeError() )
268 {
269 // store U (to be revised)
270 DestinationType Uerr( U );
271
272 // update solution
273 U *= delta_;
274 for( int s = 0; s < stages(); ++s )
275 U.axpy( beta_[ s ], *update_[ s ] );
276
277 //error = infNorm( U, Uerr );
278 Uerr.axpy( -1.0, U );
279 const double errorU = Uerr.scalarProductDofs( Uerr );
280 const double normU = U.scalarProductDofs( U );
281
282 if( normU > 0 && errorU > 0 )
283 {
284 error = std::sqrt( errorU / normU );
285 }
286 std::cout << std::scientific << "Error in RK = " << error << " norm " << errorU << " " << normU << std::endl;
287 //std::cout << std::scientific << "Error in RK = " << error << std::endl;
288 }
289 else
290 {
291 // update solution
292 for( int s = 0; s < stages(); ++s )
293 U.axpy( beta_[ s ], *update_[ s ] );
294 }
295 // set error to monitor
296 monitor.error_ = error;
297
298 // update time step size
299 timeStepControl_.timeStepEstimate( helmholtzOp_.timeStepEstimate(), sourceTerm_.timeStepEstimate(), monitor );
300 }
301
302 int stages () const { return stages_; }
303
304 int order () const { return ord_; }
305
306 void description ( std::ostream &out ) const
307 {
308 out << "Generic " << stages() << "-stage implicit Runge-Kutta solver.\\\\" << std::endl;
309 }
310
311 protected:
312
313 double infNorm(const DestinationType& U, const DestinationType& Uerr ) const
314 {
315 typedef typename DestinationType :: ConstDofIteratorType ConstDofIteratorType ;
316 const ConstDofIteratorType uend = U.dend();
317 double res = 0;
318 for( ConstDofIteratorType u = U.dbegin(), uerr = Uerr.dbegin(); u != uend; ++u, ++uerr )
319 {
320 double uval = *u;
321 double uerrval = *uerr ;
322 double div = std::abs( std::max( uval, uerrval ) );
323
324 double norm = std::abs( uval - uerrval );
325 if( std::abs(div) > 1e-12 )
326 norm /= div;
327 res = std::max( res, norm );
328 }
329 return res;
330 }
331
332 HelmholtzOperatorType& helmholtzOp_;
333 LinearInverseOperatorType linearSolver_;
334
335 TimeStepControl timeStepControl_;
336 SourceTerm sourceTerm_;
337
338 int stages_;
339 int ord_;
340
341 double delta_;
342 Dune::DynamicMatrix< double > alpha_, alpha2_;
343 Dune::DynamicVector< double > gamma_, beta_, c_;
344
345 DestinationType rhs_,temp_;
346 std::vector< DestinationType * > update_;
347
348 const int maxLinearIterations_;
349
350 const PreconditionOperatorType *preconditioner_;
351 };
352
353} // namespace DuneODE
354
355#endif // #ifndef DUNE_FEM_SOLVER_RUNGEKUTTA_BASICROW_HH
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 &parameter)
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 &parameter)
constructor
Definition: basicrow.hh:150
BasicROWRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, TimeProviderType &timeProvider, const ButcherTable &butcherTable, const TimeStepControlType &timeStepControl, const NonlinearSolverParameterType &parameter)
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
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Jan 10, 23:35, 2026)