oomph::TR Class Reference

#include <trapezoid_rule.h>

+ Inheritance diagram for oomph::TR:

Public Member Functions

 TR (const bool &adaptive=false)
 
virtual ~TR ()
 Virtual destructor. More...
 
unsigned order () const
 Return the actual order of the scheme. More...
 
void set_weights ()
 Set the weights. More...
 
void set_error_weights ()
 
void set_predictor_weights ()
 Function to set the predictor weights. More...
 
unsigned nprev_values () const
 Number of previous values available. More...
 
unsigned derivative_index (const unsigned &t) const
 Location in data of derivatives. More...
 
unsigned predicted_value_index () const
 Location of predicted value. More...
 
unsigned ndt () const
 Number of timestep increments that need to be stored by the scheme. More...
 
void assign_initial_values_impulsive (Data *const &data_pt)
 
void assign_initial_positions_impulsive (Node *const &node_pt)
 
void actions_after_timestep (Problem *problem_pt)
 
void actions_before_timestep (Problem *problem_pt)
 
void setup_initial_derivative (Problem *problem_pt)
 
void shift_time_values (Data *const &data_pt)
 
void shift_time_positions (Node *const &node_pt)
 
void calculate_predicted_positions (Node *const &node_pt)
 Function to calculate predicted positions at a node. More...
 
void calculate_predicted_values (Data *const &data_pt)
 Function to calculate predicted data values in a Data object. More...
 
double temporal_error_in_position (Node *const &node_pt, const unsigned &i)
 Compute the error in the position i at a node. More...
 
double temporal_error_in_value (Data *const &data_pt, const unsigned &i)
 Compute the error in the value i in a Data structure. More...
 
- Public Member Functions inherited from oomph::TimeStepper
 TimeStepper (const unsigned &tstorage, const unsigned &max_deriv)
 
 TimeStepper ()
 Broken empty constructor. More...
 
 TimeStepper (const TimeStepper &)=delete
 Broken copy constructor. More...
 
void operator= (const TimeStepper &)=delete
 Broken assignment operator. More...
 
virtual ~TimeStepper ()
 virtual destructor More...
 
unsigned highest_derivative () const
 Highest order derivative that the scheme can compute. More...
 
doubletime ()
 Return current value of continous time. More...
 
double time () const
 Return current value of continous time. More...
 
virtual unsigned nprev_values_for_value_at_evaluation_time () const
 
void make_steady ()
 
bool is_steady () const
 
bool predict_by_explicit_step () const
 
ExplicitTimeStepperexplicit_predictor_pt ()
 
void set_predictor_pt (ExplicitTimeStepper *_pred_pt)
 
void update_predicted_time (const double &new_time)
 
void check_predicted_values_up_to_date () const
 Check that the predicted values are the ones we want. More...
 
unsigned predictor_storage_index () const
 
void enable_warning_in_assign_initial_data_values ()
 
void disable_warning_in_assign_initial_data_values ()
 
const DenseMatrix< double > * weights_pt () const
 Get a (const) pointer to the weights. More...
 
virtual void undo_make_steady ()
 
std::string type () const
 
void time_derivative (const unsigned &i, Data *const &data_pt, Vector< double > &deriv)
 
double time_derivative (const unsigned &i, Data *const &data_pt, const unsigned &j)
 Evaluate i-th derivative of j-th value in Data. More...
 
void time_derivative (const unsigned &i, Node *const &node_pt, Vector< double > &deriv)
 
double time_derivative (const unsigned &i, Node *const &node_pt, const unsigned &j)
 
Time *const & time_pt () const
 Access function for the pointer to time (const version) More...
 
Time *& time_pt ()
 
virtual double weight (const unsigned &i, const unsigned &j) const
 Access function for j-th weight for the i-th derivative. More...
 
unsigned ntstorage () const
 
bool adaptive_flag () const
 Function to indicate whether the scheme is adaptive (false by default) More...
 

Public Attributes

bool Initial_derivative_set
 
bool Shift_f
 

Private Member Functions

 TR (const TR &dummy)=delete
 Broken copy constructor. More...
 
void operator= (const TR &dummy)=delete
 Broken assignment operator. More...
 

Private Attributes

Vector< doublePredictor_weight
 Private data for the predictor weights. More...
 
double Error_weight
 Private data for the error weight. More...
 

Additional Inherited Members

- Protected Attributes inherited from oomph::TimeStepper
TimeTime_pt
 Pointer to discrete time storage scheme. More...
 
DenseMatrix< doubleWeight
 Storage for the weights associated with the timestepper. More...
 
std::string Type
 
bool Adaptive_Flag
 
bool Is_steady
 
bool Shut_up_in_assign_initial_data_values
 
bool Predict_by_explicit_step
 
