explicit_timesteppers.cc
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 // Non-inline member functions for the explicit timesteppers
27 #include "explicit_timesteppers.h"
28 #include "timesteppers.h"
29 
30 namespace oomph
31 {
32  //===================================================================
33  /// Dummy value of time always set to zero
34  //==================================================================
36 
37  //==================================================================
38  /// A single virtual function that returns the residuals
39  /// vector multiplied by the inverse mass matrix
40  //=================================================================
42  {
43  std::ostringstream error_stream;
44  error_stream
45  << "Empty default function called.\n"
46  << "The function must return the solution x of the linear system\n"
47  << " M x = R\n"
48  << "in order for the object to be used by an ExplicitTimeStepper.\n"
49  << "NOTE: It is the responsibility of the object to set the size \n"
50  << " of the vector x\n";
51 
52 
53  throw OomphLibError(
54  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
55  }
56 
57  //=======================================================================
58  /// Function that should get the values of the dofs in the object
59  //=======================================================================
61  {
62  std::ostringstream error_stream;
63  error_stream
64  << "Empty default function called.\n"
65  << "The function must return the current values of the degrees of \n"
66  << "freedom in the object.\n"
67  << "Note: It is the responsibility of the object to set the size\n"
68  << "of the vector\n";
69 
70  throw OomphLibError(
71  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
72  }
73 
74  //=======================================================================
75  /// Function that should get the values of the dofs in the object
76  //=======================================================================
78  DoubleVector& dofs) const
79  {
80  std::ostringstream error_stream;
81  error_stream
82  << "Empty default function called.\n"
83  << "The function must return the t'th history values of the degrees of \n"
84  << "freedom in the object.\n"
85  << "Note: It is the responsibility of the object to set the size\n"
86  << "of the vector\n";
87 
88  throw OomphLibError(
89  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
90  }
91 
92  //=====================================================================
93  /// Function that sets the values of the dofs in the object
94  //====================================================================
96  {
97  std::ostringstream error_stream;
98  error_stream
99  << "Empty default function called.\n"
100  << "The function must set the current values of the degrees of \n"
101  << "freedom in the object.\n";
102 
103  throw OomphLibError(
104  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
105  }
106 
107 
108  //====================================================================
109  /// Function that adds the lambda multiplied by the increment_dofs
110  /// vector to the dofs in the object
111  //====================================================================
113  const double& lambda, const DoubleVector& increment_dofs)
114  {
115  std::ostringstream error_stream;
116  error_stream
117  << "Empty default function called.\n"
118  << "The function must add lambda multiplied by the contents of the\n"
119  << "input vector to the degrees of freedom in the object.\n"
120  << "Note: It is the responsibility of the object to ensure that the\n"
121  << " the degrees of freedom are in the same order as those \n"
122  << " returned by get_dvaluesdt()\n";
123 
124  throw OomphLibError(
125  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
126  }
127 
128  //==================================================================
129  /// Virtual function that should be overloaded to return access
130  /// to the local time in the object
131  //=================================================================
133  {
134  std::ostringstream error_stream;
135  error_stream << "Empty default function called.\n"
136  << "The function must return a reference to the local time in "
137  "the object\n";
138 
139  throw OomphLibError(
140  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
141 
142  return Dummy_time_value;
143  }
144 
145  /// Virtual function that should be overloaded to return a pointer to a
146  /// Time object.
148  {
149  std::ostringstream error_stream;
150  error_stream
151  << "Empty default function called.\n"
152  << "The function must return a pointer to an oomph-lib Time object.\n";
153  throw OomphLibError(
154  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
155 
156  return 0;
157  }
158 
159 
160  //================================================================
161  /// Euler timestepping x^{t+1} = x^{t} + dt M^{-1} L(x^{t})
162  //=================================================================
164  const double& dt)
165  {
167  object_pt->actions_before_explicit_stage();
168 
169  // Vector that will hold the inverse mass matrix multiplied by the
170  // residuals
171  DoubleVector minv_res;
172  // Get M^{-1} R
173  object_pt->get_dvaluesdt(minv_res);
174 
175  // Add the result to the unknowns
176  object_pt->add_to_dofs(dt, minv_res);
177  // Increment the time by the appropriate amount
178  object_pt->time() += dt;
179 
180  // Call any actions required after the change in the unknowns
181  object_pt->actions_after_explicit_stage();
182  object_pt->actions_after_explicit_timestep();
183  }
184 
185 
186  //====================================================================
187  // Broken default timestep function for RungeKutta schemes
188  //====================================================================
189  template<unsigned ORDER>
191  ExplicitTimeSteppableObject* const& object_pt, const double& dt)
192  {
193  std::ostringstream error_stream;
194  error_stream << "Timestep not implemented for order " << ORDER << "\n";
195  throw OomphLibError(
196  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
197  }
198 
199  //===================================================================
200  /// Explicit specialisation for fourth-order RK scheme
201  //==================================================================
202  template<>
204  const double& dt)
205  {
207 
208  // Stage 1
209  // ============================================================
210  object_pt->actions_before_explicit_stage();
211 
212  // Store the initial values and initial time
213  DoubleVector u;
214  object_pt->get_dofs(u);
215 
216  // Now get the first unknowns
217  DoubleVector k1;
218  object_pt->get_dvaluesdt(k1);
219 
220  // Add to the residuals
221  object_pt->add_to_dofs(0.5 * dt, k1);
222  // Increment the time
223  object_pt->time() += 0.5 * dt;
224  object_pt->actions_after_explicit_stage();
225 
226 
227  // Stage 2
228  // ============================================================
229  object_pt->actions_before_explicit_stage();
230 
231  // Get the next unknowns
232  DoubleVector k2;
233  object_pt->get_dvaluesdt(k2);
234 
235  // Now reset the residuals
236  object_pt->set_dofs(u);
237  object_pt->add_to_dofs(0.5 * dt, k2);
238  // Time remains the same
239  object_pt->actions_after_explicit_stage();
240 
241  // Stage 3
242  // ============================================================
243  object_pt->actions_before_explicit_stage();
244 
245  // Get the next unknowns
246  DoubleVector k3;
247  object_pt->get_dvaluesdt(k3);
248 
249  // Now reset the residuals
250  object_pt->set_dofs(u);
251  object_pt->add_to_dofs(dt, k3);
252  // Increment the time (now at initial_time + dt)
253  object_pt->time() += 0.5 * dt;
254  object_pt->actions_after_explicit_stage();
255 
256  // Stage 4
257  // ============================================================
258  object_pt->actions_before_explicit_stage();
259 
260  // Get the final unknowns
261  DoubleVector k4;
262  object_pt->get_dvaluesdt(k4);
263 
264  // Set the final values of the unknowns
265  object_pt->set_dofs(u);
266  object_pt->add_to_dofs((dt / 6.0), k1);
267  object_pt->add_to_dofs((dt / 3.0), k2);
268  object_pt->add_to_dofs((dt / 3.0), k3);
269  object_pt->add_to_dofs((dt / 6.0), k4);
270  object_pt->actions_after_explicit_stage();
271 
272  object_pt->actions_after_explicit_timestep();
273  }
274 
275  //===================================================================
276  /// Explicit specialisation for second-order RK scheme
277  //==================================================================
278  template<>
280  const double& dt)
281  {
283 
284  // Stage 1
285  // ============================================================
286  object_pt->actions_before_explicit_stage();
287 
288  // Store the initial values
289  DoubleVector u;
290  object_pt->get_dofs(u);
291 
292  // Get f1 (time derivative at t0, y0) and add to dofs
293  DoubleVector f1;
294  object_pt->get_dvaluesdt(f1);
295  object_pt->add_to_dofs(dt, f1);
296 
297  // Advance time to t1 = t0 + dt
298  object_pt->time() += dt;
299 
300  object_pt->actions_after_explicit_stage();
301 
302 
303  // Stage 2
304  // ============================================================
305  object_pt->actions_before_explicit_stage();
306 
307  // get f2 (with t=t1, y = y0 + h f0)
308  DoubleVector f2;
309  object_pt->get_dvaluesdt(f2);
310 
311  // Final answer is starting dofs + h/2 * (f1 + f2)
312  object_pt->set_dofs(u);
313  object_pt->add_to_dofs(0.5 * dt, f1);
314  object_pt->add_to_dofs(0.5 * dt, f2);
315 
316  object_pt->actions_after_explicit_stage();
317 
318  // Done, do actions
319  object_pt->actions_after_explicit_timestep();
320  }
321 
322 
323  //=================================================================
324  // General constructor for LowOrder RK schemes
325  //==================================================================
326  template<unsigned ORDER>
328  {
329  Type = "LowStorageRungeKutta";
330  }
331 
332  //======================================================================
333  /// Broken default timestep for LowStorageRungeKutta
334  //======================================================================
335  template<unsigned ORDER>
337  ExplicitTimeSteppableObject* const& object_pt, const double& dt)
338  {
339  std::ostringstream error_stream;
340  error_stream << "Timestep not implemented for order " << ORDER << "\n";
341  throw OomphLibError(
342  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
343  }
344 
345  //================================================================
346  // Specialised constructor for fourth-order RK scheme
347  //================================================================
348  template<>
350  {
351  Type = "LowStorageRungeKutta";
352 
353  A.resize(5);
354  A[0] = 0.0;
355  A[1] = -567301805773.0 / 1357537059087.0;
356  A[2] = -2404267990393.0 / 2016746695238.0;
357  A[3] = -3550918686646.0 / 2091501179385.0;
358  A[4] = -1275806237668.0 / 842570457699.0;
359 
360  B.resize(5);
361  B[0] = 1432997174477.0 / 9575080441755.0;
362  B[1] = 5161836677717.0 / 13612068292357.0;
363  B[2] = 1720146321549.0 / 2090206949498.0;
364  B[3] = 3134564353537.0 / 4481467310338.0;
365  B[4] = 2277821191437.0 / 14882151754819.0;
366 
367  C.resize(5);
368  C[0] = B[0];
369  C[1] = 2526269341429.0 / 6820363962896.0;
370  C[2] = 2006345519317.0 / 3224310063776.0;
371  C[3] = 2802321613138.0 / 2924317926251.0;
372  C[4] = 1.0;
373  }
374 
375 
376  // Explicit specialisation for fourth-order RK scheme
377  template<>
379  ExplicitTimeSteppableObject* const& object_pt, const double& dt)
380  {
382 
383  // Store the initial time
384  const double initial_time = object_pt->time();
385  // Temporary storage
386  DoubleVector k;
387 
388  // Temporary storage for the inverse mass matrix multiplied by the residuals
389  DoubleVector minv_res;
390 
391  // Loop over the number of steps in the scheme
392  for (unsigned i = 0; i < 5; i++)
393  {
394  object_pt->actions_before_explicit_stage();
395 
396  // Get the inverse mass matrix multiplied by the current value
397  // of the residuals
398  object_pt->get_dvaluesdt(minv_res);
399  // Get the values of k
400  const unsigned n_dof = minv_res.nrow();
401 
402  // First time round resize k and initialise to zero
403  if (i == 0)
404  {
405  k.build(minv_res.distribution_pt(), 0.0);
406  }
407  // Now construct the next value of k
408  for (unsigned n = 0; n < n_dof; n++)
409  {
410  k[n] *= A[i];
411  k[n] += dt * minv_res[n];
412  }
413 
414  // Add to the residuals
415  object_pt->add_to_dofs(B[i], k);
416  // Set the new time
417  object_pt->time() = initial_time + C[i] * dt;
418 
419  object_pt->actions_after_explicit_stage();
420  }
421 
422  object_pt->actions_after_explicit_timestep();
423  }
424 
425  // ??ds this could be heavily optimised if needed. Keeping it simple for
426  // now
428  const double& dt)
429  {
430  using namespace StringConversion;
431 
433  object_pt->actions_before_explicit_stage();
434 
435  // Storage indicies for the history values that we need
436  unsigned tn = 1;
437  unsigned tnm1 = tn + 1;
438  unsigned tnm2 = tnm1 + 1;
439 
440  // Check dts are the same, this will need to be removed if ebdf3 is being
441  // used as something other than a predictor... But seeing as it isn't
442  // stable that isn't likely.
443 #ifdef PARANOID
444  if (std::abs(dt - object_pt->time_pt()->dt(0)) > 1e-15)
445  {
446  std::string err =
447  "Inconsistent dts! Predictor is stepping by " + to_string(dt);
448  err += " but dt(0) = " + to_string(object_pt->time_pt()->dt(0));
449  throw OomphLibError(
450  err, OOMPH_EXCEPTION_LOCATION, OOMPH_CURRENT_FUNCTION);
451  }
452 #endif
453 
454  // Get older dt values
455  double dtnm1 = object_pt->time_pt()->dt(1);
456  double dtnm2 = object_pt->time_pt()->dt(2);
457 
458  // Calculate weights for these dts
459  set_weights(dt, dtnm1, dtnm2);
460 
461  // Get derivative value at step n (even though this uses values from t=0 =
462  // step n+1, it's ok because we haven't changed the values in that slot
463  // yet).
464  DoubleVector fn;
465  object_pt->get_dvaluesdt(fn);
466  fn *= Fn_weight;
467 
468  // Extract history values and multiply by their weights
469  DoubleVector ynp1, yn, ynm1, ynm2;
470  object_pt->get_dofs(tn, yn);
471  yn *= Yn_weight;
472  object_pt->get_dofs(tnm1, ynm1);
473  ynm1 *= Ynm1_weight;
474  object_pt->get_dofs(tnm2, ynm2);
475  ynm2 *= Ynm2_weight;
476 
477 
478  // Add everything together
479  ynp1 = yn;
480  ynp1 += ynm1;
481  ynp1 += ynm2;
482  ynp1 += fn;
483 
484  // Done, update things in the object
485  object_pt->set_dofs(ynp1);
486  object_pt->time() += dt;
487 
488  // Actions functions
489  object_pt->actions_after_explicit_stage();
490  object_pt->actions_after_explicit_timestep();
491  }
492 
493  /// Calculate the weights for this set of step sizes.
494  void EBDF3::set_weights(const double& dtn,
495  const double& dtnm1,
496  const double& dtnm2)
497  {
498  using namespace std;
499 
500  // If this is slow we can probably optimise by doing direct
501  // multiplication instead of using pow.
502 
503  // Generated using sympy from my python ode code.
504 
505  double denom = pow(dtnm1, 4) * dtnm2 + 2 * pow(dtnm1, 3) * pow(dtnm2, 2) +
506  pow(dtnm1, 2) * pow(dtnm2, 3);
507 
508  Yn_weight =
509  -(2 * pow(dtn, 3) * dtnm1 * dtnm2 + pow(dtn, 3) * pow(dtnm2, 2) +
510  3 * pow(dtn, 2) * pow(dtnm1, 2) * dtnm2 +
511  3 * pow(dtn, 2) * dtnm1 * pow(dtnm2, 2) + pow(dtn, 2) * pow(dtnm2, 3) -
512  pow(dtnm1, 4) * dtnm2 - 2 * pow(dtnm1, 3) * pow(dtnm2, 2) -
513  pow(dtnm1, 2) * pow(dtnm2, 3)) /
514  denom;
515 
516  Ynm1_weight =
517  -(-pow(dtn, 3) * pow(dtnm1, 2) - 2 * pow(dtn, 3) * dtnm1 * dtnm2 -
518  pow(dtn, 3) * pow(dtnm2, 2) - pow(dtn, 2) * pow(dtnm1, 3) -
519  3 * pow(dtn, 2) * pow(dtnm1, 2) * dtnm2 -
520  3 * pow(dtn, 2) * dtnm1 * pow(dtnm2, 2) - pow(dtn, 2) * pow(dtnm2, 3)) /
521  denom;
522 
523  Ynm2_weight =
524  -(pow(dtn, 3) * pow(dtnm1, 2) + pow(dtn, 2) * pow(dtnm1, 3)) / denom;
525 
526  Fn_weight =
527  -(-pow(dtn, 3) * pow(dtnm1, 2) * dtnm2 -
528  pow(dtn, 3) * dtnm1 * pow(dtnm2, 2) -
529  2 * pow(dtn, 2) * pow(dtnm1, 3) * dtnm2 -
530  3 * pow(dtn, 2) * pow(dtnm1, 2) * pow(dtnm2, 2) -
531  pow(dtn, 2) * dtnm1 * pow(dtnm2, 3) - dtn * pow(dtnm1, 4) * dtnm2 -
532  2 * dtn * pow(dtnm1, 3) * pow(dtnm2, 2) -
533  dtn * pow(dtnm1, 2) * pow(dtnm2, 3)) /
534  denom;
535  }
536 
537 
538  // Force build of templates
539  template class RungeKutta<4>;
540  template class LowStorageRungeKutta<4>;
541 } // namespace oomph
e
Definition: cfortran.h:571
cstr elem_len * i
Definition: cfortran.h:603
char t
Definition: cfortran.h:568
LinearAlgebraDistribution * distribution_pt() const
access to the LinearAlgebraDistribution
unsigned nrow() const
access function to the number of global rows.
A vector in the mathematical sense, initially developed for linear algebra type applications....
Definition: double_vector.h:58
void build(const DoubleVector &old_vector)
Just copys the argument DoubleVector.
void timestep(ExplicitTimeSteppableObject *const &object_pt, const double &dt)
Function that is used to advance the solution by time dt.
void set_weights(const double &dtn, const double &dtnm1, const double &dtnm2)
Calculate the weights for this set of step sizes.
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 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.
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.
=========================================================== Runge Kutta Timestepping that uses low st...
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.
An OomphLibError object which should be thrown when an run-time error is encountered....
=========================================================== Standard Runge Kutta Timestepping
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
double & dt(const unsigned &t=0)
Return the value of the t-th stored timestep (t=0: present; t>0: previous).
Definition: timesteppers.h:136
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
std::string to_string(T object, unsigned float_precision=8)
Conversion function that should work for anything with operator<< defined (at least all basic types).
//////////////////////////////////////////////////////////////////// ////////////////////////////////...