explicit_timesteppers.h
Go to the documentation of this file.
1 // LIC// ====================================================================
2 // LIC// This file forms part of oomph-lib, the object-oriented,
3 // LIC// multi-physics finite-element library, available
4 // LIC// at http://www.oomph-lib.org.
5 // LIC//
6 // LIC// Copyright (C) 2006-2023 Matthias Heil and Andrew Hazel
7 // LIC//
8 // LIC// This library is free software; you can redistribute it and/or
9 // LIC// modify it under the terms of the GNU Lesser General Public
10 // LIC// License as published by the Free Software Foundation; either
11 // LIC// version 2.1 of the License, or (at your option) any later version.
12 // LIC//
13 // LIC// This library is distributed in the hope that it will be useful,
14 // LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // LIC// Lesser General Public License for more details.
17 // LIC//
18 // LIC// You should have received a copy of the GNU Lesser General Public
19 // LIC// License along with this library; if not, write to the Free Software
20 // LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21 // LIC// 02110-1301 USA.
22 // LIC//
23 // LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24 // LIC//
25 // LIC//====================================================================
26 // Header functions for classes that represent explicit time-stepping schemes
27 
28 // Include guards to prevent multiple inclusion of this header
29 #ifndef OOMPH_EXPLICIT_TIMESTEPPERS
30 #define OOMPH_EXPLICIT_TIMESTEPPERS
31 
32 // Config header generated by autoconfig
33 #ifdef HAVE_CONFIG_H
34 #include <oomph-lib-config.h>
35 #endif
36 
37 
38 #ifdef OOMPH_HAS_MPI
39 #include "mpi.h"
40 #endif
41 
42 #include "Vector.h"
43 #include "double_vector.h"
44 #include "oomph_utilities.h"
45 
46 namespace oomph
47 {
48  class Time;
49 
50  //===============================================================
51  /// Class for objects than can be advanced in time by an Explicit
52  /// Timestepper.
53  /// WARNING: For explicit time stepping to work the object's residual
54  /// function (as used by get_inverse_mass_matrix_times_residuals(..)) MUST
55  /// be in the form r = f(t, u) - [timestepper approximation to dudt]!
56  /// Standard implicit time stepping will work with plenty of residuals that
57  /// don't fit into this form. Some examples where implicit time stepping
58  /// will work fine but explicit will fail:
59  /// 1) The negation of the above formula, this implementation will end up
60  /// using dudt = - f(u,t).
61  /// 2) A residual which is implicit or non-linear in dudt, such as r = dudt
62  /// - u x dudt.
63  //===============================================================
65  {
66  // Dummy double value for time
67  static double Dummy_time_value;
68 
69  public:
70  /// Empty constructor
72 
73  /// Broken copy constructor
75 
76  /// Broken assignment operator
77  void operator=(const ExplicitTimeSteppableObject&) = delete;
78 
79  /// Empty destructor
81 
82  /// A single virtual function that returns the residuals
83  /// vector multiplied by the inverse mass matrix
84  virtual void get_dvaluesdt(DoubleVector& minv_res);
85 
86  /// Function that gets the values of the dofs in the object
87  virtual void get_dofs(DoubleVector& dofs) const;
88 
89  /// Function that gets the history values of the dofs in the object
90  virtual void get_dofs(const unsigned& t, DoubleVector& dofs) const;
91 
92  /// Function that sets the values of the dofs in the object
93  virtual void set_dofs(const DoubleVector& dofs);
94 
95  /// Function that adds the values to the dofs
96  virtual void add_to_dofs(const double& lambda,
97  const DoubleVector& increment_dofs);
98 
99  /// Empty virtual function to do anything needed before a stage of
100  /// an explicit time step (Runge-Kutta steps contain multiple stages per
101  /// time step, most others only contain one).
103 
104  /// Empty virtual function that should be overloaded to update any
105  /// dependent data or boundary conditions that should be advanced after each
106  /// stage of an explicit time step (Runge-Kutta steps contain multiple
107  /// stages per time step, most others only contain one).
109 
110  /// Empty virtual function that can be overloaded to do anything
111  /// needed before an explicit step.
113 
114  /// Empty virtual function that can be overloaded to do anything
115  /// needed after an explicit step.
117 
118  /// Broken virtual function that should be overloaded to
119  /// return access to the local time in the object
120  virtual double& time();
121 
122  /// Virtual function that should be overloaded to return a pointer to
123  /// a Time object.
124  virtual Time* time_pt() const;
125  };
126 
127 
128  //=====================================================================
129  /// A Base class for explicit timesteppers
130  //=====================================================================
132  {
133  protected:
134  /// String that indicates the type of the timestepper
135  /// (e.g. "RungeKutta", etc.)
137 
138  public:
139  /// Empty Constructor.
141 
142  /// Broken copy constructor
144 
145  /// Broken assignment operator
146  void operator=(const ExplicitTimeStepper&) = delete;
147 
148  /// Empty virtual destructor --- no memory is allocated in this class
149  virtual ~ExplicitTimeStepper() {}
150 
151  /// Pure virtual function that is used to advance time in the object
152  // referenced by object_pt by an amount dt
153  virtual void timestep(ExplicitTimeSteppableObject* const& object_pt,
154  const double& dt) = 0;
155  };
156 
157 
158  /// ===========================================================
159  /// Simple first-order Euler Timestepping
160  //============================================================
161  class Euler : public ExplicitTimeStepper
162  {
163  public:
164  /// Constructor, set the type
166  {
167  Type = "Euler";
168  }
169 
170  /// Broken copy constructor
171  Euler(const Euler&) = delete;
172 
173  /// Broken assignment operator
174  void operator=(const Euler&) = delete;
175 
176  /// Overload function that is used to advance time in the object
177  /// reference by object_pt by an amount dt
178  void timestep(ExplicitTimeSteppableObject* const& object_pt,
179  const double& dt);
180  };
181 
182  /// ===========================================================
183  /// Standard Runge Kutta Timestepping
184  //============================================================
185  template<unsigned ORDER>
187  {
188  public:
189  /// Constructor, set the type
191  {
192  Type = "RungeKutta";
193  }
194 
195  /// Broken copy constructor
196  RungeKutta(const RungeKutta&) = delete;
197 
198  /// Broken assignment operator
199  void operator=(const RungeKutta&) = delete;
200 
201  /// Function that is used to advance time in the object
202  // reference by object_pt by an amount dt
203  void timestep(ExplicitTimeSteppableObject* const& object_pt,
204  const double& dt);
205  };
206 
207 
208  /// ===========================================================
209  /// Runge Kutta Timestepping that uses low storage
210  //============================================================
211  template<unsigned ORDER>
213  {
214  // Storage for the coefficients
216 
217  public:
218  /// Constructor, set the type
220 
221  /// Broken copy constructor
223 
224  /// Broken assignment operator
225  void operator=(const LowStorageRungeKutta&) = delete;
226 
227  /// Function that is used to advance the solution by time dt
228  void timestep(ExplicitTimeSteppableObject* const& object_pt,
229  const double& dt);
230  };
231 
232 
233  /// ===========================================================
234  /// An explicit version of BDF3 (i.e. uses derivative evaluation at y_n
235  /// instead of y_{n+1}). Useful as a predictor because it is third order
236  /// accurate but requires only one function evaluation (i.e. only one mass
237  /// matrix inversion + residual calculation).
238  //============================================================
239  class EBDF3 : public ExplicitTimeStepper
240  {
241  double Yn_weight;
242  double Ynm1_weight;
243  double Ynm2_weight;
244  double Fn_weight;
245 
246  public:
247  /// Constructor, set the type
248  EBDF3() {}
249 
250  /// Broken copy constructor
251  EBDF3(const EBDF3&) = delete;
252 
253  /// Broken assignment operator
254  void operator=(const EBDF3&) = delete;
255 
256  void set_weights(const double& dtn,
257  const double& dtnm1,
258  const double& dtnm2);
259 
260  /// Function that is used to advance the solution by time dt
261  void timestep(ExplicitTimeSteppableObject* const& object_pt,
262  const double& dt);
263  };
264 
265 } // namespace oomph
266 
267 #endif
char t
Definition: cfortran.h:568
A vector in the mathematical sense, initially developed for linear algebra type applications....
Definition: double_vector.h:58
=========================================================== An explicit version of BDF3 (i....
void timestep(ExplicitTimeSteppableObject *const &object_pt, const double &dt)
Function that is used to advance the solution by time dt.
EBDF3()
Constructor, set the type.
void operator=(const EBDF3 &)=delete
Broken assignment operator.
void set_weights(const double &dtn, const double &dtnm1, const double &dtnm2)
Calculate the weights for this set of step sizes.
EBDF3(const EBDF3 &)=delete
Broken copy constructor.
=========================================================== Simple first-order Euler Timestepping
void operator=(const Euler &)=delete
Broken assignment operator.
Euler()
Constructor, set the type.
Euler(const Euler &)=delete
Broken copy constructor.
void timestep(ExplicitTimeSteppableObject *const &object_pt, const double &dt)
Overload function that is used to advance time in the object reference by object_pt by an amount dt.
Class for objects than can be advanced in time by an Explicit Timestepper. WARNING: For explicit time...
virtual void add_to_dofs(const double &lambda, const DoubleVector &increment_dofs)
Function that adds the values to the dofs.
virtual ~ExplicitTimeSteppableObject()
Empty destructor.
virtual void get_dofs(DoubleVector &dofs) const
Function that gets the values of the dofs in the object.
virtual void actions_after_explicit_timestep()
Empty virtual function that can be overloaded to do anything needed after an explicit step.
virtual double & time()
Broken virtual function that should be overloaded to return access to the local time in the object.
virtual Time * time_pt() const
Virtual function that should be overloaded to return a pointer to a Time object.
virtual void actions_before_explicit_stage()
Empty virtual function to do anything needed before a stage of an explicit time step (Runge-Kutta ste...
virtual void set_dofs(const DoubleVector &dofs)
Function that sets the values of the dofs in the object.
void operator=(const ExplicitTimeSteppableObject &)=delete
Broken assignment operator.
virtual void get_dvaluesdt(DoubleVector &minv_res)
A single virtual function that returns the residuals vector multiplied by the inverse mass matrix.
virtual void actions_after_explicit_stage()
Empty virtual function that should be overloaded to update any dependent data or boundary conditions ...
virtual void actions_before_explicit_timestep()
Empty virtual function that can be overloaded to do anything needed before an explicit step.
static double Dummy_time_value
Dummy value of time always set to zero.
ExplicitTimeSteppableObject(const ExplicitTimeSteppableObject &)=delete
Broken copy constructor.
A Base class for explicit timesteppers.
ExplicitTimeStepper(const ExplicitTimeStepper &)=delete
Broken copy constructor.
ExplicitTimeStepper()
Empty Constructor.
virtual ~ExplicitTimeStepper()
Empty virtual destructor — no memory is allocated in this class.
std::string Type
String that indicates the type of the timestepper (e.g. "RungeKutta", etc.)
virtual void timestep(ExplicitTimeSteppableObject *const &object_pt, const double &dt)=0
Pure virtual function that is used to advance time in the object.
void operator=(const ExplicitTimeStepper &)=delete
Broken assignment operator.
=========================================================== Runge Kutta Timestepping that uses low st...
LowStorageRungeKutta(const LowStorageRungeKutta &)=delete
Broken copy constructor.
void operator=(const LowStorageRungeKutta &)=delete
Broken assignment operator.
LowStorageRungeKutta()
Constructor, set the type.
void timestep(ExplicitTimeSteppableObject *const &object_pt, const double &dt)
Function that is used to advance the solution by time dt.
=========================================================== Standard Runge Kutta Timestepping
RungeKutta(const RungeKutta &)=delete
Broken copy constructor.
void operator=(const RungeKutta &)=delete
Broken assignment operator.
RungeKutta()
Constructor, set the type.
void timestep(ExplicitTimeSteppableObject *const &object_pt, const double &dt)
Function that is used to advance time in the object.
Class to keep track of discrete/continous time. It is essential to have a single Time object when usi...
Definition: timesteppers.h:63
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...