ExplicitTimeStepperExplicit_predictor_pt
 
double Predicted_time
 
int Predictor_storage_index
 

Detailed Description

Trapezoid rule time stepping scheme.

This method requires a value of dy/dt at the initial time. The implementation of this calculation is exactly the same as is used for explicit time stepping.

The function setup_initial_derivative(Problem* problem_pt) should be called after the initial conditions have been set, but before beginning time stepping, to compute this initial value of dy/dt.

Warning: moving nodes not implemented (I have no test case).

Constructor & Destructor Documentation

◆ TR() [1/2]

oomph::TR::TR ( const bool adaptive = false)
inline

Constructor, storage for two history derivatives (one for TR and one for the predictor step), one history value, present value and predicted value.

63  : TimeStepper(2 + 2 + 1, 1)
64  {
65  // Set the weight for the zero-th derivative
66  Weight(0, 0) = 1.0;
67 
68  Initial_derivative_set = false;
69 
71 
72  // Initialise adaptivity stuff
73  Predictor_weight.assign(4, 0.0);
74  Error_weight = 0.0;
75 
76  Shift_f = true;
77  }
Vector< double > Predictor_weight
Private data for the predictor weights.
Definition: trapezoid_rule.h:331
bool Shift_f
Definition: trapezoid_rule.h:253
double Error_weight
Private data for the error weight.
Definition: trapezoid_rule.h:334
bool Initial_derivative_set
Definition: trapezoid_rule.h:251
DenseMatrix< double > Weight
Storage for the weights associated with the timestepper.
Definition: timesteppers.h:237
TimeStepper()
Broken empty constructor.
Definition: timesteppers.h:301
bool Adaptive_Flag
Definition: timesteppers.h:245
int adaptive
Definition: jeffery_hamel.cc:106

References Global_Physical_Variables::adaptive, oomph::TimeStepper::Adaptive_Flag, Error_weight, Initial_derivative_set, Predictor_weight, Shift_f, and oomph::TimeStepper::Weight.

◆ ~TR()

virtual oomph::TR::~TR ( )
inlinevirtual

Virtual destructor.

80 {}

◆ TR() [2/2]

oomph::TR::TR ( const TR dummy)
privatedelete

Broken copy constructor.

Member Function Documentation

◆ actions_after_timestep()

void oomph::TR::actions_after_timestep ( Problem problem_pt)
inlinevirtual

Interface for any actions that need to be performed after a time step.

Reimplemented from oomph::TimeStepper.

172  {
173  // only ever want this to be true for one step
174  Shift_f = true;
175  }

References Shift_f.

◆ actions_before_timestep()

void oomph::TR::actions_before_timestep ( Problem problem_pt)
inlinevirtual

Interface for any actions that need to be performed before a time step.

Reimplemented from oomph::TimeStepper.

