oomph::EBDF3 Class Reference

#include <explicit_timesteppers.h>

+ Inheritance diagram for oomph::EBDF3:

Public Member Functions

 EBDF3 ()
 Constructor, set the type. More...
 
 EBDF3 (const EBDF3 &)=delete
 Broken copy constructor. More...
 
void operator= (const EBDF3 &)=delete
 Broken assignment operator. More...
 
void set_weights (const double &dtn, const double &dtnm1, const double &dtnm2)
 Calculate the weights for this set of step sizes. More...
 
void timestep (ExplicitTimeSteppableObject *const &object_pt, const double &dt)
 Function that is used to advance the solution by time dt. More...
 
- Public Member Functions inherited from oomph::ExplicitTimeStepper
 ExplicitTimeStepper ()
 Empty Constructor. More...
 
 ExplicitTimeStepper (const ExplicitTimeStepper &)=delete
 Broken copy constructor. More...
 
void operator= (const ExplicitTimeStepper &)=delete
 Broken assignment operator. More...
 
virtual ~ExplicitTimeStepper ()
 Empty virtual destructor — no memory is allocated in this class. More...
 

Private Attributes

double Yn_weight
 
double Ynm1_weight
 
double Ynm2_weight
 
double Fn_weight
 

Additional Inherited Members

- Protected Attributes inherited from oomph::ExplicitTimeStepper
std::string Type
 

Detailed Description

=========================================================== An explicit version of BDF3 (i.e. uses derivative evaluation at y_n instead of y_{n+1}). Useful as a predictor because it is third order accurate but requires only one function evaluation (i.e. only one mass matrix inversion + residual calculation).

Constructor & Destructor Documentation

◆ EBDF3() [1/2]

oomph::EBDF3::EBDF3 ( )
inline

Constructor, set the type.

248 {}

◆ EBDF3() [2/2]

oomph::EBDF3::EBDF3 ( const EBDF3 )
delete

Broken copy constructor.

Member Function Documentation

◆ operator=()

void oomph::EBDF3::operator= ( const EBDF3 )
delete

Broken assignment operator.

◆ set_weights()

void oomph::EBDF3::set_weights ( const double dtn,
const double dtnm1,
const double dtnm2 
)

Calculate the weights for this set of step sizes.

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  }
double Yn_weight
Definition: explicit_timesteppers.h:241
double Ynm2_weight
Definition: explicit_timesteppers.h:243
double Ynm1_weight
Definition: explicit_timesteppers.h:242
double Fn_weight
Definition: explicit_timesteppers.h:244
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 pow(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:625

References Fn_weight, Eigen::bfloat16_impl::pow(), Yn_weight, Ynm1_weight, and Ynm2_weight.

Referenced by timestep().

◆ timestep()

void oomph::EBDF3::timestep ( ExplicitTimeSteppableObject *const &  object_pt,
const double dt 
)
virtual

Function that is used to advance the solution by time dt.

Implements oomph::ExplicitTimeStepper.

429  {
430  using namespace StringConversion;
431 
432  object_pt->actions_before_explicit_timestep();
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(
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  }
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
std::vector< double > DoubleVector
loads clump configuration
Definition: ClumpInput.h:26
Array< double, 1, 3 > e(1./3., 0.5, 2.)
void set_weights(const double &dtn, const double &dtnm1, const double &dtnm2)
Calculate the weights for this set of step sizes.
Definition: explicit_timesteppers.cc:494
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286
std::string to_string(T object, unsigned float_precision=8)
Definition: oomph_utilities.h:189
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References abs(), oomph::ExplicitTimeSteppableObject::actions_after_explicit_stage(), oomph::ExplicitTimeSteppableObject::actions_after_explicit_timestep(), oomph::ExplicitTimeSteppableObject::actions_before_explicit_stage(), oomph::ExplicitTimeSteppableObject::actions_before_explicit_timestep(), oomph::Time::dt(), e(), Fn_weight, oomph::ExplicitTimeSteppableObject::get_dofs(), oomph::ExplicitTimeSteppableObject::get_dvaluesdt(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::ExplicitTimeSteppableObject::set_dofs(), set_weights(), oomph::Global_string_for_annotation::string(), oomph::ExplicitTimeSteppableObject::time(), oomph::ExplicitTimeSteppableObject::time_pt(), oomph::StringConversion::to_string(), Yn_weight, Ynm1_weight, and Ynm2_weight.

Member Data Documentation

◆ Fn_weight

double oomph::EBDF3::Fn_weight
private

Referenced by set_weights(), and timestep().

◆ Yn_weight

double oomph::EBDF3::Yn_weight
private

Referenced by set_weights(), and timestep().

◆ Ynm1_weight

double oomph::EBDF3::Ynm1_weight
private

Referenced by set_weights(), and timestep().

◆ Ynm2_weight

double oomph::EBDF3::Ynm2_weight
private

Referenced by set_weights(), and timestep().


The documentation for this class was generated from the following files: