dune-fem  2.4.1-rc
basicrow.hh
Go to the documentation of this file.
1 #ifndef DUNE_FEM_SOLVER_RUNGEKUTTA_BASICROW_HH
2 #define DUNE_FEM_SOLVER_RUNGEKUTTA_BASICROW_HH
3 
4 //- system includes
5 #include <cassert>
6 #include <limits>
7 #include <sstream>
8 #include <vector>
9 
10 //- dune-common includes
11 #include <dune/common/dynmatrix.hh>
12 #include <dune/common/dynvector.hh>
13 #include <dune/common/nullptr.hh>
14 
15 //- dune-fem includes
18 
19 namespace DuneODE
20 {
21 
22  // NoROWRungeKuttaSourceTerm
23  // ------------------------------
24 
26  {
27  template< class T >
28  bool operator() ( double time, double timeStepSize, int stage, const T &u, const std::vector< T * > &update, T &source )
29  {
30  return false;
31  }
32 
33  template< class T >
34  void limit( T& update, const double time ) {}
35 
36  template< class T >
37  double initialTimeStepEstimate ( double time, const T &u ) const
38  {
39  // return negative value to indicate that implicit time step should be used
40  return -1.0;
41  }
42 
43  double timeStepEstimate () const
44  {
46  }
47 
48  };
49 
51  : public Dune::Fem::LocalParameter< ROWSolverParameter, ROWSolverParameter >
52  {
53  protected:
54  // key prefix, default is fem.solver.row. (can be overloaded by user)
55  const std::string keyPrefix_;
56 
58 
59  public:
60  ROWSolverParameter( const std::string keyPrefix, const Dune::Fem::ParameterReader &parameter = Dune::Fem::Parameter::container() )
61  : keyPrefix_( keyPrefix ),
62  parameter_( parameter )
63  {}
64 
66  : keyPrefix_( "fem.solver.row." ),
67  parameter_( parameter )
68  {}
69 
70  const Dune::Fem::ParameterReader &parameter () const { return parameter_; }
71 
72  virtual double linAbsTolParameter ( ) const
73  {
74  return parameter().getValue< double >( keyPrefix_ + "linabstol", 1e-6 );
75  }
76 
77  virtual double linReductionParameter ( ) const
78  {
79  return parameter().getValue< double >( keyPrefix_ + "linreduction", 1e-4 );
80  }
81 
82  virtual bool linearSolverVerbose () const
83  {
84  return parameter().getValue< bool >( keyPrefix_ + "verbose", false );
85  }
86 
87  virtual int maxLinearIterationsParameter () const
88  {
89  return parameter().getValue< int >( keyPrefix_ + "maxlineariterations", std::numeric_limits< int >::max() );
90  }
91  };
92 
93 
94 
96  template< class HelmholtzOperator, class NonlinearSolver, class TimeStepControl, class SourceTerm = NoROWRungeKuttaSourceTerm >
98  : public OdeSolverInterface< typename HelmholtzOperator::DomainFunctionType >
99  {
102 
103  public:
106 
107  typedef HelmholtzOperator HelmholtzOperatorType;
108  typedef NonlinearSolver NonlinearSolverType;
109  typedef TimeStepControl TimeStepControlType;
110  typedef SourceTerm SourceTermType;
111 
113  typedef typename NonlinearSolver::ParametersType NonlinearSolverParametersType;
114 
115  typedef typename HelmholtzOperator::SpaceOperatorType::PreconditionOperatorType PreconditionOperatorType;
116 
118 
126  template< class ButcherTable >
127  BasicROWRungeKuttaSolver ( HelmholtzOperatorType &helmholtzOp,
128  TimeProviderType& timeProvider,
129  const ButcherTable &butcherTable,
130  const TimeStepControlType &timeStepControl,
131  const SourceTermType &sourceTerm,
132  const ParametersType& parameter,
133  const NonlinearSolverParametersType& nlsParam )
134  : helmholtzOp_( helmholtzOp ),
135  nonlinearSolver_( helmholtzOp_, nlsParam ),
136  timeStepControl_( timeStepControl ),
137  sourceTerm_( sourceTerm ),
138  stages_( butcherTable.stages() ),
139  alpha_( butcherTable.A() ),
140  alpha2_( butcherTable.B() ),
141  gamma_( stages() ),
142  beta_( stages() ),
143  c_( butcherTable.c() ),
144  rhs_( "RK rhs", helmholtzOp_.space() ),
145  temp_( "RK temp", helmholtzOp_.space() ),
146  update_( stages(), nullptr ),
147  linAbsTol_( parameter.linAbsTolParameter( ) ),
148  linReduction_( parameter.linReductionParameter( ) ),
149  linVerbose_( parameter.linearSolverVerbose() ),
150  maxLinearIterations_( parameter.maxLinearIterationsParameter() ),
151  preconditioner_(helmholtzOp.spaceOperator().preconditioner())
152  {
153  setup( butcherTable );
154  }
155 
156  template< class ButcherTable >
157  BasicROWRungeKuttaSolver ( HelmholtzOperatorType &helmholtzOp,
158  TimeProviderType& timeProvider,
159  const ButcherTable &butcherTable,
160  const TimeStepControlType &timeStepControl,
161  const SourceTermType &sourceTerm,
163  : BasicROWRungeKuttaSolver( helmholtzOp, timeProvider, butcherTable, timeStepControl, sourceTerm, ParametersType( parameter ),
164  NonlinearSolverParametersType( parameter ) )
165  {}
166 
173  template< class ButcherTable >
174  BasicROWRungeKuttaSolver ( HelmholtzOperatorType &helmholtzOp,
175  TimeProviderType& timeProvider,
176  const ButcherTable &butcherTable,
177  const TimeStepControlType &timeStepControl,
178  const ParametersType& parameter,
179  const NonlinearSolverParametersType& nlsParam )
180  : BasicROWRungeKuttaSolver( helmholtzOp, timeProvider, butcherTable, timeStepControl, SourceTermType(), parameter, nlsParam )
181  {}
182 
183  template< class ButcherTable >
184  BasicROWRungeKuttaSolver ( HelmholtzOperatorType &helmholtzOp,
185  TimeProviderType& timeProvider,
186  const ButcherTable &butcherTable,
187  const TimeStepControlType &timeStepControl,
189  : BasicROWRungeKuttaSolver( helmholtzOp, timeProvider, butcherTable, timeStepControl, SourceTermType(), ParametersType( parameter ),
190  NonlinearSolverParametersType( parameter ) )
191  {}
192 
198  template< class ButcherTable >
199  BasicROWRungeKuttaSolver ( HelmholtzOperatorType &helmholtzOp,
200  TimeProviderType& timeProvider,
201  const ButcherTable &butcherTable,
202  const ParametersType& parameter,
203  const NonlinearSolverParametersType& nlsParam )
204  : BasicROWRungeKuttaSolver( helmholtzOp, timeProvider, butcherTable, TimeStepControlType(), SourceTermType(), parameter, nlsParam )
205  {}
206 
207  template< class ButcherTable >
208  BasicROWRungeKuttaSolver ( HelmholtzOperatorType &helmholtzOp,
209  TimeProviderType& timeProvider,
210  const ButcherTable &butcherTable,
212  : BasicROWRungeKuttaSolver( helmholtzOp, timeProvider, butcherTable, TimeStepControlType(), SourceTermType(), ParametersType( parameter ),
213  NonlinearSolverParametersType( parameter ) )
214  {}
215 
216  template< class ButcherTable >
217  void setup( const ButcherTable& butcherTable )
218  {
219  std::cout << "ROW method of order=" << butcherTable.order() << " with " << stages_ << " stages" << std::endl;
220  // create intermediate functions
221  for( int i = 0; i < stages(); ++i )
222  {
223  std::ostringstream name;
224  name << "RK stage " << i;
225  update_[ i ] = new DestinationType( name.str(), helmholtzOp_.space() );
226  }
227 
228  // compute coefficients
229  for( int i = 0; i < stages(); ++i )
230  gamma_[ i ] = alpha_[ i ][ i ];
231 
232  alpha_.invert();
233  alpha_.mtv( butcherTable.b(), beta_ );
234  alpha2_.rightmultiply( alpha_ );
235  }
236 
239  {
240  for( int i = 0; i < stages(); ++i )
241  delete update_[ i ];
242  }
243 
245  void initialize ( const DestinationType &U0 )
246  {
247  const double time = timeStepControl_.time();
248 
249  helmholtzOp_.setTime( time );
250  helmholtzOp_.initializeTimeStepSize( U0 );
251  const double helmholtzEstimate = helmholtzOp_.timeStepEstimate();
252 
253  double sourceTermEstimate = sourceTerm_.initialTimeStepEstimate( time, U0 );
254  // negative time step is given by the empty source term
255  if( sourceTermEstimate < 0.0 ) sourceTermEstimate = helmholtzEstimate ;
256 
257  timeStepControl_.initialTimeStepSize( helmholtzEstimate, sourceTermEstimate );
258  }
259 
260  using BaseType::solve;
261 
263  void solve ( DestinationType &U, MonitorType &monitor )
264  {
265  monitor.reset();
266 
267  typename HelmholtzOperatorType::JacobianOperatorType jOp( "jacobianOperator", U.space(), U.space() );
268 
269  const double time = timeStepControl_.time();
270  const double timeStepSize = timeStepControl_.timeStepSize();
271  assert( timeStepSize > 0.0 );
272  // std::cout << "solving... " << time << " : " << U << std::endl;
273  for( int s = 0; s < stages(); ++s )
274  {
275  // update for stage s
276  DestinationType& updateStage = *update_[ s ];
277 
278  rhs_.assign( U );
279  for( int k = 0; k < s; ++k )
280  rhs_.axpy( alpha2_[ s ][ k ], *update_[ k ] );
281  helmholtzOp_.spaceOperator()(rhs_,updateStage);
282  updateStage *= (gamma_[s]*timeStepSize);
283  for( int k = 0; k < s; ++k )
284  updateStage.axpy( -gamma_[s]*alpha_[ s ][ k ], *update_[ k ] );
285 
286  rhs_.assign( updateStage );
287 
288  // solve the system
289  const double stageTime = time + c_[ s ]*timeStepSize;
290  helmholtzOp_.setTime( stageTime );
291  // compute jacobian if the diagonal entry in the butcher tableau changes
292  // if ( s==0 || (gamma_[s-1] != gamma_[s]) )
293  {
294  // std::cout << "lambda=" << gamma_[ s ]*timeStepSize << std::endl;
295  helmholtzOp_.setLambda( gamma_[ s ]*timeStepSize );
296  helmholtzOp_.jacobian( U, jOp );
297  }
298  const int remLinearIts = maxLinearIterations_;
299  if (preconditioner_)
300  {
301  typename NonlinearSolverType::LinearInverseOperatorType jInv( jOp, *preconditioner_, linReduction_, linAbsTol_, remLinearIts, linVerbose_ );
302  jInv( rhs_, updateStage );
303  monitor.linearSolverIterations_ += jInv.iterations();
304  }
305  else
306  {
307  typename NonlinearSolverType::LinearInverseOperatorType jInv( jOp, linReduction_, linAbsTol_, remLinearIts, linVerbose_ );
308  jInv( rhs_, updateStage );
309  monitor.linearSolverIterations_ += jInv.iterations();
310  }
311  }
312 
313  double error = 0.0;
314  if(0 && timeStepControl_.computeError() )
315  {
316  // store U (to be revised)
317  DestinationType Uerr( U );
318 
319  // update solution
320  U *= delta_;
321  for( int s = 0; s < stages(); ++s )
322  U.axpy( beta_[ s ], *update_[ s ] );
323 
324  //error = infNorm( U, Uerr );
325  Uerr.axpy( -1.0, U );
326  const double errorU = Uerr.scalarProductDofs( Uerr );
327  const double normU = U.scalarProductDofs( U );
328 
329  if( normU > 0 && errorU > 0 )
330  {
331  error = std::sqrt( errorU / normU );
332  }
333  std::cout << std::scientific << "Error in RK = " << error << " norm " << errorU << " " << normU << std::endl;
334  //std::cout << std::scientific << "Error in RK = " << error << std::endl;
335  }
336  else
337  {
338  // update solution
339  for( int s = 0; s < stages(); ++s )
340  U.axpy( beta_[ s ], *update_[ s ] );
341  }
342  // set error to monitor
343  monitor.error_ = error;
344 
345  // update time step size
346  timeStepControl_.timeStepEstimate( helmholtzOp_.timeStepEstimate(), sourceTerm_.timeStepEstimate(), monitor );
347  }
348 
349  int stages () const { return stages_; }
350 
351  void description ( std::ostream &out ) const
352  {
353  out << "Generic " << stages() << "-stage implicit Runge-Kutta solver.\\\\" << std::endl;
354  }
355 
356  protected:
357 
358  double infNorm(const DestinationType& U, const DestinationType& Uerr ) const
359  {
360  typedef typename DestinationType :: ConstDofIteratorType ConstDofIteratorType ;
361  const ConstDofIteratorType uend = U.dend();
362  double res = 0;
363  for( ConstDofIteratorType u = U.dbegin(), uerr = Uerr.dbegin(); u != uend; ++u, ++uerr )
364  {
365  double uval = *u;
366  double uerrval = *uerr ;
367  double div = std::abs( std::max( uval, uerrval ) );
368 
369  double norm = std::abs( uval - uerrval );
370  if( std::abs(div) > 1e-12 )
371  norm /= div;
372  res = std::max( res, norm );
373  }
374  return res;
375  }
376 
377  HelmholtzOperatorType &helmholtzOp_;
378  NonlinearSolverType nonlinearSolver_;
379  TimeStepControl timeStepControl_;
380  SourceTerm sourceTerm_;
381 
382  int stages_;
383  double delta_;
384  Dune::DynamicMatrix< double > alpha_, alpha2_;
385  Dune::DynamicVector< double > gamma_, beta_, c_;
386 
387  DestinationType rhs_,temp_;
388  std::vector< DestinationType * > update_;
389 
390  const double linAbsTol_, linReduction_;
391  const bool linVerbose_;
393 
394  const PreconditionOperatorType *preconditioner_;
395  };
396 
397 } // namespace DuneODE
398 
399 #endif // #ifndef DUNE_FEM_SOLVER_RUNGEKUTTA_BASICROW_HH
const Dune::Fem::ParameterReader & parameter() const
Definition: basicrow.hh:70
Dune::DynamicVector< double > gamma_
Definition: basicrow.hh:385
static double sqrt(const Double &v)
Definition: double.hh:870
SourceTerm SourceTermType
Definition: basicrow.hh:110
Interface class for ODE Solver.
Definition: odesolverinterface.hh:11
ROWSolverParameter(const std::string keyPrefix, const Dune::Fem::ParameterReader &parameter=Dune::Fem::Parameter::container())
Definition: basicrow.hh:60
BasicROWRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, TimeProviderType &timeProvider, const ButcherTable &butcherTable, const TimeStepControlType &timeStepControl, const SourceTermType &sourceTerm, const Dune::Fem::ParameterReader &parameter=Dune::Fem::Parameter::container())
Definition: basicrow.hh:157
DestinationType temp_
Definition: basicrow.hh:387
bool operator()(double time, double timeStepSize, int stage, const T &u, const std::vector< T * > &update, T &source)
Definition: basicrow.hh:28
virtual int maxLinearIterationsParameter() const
Definition: basicrow.hh:87
HelmholtzOperator HelmholtzOperatorType
Definition: basicrow.hh:107
static constexpr T max(T a)
Definition: utility.hh:65
BaseType::MonitorType MonitorType
Definition: basicrow.hh:104
void solve(DestinationType &U, MonitorType &monitor)
solve the system
Definition: basicrow.hh:263
Definition: basicrow.hh:50
ROWSolverParameter(const Dune::Fem::ParameterReader &parameter=Dune::Fem::Parameter::container())
Definition: basicrow.hh:65
std::vector< DestinationType * > update_
Definition: basicrow.hh:388
Definition: multistep.hh:16
virtual bool linearSolverVerbose() const
Definition: basicrow.hh:82
virtual double linAbsTolParameter() const
Definition: basicrow.hh:72
double timeStepEstimate() const
Definition: basicrow.hh:43
BasicROWRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, TimeProviderType &timeProvider, const ButcherTable &butcherTable, const TimeStepControlType &timeStepControl, const ParametersType &parameter, const NonlinearSolverParametersType &nlsParam)
constructor
Definition: basicrow.hh:174
const bool linVerbose_
Definition: basicrow.hh:391
TimeStepControl timeStepControl_
Definition: basicrow.hh:379
void limit(T &update, const double time)
Definition: basicrow.hh:34
Dune::DynamicMatrix< double > alpha_
Definition: basicrow.hh:384
ROWSolverParameter ParametersType
Definition: basicrow.hh:112
const std::string keyPrefix_
Definition: basicrow.hh:55
NonlinearSolverType nonlinearSolver_
Definition: basicrow.hh:378
general base for time providers
Definition: timeprovider.hh:35
Double abs(const Double &a)
Definition: double.hh:860
NonlinearSolver NonlinearSolverType
Definition: basicrow.hh:108
void initialize(const DestinationType &U0)
apply operator once to get dt estimate
Definition: basicrow.hh:245
HelmholtzOperator::SpaceOperatorType::PreconditionOperatorType PreconditionOperatorType
Definition: basicrow.hh:115
ROW RungeKutta ODE solver.
Definition: basicrow.hh:97
Definition: io/parameter.hh:549
const PreconditionOperatorType * preconditioner_
Definition: basicrow.hh:394
Dune::Fem::TimeProviderBase TimeProviderType
Definition: basicrow.hh:117
static ParameterContainer & container()
Definition: io/parameter.hh:190
virtual double linReductionParameter() const
Definition: basicrow.hh:77
int stages_
Definition: basicrow.hh:382
BaseType::DestinationType DestinationType
Definition: basicrow.hh:105
BasicROWRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, TimeProviderType &timeProvider, const ButcherTable &butcherTable, const ParametersType &parameter, const NonlinearSolverParametersType &nlsParam)
constructor
Definition: basicrow.hh:199
BasicROWRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, TimeProviderType &timeProvider, const ButcherTable &butcherTable, const TimeStepControlType &timeStepControl, const SourceTermType &sourceTerm, const ParametersType &parameter, const NonlinearSolverParametersType &nlsParam)
constructor
Definition: basicrow.hh:127
int stages() const
Definition: basicrow.hh:349
Dune::Fem::ParameterReader parameter_
Definition: basicrow.hh:57
Definition: odesolverinterface.hh:17
const double linReduction_
Definition: basicrow.hh:390
void setup(const ButcherTable &butcherTable)
Definition: basicrow.hh:217
NonlinearSolver::ParametersType NonlinearSolverParametersType
Definition: basicrow.hh:113
double initialTimeStepEstimate(double time, const T &u) const
Definition: basicrow.hh:37
SourceTerm sourceTerm_
Definition: basicrow.hh:380
double delta_
Definition: basicrow.hh:383
~BasicROWRungeKuttaSolver()
destructor
Definition: basicrow.hh:238
HelmholtzOperatorType & helmholtzOp_
Definition: basicrow.hh:377
BasicROWRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, TimeProviderType &timeProvider, const ButcherTable &butcherTable, const TimeStepControlType &timeStepControl, const Dune::Fem::ParameterReader &parameter=Dune::Fem::Parameter::container())
Definition: basicrow.hh:184
Definition: basicrow.hh:25
BasicROWRungeKuttaSolver(HelmholtzOperatorType &helmholtzOp, TimeProviderType &timeProvider, const ButcherTable &butcherTable, const Dune::Fem::ParameterReader &parameter=Dune::Fem::Parameter::container())
Definition: basicrow.hh:208
const int maxLinearIterations_
Definition: basicrow.hh:392
double infNorm(const DestinationType &U, const DestinationType &Uerr) const
Definition: basicrow.hh:358
void description(std::ostream &out) const
print description of ODE solver to out stream
Definition: basicrow.hh:351
TimeStepControl TimeStepControlType
Definition: basicrow.hh:109
DestinationImp DestinationType
type of destination
Definition: odesolverinterface.hh:53