178  {
179 #ifdef PARANOID
181  {
182  std::string err = "Initial derivative of TR has not been set";
183  throw OomphLibError(
185  }
186 #endif
187  }
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References Initial_derivative_set, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and oomph::Global_string_for_annotation::string().

◆ assign_initial_positions_impulsive()

void oomph::TR::assign_initial_positions_impulsive ( Node *const &  node_pt)
inlinevirtual

Initialise the time-history for the nodal positions corresponding to an impulsive start.

Implements oomph::TimeStepper.

165  {
166  throw OomphLibError("Function not yet implemented",
169  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ assign_initial_values_impulsive()

void oomph::TR::assign_initial_values_impulsive ( Data *const &  data_pt)
inlinevirtual

Initialise the time-history for the Data values, corresponding to an impulsive start.

Implements oomph::TimeStepper.

156  {
157  throw OomphLibError("Function not yet implemented",
160  }

References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

◆ calculate_predicted_positions()

void oomph::TR::calculate_predicted_positions ( Node *const &  node_pt)
inlinevirtual

Function to calculate predicted positions at a node.

Reimplemented from oomph::TimeStepper.

265  {
266  // Loop over the dimensions
267  unsigned n_dim = node_pt->ndim();
268  for (unsigned j = 0; j < n_dim; j++)
269  {
270  // If the node is not copied
271  if (!node_pt->position_is_a_copy(j))
272  {
273  // Initialise the predictor to zero
274  double predicted_value = 0.0;
275 
276  // Loop over all the stored data and add appropriate values
277  // to the predictor
278  for (unsigned i = 1; i < 4; i++)
279  {
280  predicted_value += node_pt->x(i, j) * Predictor_weight[i];
281  }
282 
283  // Store the predicted value
284  node_pt->x(predicted_value_index(), j) = predicted_value;
285  }
286  }
287  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
unsigned predicted_value_index() const
Location of predicted value.
Definition: trapezoid_rule.h:142
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References i, j, oomph::Node::ndim(), oomph::Node::position_is_a_copy(), predicted_value_index(), Predictor_weight, and oomph::Node::x().

◆ calculate_predicted_values()

void oomph::TR::calculate_predicted_values ( Data *const &  data_pt)
inlinevirtual

Function to calculate predicted data values in a Data object.

Reimplemented from oomph::TimeStepper.

291  {
292  // Loop over the values
293  unsigned n_value = data_pt->nvalue();
294  for (unsigned j = 0; j < n_value; j++)
295  {
296  // If the value is not copied
297  if (!data_pt->is_a_copy(j))
298  {
299  // Loop over all the stored data and add appropriate
300  // values to the predictor
301  double predicted_value = 0.0;
302  for (unsigned i = 1; i < 4; i++)
303  {
304  predicted_value += data_pt->value(i, j) * Predictor_weight[i];
305  }
306 
307  // Store the predicted value
308  data_pt->set_value(predicted_value_index(), j, predicted_value);
309  }
310  }
311  }

References i, oomph::Data::is_a_copy(), j, oomph::Data::nvalue(), predicted_value_index(), Predictor_weight, oomph::Data::set_value(), and oomph::Data::value().

◆ derivative_index()

unsigned oomph::TR::derivative_index ( const unsigned t) const
inline

Location in data of derivatives.

129  {
130 #ifdef PARANOID
131  if (t >= 2)
132  {
133  std::string err = "Only storing two derivatives.";
134  throw OomphLibError(
136  }
137 #endif
138  return t + 2;
139  }
t
Definition: plotPSD.py:36

References OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::Global_string_for_annotation::string(), and plotPSD::t.

Referenced by predicted_value_index(), set_predictor_weights(), set_weights(), setup_initial_derivative(), and shift_time_values().

◆ ndt()

unsigned oomph::TR::ndt ( ) const
inlinevirtual

Number of timestep increments that need to be stored by the scheme.

Implements oomph::TimeStepper.

149  {
150  return 2;
151  }

◆ nprev_values()

unsigned oomph::TR::nprev_values ( ) const
inlinevirtual

Number of previous values available.

Implements oomph::TimeStepper.

123  {
124  return 1;
125  }

◆ operator=()

void oomph::TR::operator= ( const TR dummy)
privatedelete

Broken assignment operator.

◆ order()

unsigned oomph::TR::order ( ) const
inlinevirtual

Return the actual order of the scheme.

Reimplemented from oomph::TimeStepper.

84  {
85  return 2;
86  }

◆ predicted_value_index()

unsigned oomph::TR::predicted_value_index ( ) const
inline

Location of predicted value.

143  {
144  return derivative_index(1) + 1;
145  }
unsigned derivative_index(const unsigned &t) const
Location in data of derivatives.
Definition: trapezoid_rule.h:128

References derivative_index().

Referenced by calculate_predicted_positions(), calculate_predicted_values(), temporal_error_in_position(), and temporal_error_in_value().

◆ set_error_weights()

void oomph::TR::set_error_weights ( )
inlinevirtual

Set the weights for the error computation, (currently empty – overwrite for specific scheme)

Reimplemented from oomph::TimeStepper.

98  {
99  double dt = Time_pt->dt(0);
100  double dtprev = Time_pt->dt(1);
101  Error_weight = 1 / (3 * (1 + (dtprev / dt)));
102  }
Time * Time_pt
Pointer to discrete time storage scheme.
Definition: timesteppers.h:234
double & dt(const unsigned &t=0)
Definition: timesteppers.h:136

References oomph::Time::dt(), Error_weight, and oomph::TimeStepper::Time_pt.

◆ set_predictor_weights()

void oomph::TR::set_predictor_weights ( )
inlinevirtual

Function to set the predictor weights.

Reimplemented from oomph::TimeStepper.

106  {
107  // Read the value of the time steps
108  double dt = Time_pt->dt(0);
109  double dtprev = Time_pt->dt(1);
110  double dtr = dt / dtprev;
111 
112  // y weights
113  Predictor_weight[0] = 0.0;
114  Predictor_weight[1] = 1.0;
115 
116  // dy weights
117  Predictor_weight[derivative_index(0)] = (dt / 2) * (2 + dtr);
118  Predictor_weight[derivative_index(1)] = -(dt / 2) * dtr;
119  }

References derivative_index(), oomph::Time::dt(), Predictor_weight, and oomph::TimeStepper::Time_pt.

◆ set_weights()

void oomph::TR::set_weights ( )
inlinevirtual

Set the weights.

Implements oomph::TimeStepper.

90  {
91  double dt = Time_pt->dt(0);
92  Weight(1, 0) = 2.0 / dt;
93  Weight(1, 1) = -2.0 / dt;
94  Weight(1, derivative_index(0)) = -1.0; // old derivative
95  }

References derivative_index(), oomph::Time::dt(), oomph::TimeStepper::Time_pt, and oomph::TimeStepper::Weight.

◆ setup_initial_derivative()

void oomph::TR::setup_initial_derivative ( Problem problem_pt)
inline
190  {
191  oomph_info
192  << "Solving for derivative at initial time."
193  << " Warning: if residual is not in the correct form this may fail."
194  << std::endl;
195 
196  // Get the derivative at initial time
197  DoubleVector f0;
198  problem_pt->get_dvaluesdt(f0);
199 
200  // store in derivatives slot ready for use in timestepping.
201  problem_pt->set_dofs(this->derivative_index(0), f0);
202 
203  Initial_derivative_set = true;
204  Shift_f = false;
205  }
std::vector< double > DoubleVector
loads clump configuration
Definition: ClumpInput.h:26
OomphInfo oomph_info
Definition: oomph_definitions.cc:319

References derivative_index(), oomph::Problem::get_dvaluesdt(), Initial_derivative_set, oomph::oomph_info, oomph::Problem::set_dofs(), and Shift_f.

Referenced by oomph::MyProblem::actions_after_set_initial_condition().

◆ shift_time_positions()

void oomph::TR::shift_time_positions ( Node *const &  node_pt)
inlinevirtual

This function advances the time history of the positions at a node.

Implements oomph::TimeStepper.

258  {
259  // do nothing, not implemented for moving nodes
260  }

◆ shift_time_values()

void oomph::TR::shift_time_values ( Data *const &  data_pt)
inlinevirtual

This function updates the Data's time history so that we can advance to the next timestep.

Implements oomph::TimeStepper.

211  {
212  // Find number of values stored
213  unsigned n_value = data_pt->nvalue();
214 
215  // Get previous step dt, time_pt has already been shifted so it's in
216  // slot 1.
217  double dtn = time_pt()->dt(1);
218 
219  // Loop over the values
220  for (unsigned j = 0; j < n_value; j++)
221  {
222  // Set previous values to the previous value, if not a copy
223  if (data_pt->is_a_copy(j) == false)
224  {
225  if (Shift_f)
226  {
227  // Calculate time derivative at step n by rearranging the TR
228  // formula.
229  double fnm1 = data_pt->value(derivative_index(0), j);
230  double ynm1 = data_pt->value(1, j);
231  double yn = data_pt->value(0, j);
232  double fn = (2 / dtn) * (yn - ynm1) - fnm1;
233 
234  data_pt->set_value(derivative_index(0), j, fn);
235 
236  // Shift time derivative
237  data_pt->set_value(derivative_index(1), j, fnm1);
238  }
239  else
240  {
241  std::cout << "didn't shift derivatives" << std::endl;
242  }
243 
244  // Shift value
245  data_pt->set_value(1, j, data_pt->value(0, j));
246  }
247  }
248  }
Time *const & time_pt() const
Access function for the pointer to time (const version)
Definition: timesteppers.h:572

References derivative_index(), oomph::Time::dt(), oomph::Data::is_a_copy(), j, oomph::Data::nvalue(), oomph::Data::set_value(), Shift_f, oomph::TimeStepper::time_pt(), and oomph::Data::value().

◆ temporal_error_in_position()

double oomph::TR::temporal_error_in_position ( Node *const &  node_pt,
const unsigned i 
)
inlinevirtual

Compute the error in the position i at a node.

Reimplemented from oomph::TimeStepper.

316  {
317  return Error_weight *
318  (node_pt->x(i) - node_pt->x(predicted_value_index(), i));
319  }

References Error_weight, i, predicted_value_index(), and oomph::Node::x().

◆ temporal_error_in_value()

double oomph::TR::temporal_error_in_value ( Data *const &  data_pt,
const unsigned i 
)
inlinevirtual

Compute the error in the value i in a Data structure.

Reimplemented from oomph::TimeStepper.

323  {
324  return Error_weight *
325  (data_pt->value(i) - data_pt->value(predicted_value_index(), i));
326  }

References Error_weight, i, predicted_value_index(), and oomph::Data::value().

Member Data Documentation

◆ Error_weight

double oomph::TR::Error_weight
private

Private data for the error weight.

Referenced by set_error_weights(), temporal_error_in_position(), temporal_error_in_value(), and TR().

◆ Initial_derivative_set

bool oomph::TR::Initial_derivative_set

◆ Predictor_weight

Vector<double> oomph::TR::Predictor_weight
private

Private data for the predictor weights.

Referenced by calculate_predicted_positions(), calculate_predicted_values(), set_predictor_weights(), and TR().

◆ Shift_f


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