Public Member Functions | Protected Attributes | List of all members
oomph::TimeStepper Class Referenceabstract

Base class for time-stepping schemes. Timestepper provides an approximation of the temporal derivatives of Data such that the i-th derivative of the j-th value in Data is represented as. More...

#include <timesteppers.h>

+ Inheritance diagram for oomph::TimeStepper:

Public Member Functions

 TimeStepper (const unsigned &tstorage, const unsigned &max_deriv)
 Constructor. Pass the amount of storage required by timestepper (present value + history values) and the order of highest time-derivative. More...
 
 TimeStepper ()
 Broken empty constructor. More...
 
 TimeStepper (const TimeStepper &)
 Broken copy constructor. More...
 
void operator= (const TimeStepper &)
 Broken assignment operator. More...
 
virtual ~TimeStepper ()
 virtual destructor More...
 
unsigned highest_derivative () const
 Highest order derivative that the scheme can compute. More...
 
virtual unsigned order () const
 Actual order (accuracy) of the scheme. More...
 
double & time ()
 Return current value of continous time. More...
 
double time () const
 Return current value of continous time. More...
 
virtual unsigned ndt () const =0
 Number of timestep increments that are required by the scheme. More...
 
virtual unsigned nprev_values_for_value_at_evaluation_time () const
 Number of previous values needed to calculate the value at the current time. i.e. how many previous values must we loop over to calculate the values at the evaluation time. For most methods this is 1, i.e. just use the value at t_{n+1}. See midpoint method for a counter-example. More...
 
virtual unsigned nprev_values () const =0
 Number of previous values available: 0 for static, 1 for BDF<1>,... More...
 
virtual void set_weights ()=0
 Function to set the weights for present timestep (don't need to pass present timestep or previous timesteps as they are available via Time_pt) More...
 
void make_steady ()
 Function to make the time stepper temporarily steady. This is trivially achieved by setting all the weights to zero. More...
 
bool is_steady () const
 Flag to indicate if a timestepper has been made steady (possibly temporarily to switch off time-dependence) More...
 
bool predict_by_explicit_step () const
 Flag: is adaptivity done by taking a separate step using an ExplicitTimeStepper object? More...
 
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
 Return the time-index in each Data where predicted values are stored if the timestepper is adaptive. More...
 
void enable_warning_in_assign_initial_data_values ()
 Enable the output of warnings due to possible fct pointer vector size mismatch in assign_initial_data_values (Default) More...
 
void disable_warning_in_assign_initial_data_values ()
 Disable the output of warnings due to possible fct pointer vector size mismatch in assign_initial_data_values. More...
 
const DenseMatrix< double > * weights_pt () const
 Get a (const) pointer to the weights. More...
 
virtual void undo_make_steady ()
 Reset the is_steady status of a specific TimeStepper to its default and re-assign the weights. More...
 
std::string type () const
 Return string that indicates the type of the timestepper (e.g. "BDF", "Newmark", etc.) More...
 
void time_derivative (const unsigned &i, Data *const &data_pt, Vector< double > &deriv)
 Evaluate i-th derivative of all values in Data and return in Vector deriv[]. More...
 
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)
 Evaluate i-th derivative of all values in Node and return in Vector deriv[] (this can't be simply combined with time_derivative(.., Data, ...) because of differences with haning nodes). More...
 
double time_derivative (const unsigned &i, Node *const &node_pt, const unsigned &j)
 Evaluate i-th derivative of j-th value in Node. Note the use of the node's value() function so that hanging nodes are taken into account (this is why the functions for Data and Node cannot be combined through simple polymorphism: value is not virtual). More...
 
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
 Return the number of doubles required to represent history (one for steady) More...
 
virtual void assign_initial_values_impulsive (Data *const &data_pt)=0
 Initialise the time-history for the Data values corresponding to an impulsive start. More...
 
virtual void assign_initial_positions_impulsive (Node *const &node_pt)=0
 Initialiset the positions for the node corresponding to an impulsive start. More...
 
virtual void shift_time_values (Data *const &data_pt)=0
 This function advances the Data's time history so that we can move on to the next timestep. More...
 
virtual void shift_time_positions (Node *const &node_pt)=0
 This function advances the time history of the positions at a node. The default should be OK, but would need to be overloaded. More...
 
bool adaptive_flag () const
 Function to indicate whether the scheme is adaptive (false by default) More...
 
virtual void set_predictor_weights ()
 Set the weights for the predictor previous timestep (currently empty – overwrite for specific scheme) More...
 
virtual void calculate_predicted_values (Data *const &data_pt)
 Do the predictor step for data stored in a Data object (currently empty – overwrite for specific scheme) More...
 
virtual void calculate_predicted_positions (Node *const &node_pt)
 Do the predictor step for the positions at a node (currently empty — overwrite for a specific scheme) More...
 
virtual void set_error_weights ()
 Set the weights for the error computation, (currently empty – overwrite for specific scheme) More...
 
virtual double temporal_error_in_position (Node *const &node_pt, const unsigned &i)
 
virtual double temporal_error_in_value (Data *const &data_pt, const unsigned &i)
 
virtual void actions_before_timestep (Problem *problem_pt)
 
virtual void actions_after_timestep (Problem *problem_pt)
 

Protected Attributes

TimeTime_pt
 Pointer to discrete time storage scheme. More...
 
DenseMatrix< double > Weight
 Storage for the weights associated with the timestepper. More...
 
std::string Type
 String that indicates the type of the timestepper (e.g. "BDF", "Newmark", etc.) More...
 
bool Adaptive_Flag
 Boolean variable to indicate whether the timestepping scheme can be adaptive. More...
 
bool Is_steady
 Bool to indicate if the timestepper is steady, i.e. its time-derivatives evaluate to zero. This status may be achieved temporarily by calling make_steady(). It can be reset to the appropriate default by the function undo_make_steady(). More...
 
bool Shut_up_in_assign_initial_data_values
 Boolean to indicate if the timestepper will output warnings when setting possibly an incorrect number of initial data values from function pointers. More...
 
bool Predict_by_explicit_step
 Flag: is adaptivity done by taking a separate step using an ExplicitTimeStepper object? More...
 
ExplicitTimeStepperExplicit_predictor_pt
 
double Predicted_time
 
int Predictor_storage_index
 The time-index in each Data object where predicted values are stored. -1 if not set. More...
 

Detailed Description

Base class for time-stepping schemes. Timestepper provides an approximation of the temporal derivatives of Data such that the i-th derivative of the j-th value in Data is represented as.

unsigned n_value=data_pt->nvalue();
Vector<double> time_derivative(n_value);
for (unsigned j=0;j<n_value;j++)
{
for (unsigned t=0;t<ntstorage();t++)
{
time_derivative[j] += data_pt->value(t,j)*weight(i,t);
}
}

where the time index t is such that

TimeSteppers can evaluate all derivatives up to their order().

The first nprev_values() `history values' represent the values at the previous timesteps. For BDF schemes, nprev_values() = ntstorage()-1 , i.e. all `history values' represent previous values. For Newmark<NSTEPS>, only the first NSTEPS `history values' represent previous values (NSTEPS=1 gives the normal Newmark scheme).

Definition at line 219 of file timesteppers.h.

Constructor & Destructor Documentation

◆ TimeStepper() [1/3]

oomph::TimeStepper::TimeStepper ( const unsigned &  tstorage,
const unsigned &  max_deriv 
)
inline

Constructor. Pass the amount of storage required by timestepper (present value + history values) and the order of highest time-derivative.

Definition at line 270 of file timesteppers.h.

References oomph::DenseMatrix< T >::resize().

◆ TimeStepper() [2/3]

oomph::TimeStepper::TimeStepper ( )
inline

Broken empty constructor.

Definition at line 292 of file timesteppers.h.

◆ TimeStepper() [3/3]

oomph::TimeStepper::TimeStepper ( const TimeStepper )
inline

Broken copy constructor.

Definition at line 300 of file timesteppers.h.

References oomph::BrokenCopy::broken_copy().

◆ ~TimeStepper()

oomph::TimeStepper::~TimeStepper ( )
virtual

virtual destructor

Destructor for time stepper, in .cc file so that we don't have to include explicit_timesteppers.h in header.

Definition at line 42 of file timesteppers.cc.

References oomph::Steady< NSTEPS >::Dummy_time, and Explicit_predictor_pt.

Member Function Documentation

◆ actions_after_timestep()

virtual void oomph::TimeStepper::actions_after_timestep ( Problem problem_pt)
inlinevirtual

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

Reimplemented in oomph::IMRByBDF, and oomph::TR.

Definition at line 618 of file timesteppers.h.

Referenced by oomph::Problem::adaptive_unsteady_newton_solve(), oomph::IMRByBDF::nprev_values_for_value_at_evaluation_time(), and oomph::Problem::unsteady_newton_solve().

◆ actions_before_timestep()

virtual void oomph::TimeStepper::actions_before_timestep ( Problem problem_pt)
inlinevirtual

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

Reimplemented in oomph::IMRByBDF, and oomph::TR.

Definition at line 614 of file timesteppers.h.

Referenced by oomph::Problem::adaptive_unsteady_newton_solve(), oomph::IMRByBDF::nprev_values_for_value_at_evaluation_time(), and oomph::Problem::unsteady_newton_solve().

◆ adaptive_flag()

bool oomph::TimeStepper::adaptive_flag ( ) const
inline

◆ assign_initial_positions_impulsive()

virtual void oomph::TimeStepper::assign_initial_positions_impulsive ( Node *const &  node_pt)
pure virtual

◆ assign_initial_values_impulsive()

virtual void oomph::TimeStepper::assign_initial_values_impulsive ( Data *const &  data_pt)
pure virtual

◆ calculate_predicted_positions()

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

Do the predictor step for the positions at a node (currently empty — overwrite for a specific scheme)

Reimplemented in oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::TR, oomph::BDF< NSTEPS >, and oomph::IMRBase.

Definition at line 594 of file timesteppers.h.

◆ calculate_predicted_values()

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

Do the predictor step for data stored in a Data object (currently empty – overwrite for specific scheme)

Reimplemented in oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::TR, oomph::BDF< NSTEPS >, oomph::IMRBase, and oomph::BDF< NSTEPS >.

Definition at line 590 of file timesteppers.h.

◆ check_predicted_values_up_to_date()

void oomph::TimeStepper::check_predicted_values_up_to_date ( ) const
inline

Check that the predicted values are the ones we want.

Definition at line 406 of file timesteppers.h.

References e, and oomph::Time::time().

Referenced by oomph::IMRBase::calculate_predicted_values().

◆ disable_warning_in_assign_initial_data_values()

void oomph::TimeStepper::disable_warning_in_assign_initial_data_values ( )
inline

Disable the output of warnings due to possible fct pointer vector size mismatch in assign_initial_data_values.

Definition at line 449 of file timesteppers.h.

◆ enable_warning_in_assign_initial_data_values()

void oomph::TimeStepper::enable_warning_in_assign_initial_data_values ( )
inline

Enable the output of warnings due to possible fct pointer vector size mismatch in assign_initial_data_values (Default)

Definition at line 444 of file timesteppers.h.

◆ explicit_predictor_pt()

ExplicitTimeStepper* oomph::TimeStepper::explicit_predictor_pt ( )
inline

Get the pointer to the explicit timestepper to use as a predictor in adaptivity if Predict_by_explicit_step is set.

Definition at line 389 of file timesteppers.h.

Referenced by oomph::Problem::calculate_predictions().

◆ highest_derivative()

unsigned oomph::TimeStepper::highest_derivative ( ) const
inline

Highest order derivative that the scheme can compute.

Definition at line 315 of file timesteppers.h.

References oomph::DenseMatrix< T >::nrow().

◆ is_steady()

bool oomph::TimeStepper::is_steady ( ) const
inline

Flag to indicate if a timestepper has been made steady (possibly temporarily to switch off time-dependence)

Definition at line 375 of file timesteppers.h.

Referenced by oomph::Problem::arc_length_step_solve_helper(), oomph::PoroelasticityEquations< 2 >::d2u_dt2(), oomph::AxisymmetricPoroelasticityEquations::d2u_dt2(), oomph::AxisymmetricLinearElasticityEquationsBase::d2u_dt2_axisymmetric_linear_elasticity(), oomph::LinearWaveEquations< DIM >::d2u_dt2_lin_wave(), oomph::LinearElasticityEquationsBase< DIM >::d2u_dt2_linear_elasticity(), oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::dc_dt_adv_diff_react(), oomph::FSIAxisymFoepplvonKarmanElement< NNODE_1D, FLUID_ELEMENT >::dposition_dt(), oomph::Node::dposition_dt(), oomph::Node::dposition_gen_dt(), oomph::PoroelasticityEquations< 2 >::dq_edge_dt(), oomph::AxisymmetricPoroelasticityEquations::dq_edge_dt(), oomph::PoroelasticityEquations< 2 >::dq_internal_dt(), oomph::AxisymmetricPoroelasticityEquations::dq_internal_dt(), oomph::PoroelasticityEquations< 2 >::du_dt(), oomph::AxisymmetricPoroelasticityEquations::du_dt(), oomph::AdvectionDiffusionEquations< DIM >::du_dt_adv_diff(), oomph::AxisymAdvectionDiffusionEquations::du_dt_axi_adv_diff(), oomph::AxisymmetricNavierStokesEquations::du_dt_axi_nst(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::du_dt_axi_nst(), oomph::AxisymmetricLinearElasticityEquationsBase::du_dt_axisymmetric_linear_elasticity(), oomph::GeneralisedAdvectionDiffusionEquations< DIM >::du_dt_cons_adv_diff(), oomph::du_dt_cons_axisym_adv_diff(), oomph::FluxTransportEquations< DIM >::du_dt_flux_transport(), oomph::LinearWaveEquations< DIM >::du_dt_lin_wave(), oomph::LinearisedAxisymmetricNavierStokesEquations::du_dt_linearised_axi_nst(), oomph::LinearisedNavierStokesEquations::du_dt_linearised_nst(), oomph::GeneralisedNewtonianNavierStokesEquations< DIM >::du_dt_nst(), oomph::NavierStokesEquations< DIM >::du_dt_nst(), oomph::SphericalAdvectionDiffusionEquations::du_dt_spherical_adv_diff(), oomph::UnsteadyHeatEquations< DIM >::du_dt_ust_heat(), oomph::WomersleyEquations< DIM >::du_dt_womersley(), oomph::Node::dx_dt(), oomph::Node::dx_gen_dt(), oomph::Problem::get_dvaluesdt(), oomph::Problem::solve_eigenproblem(), oomph::Problem::steady_newton_solve(), and oomph::SegregatableFSIProblem::steady_segregated_solve().

◆ make_steady()

void oomph::TimeStepper::make_steady ( )
inline

Function to make the time stepper temporarily steady. This is trivially achieved by setting all the weights to zero.

Definition at line 360 of file timesteppers.h.

References oomph::DenseMatrix< T >::initialise().

Referenced by oomph::Problem::arc_length_step_solve_helper(), oomph::Problem::get_dvaluesdt(), oomph::Problem::solve_eigenproblem(), oomph::Problem::steady_newton_solve(), and oomph::SegregatableFSIProblem::steady_segregated_solve().

◆ ndt()

virtual unsigned oomph::TimeStepper::ndt ( ) const
pure virtual

◆ nprev_values()

virtual unsigned oomph::TimeStepper::nprev_values ( ) const
pure virtual

Number of previous values available: 0 for static, 1 for BDF<1>,...

Implemented in oomph::BDF< NSTEPS >, oomph::Newmark< NSTEPS >, oomph::Steady< NSTEPS >, oomph::Steady< 0 >, oomph::ContinuationStorageScheme, oomph::PeriodicOrbitTimeDiscretisation, oomph::TR, and oomph::IMRBase.

Referenced by oomph::AlgebraicFSIDrivenCavityMesh< ELEMENT >::algebraic_node_update(), oomph::AlgebraicCollapsibleChannelMesh< ELEMENT >::algebraic_node_update(), oomph::Circle::Circle(), oomph::RefineablePolarTaylorHoodElement::get_interpolated_values(), oomph::RefineablePolarCrouzeixRaviartElement::get_interpolated_values(), oomph::FiniteElement::get_x(), oomph::MacroElementNodeUpdateNode::node_update(), oomph::AlgebraicRefineableQuarterTubeMesh< ELEMENT >::node_update_central_region(), oomph::AlgebraicRefineableQuarterCircleSectorMesh< ELEMENT >::node_update_in_central_box(), oomph::AlgebraicRefineableQuarterCircleSectorMesh< ELEMENT >::node_update_in_lower_right_box(), oomph::AlgebraicRefineableQuarterCircleSectorMesh< ELEMENT >::node_update_in_upper_left_box(), oomph::AlgebraicRefineableQuarterTubeMesh< ELEMENT >::node_update_lower_right_region(), oomph::AlgebraicRefineableQuarterTubeMesh< ELEMENT >::node_update_upper_left_region(), oomph::PseudoBucklingRing::position(), oomph::PseudoBucklingRing::PseudoBucklingRing(), oomph::ImmersedRigidBodyTriangleMeshPolygon::reset_reference_configuration(), oomph::TetMeshBase::snap_nodes_onto_geometric_objects(), and oomph::UnstructuredTwoDMeshGeometryBase::snap_nodes_onto_geometric_objects().

◆ nprev_values_for_value_at_evaluation_time()

virtual unsigned oomph::TimeStepper::nprev_values_for_value_at_evaluation_time ( ) const
inlinevirtual

Number of previous values needed to calculate the value at the current time. i.e. how many previous values must we loop over to calculate the values at the evaluation time. For most methods this is 1, i.e. just use the value at t_{n+1}. See midpoint method for a counter-example.

Reimplemented in oomph::IMRByBDF, oomph::IMR, and oomph::IMRBase.

Definition at line 348 of file timesteppers.h.

◆ ntstorage()

unsigned oomph::TimeStepper::ntstorage ( ) const
inline

Return the number of doubles required to represent history (one for steady)

Definition at line 562 of file timesteppers.h.

References oomph::DenseMatrix< T >::ncol().

Referenced by oomph::Node::add_values_to_vector(), oomph::PeriodicOrbitTimeDiscretisation::assign_initial_data_values(), oomph::RefineableSolidQElement< 2 >::build(), oomph::RefineableSolidQElement< 3 >::build(), oomph::Multi_domain_functions::construct_new_external_halo_master_node_helper(), oomph::Multi_domain_functions::construct_new_external_halo_node_helper(), oomph::RefineableTriangleMesh< ELEMENT >::construct_new_halo_node_helper(), oomph::RefineableTriangleMesh< ELEMENT >::construct_new_node_load_balance_helper(), oomph::Node::copy(), oomph::PoroelasticityEquations< 2 >::d2u_dt2(), oomph::AxisymmetricPoroelasticityEquations::d2u_dt2(), oomph::AxisymmetricLinearElasticityEquationsBase::d2u_dt2_axisymmetric_linear_elasticity(), oomph::LinearWaveEquations< DIM >::d2u_dt2_lin_wave(), oomph::LinearElasticityEquationsBase< DIM >::d2u_dt2_linear_elasticity(), oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::dc_dt_adv_diff_react(), oomph::SurfactantTransportInterfaceElement::dcdt_surface(), oomph::FSIAxisymFoepplvonKarmanElement< NNODE_1D, FLUID_ELEMENT >::dposition_dt(), oomph::Node::dposition_dt(), oomph::Node::dposition_gen_dt(), oomph::PoroelasticityEquations< 2 >::dq_edge_dt(), oomph::AxisymmetricPoroelasticityEquations::dq_edge_dt(), oomph::PoroelasticityEquations< 2 >::dq_internal_dt(), oomph::AxisymmetricPoroelasticityEquations::dq_internal_dt(), oomph::PoroelasticityEquations< 2 >::du_dt(), oomph::AxisymmetricPoroelasticityEquations::du_dt(), oomph::AdvectionDiffusionEquations< DIM >::du_dt_adv_diff(), oomph::AxisymAdvectionDiffusionEquations::du_dt_axi_adv_diff(), oomph::AxisymmetricNavierStokesEquations::du_dt_axi_nst(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::du_dt_axi_nst(), oomph::AxisymmetricLinearElasticityEquationsBase::du_dt_axisymmetric_linear_elasticity(), oomph::GeneralisedAdvectionDiffusionEquations< DIM >::du_dt_cons_adv_diff(), oomph::du_dt_cons_axisym_adv_diff(), oomph::FluxTransportEquations< DIM >::du_dt_flux_transport(), oomph::LinearWaveEquations< DIM >::du_dt_lin_wave(), oomph::LinearisedAxisymmetricNavierStokesEquations::du_dt_linearised_axi_nst(), oomph::LinearisedNavierStokesEquations::du_dt_linearised_nst(), oomph::GeneralisedNewtonianNavierStokesEquations< DIM >::du_dt_nst(), oomph::NavierStokesEquations< DIM >::du_dt_nst(), oomph::PolarNavierStokesEquations::du_dt_pnst(), oomph::SphericalAdvectionDiffusionEquations::du_dt_spherical_adv_diff(), oomph::SphericalNavierStokesEquations::du_dt_spherical_nst(), oomph::UnsteadyHeatEquations< DIM >::du_dt_ust_heat(), oomph::WomersleyEquations< DIM >::du_dt_womersley(), oomph::Node::dump(), oomph::Node::dx_dt(), oomph::Node::dx_gen_dt(), oomph::SolidFiniteElement::fill_in_jacobian_for_newmark_accel(), oomph::Missing_masters_functions::get_required_master_nodal_information_helper(), oomph::Multi_domain_functions::get_required_master_nodal_information_helper(), oomph::Missing_masters_functions::get_required_nodal_information_helper(), oomph::Multi_domain_functions::get_required_nodal_information_helper(), oomph::ContinuationStorageScheme::modify_storage(), oomph::PeriodicOrbitTimeDiscretisation::ndt(), oomph::GenericLagrangeInterpolatedProjectableElement< ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableTimeHarmonicLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableDarcyElement< DARCY_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableLinearElasticityElement< LINEAR_ELAST_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectablePMLTimeHarmonicLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableTimeHarmonicFourierDecomposedLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableAxisymLinearElasticityElement< AXISYM_LINEAR_ELAST_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectablePVDElement< PVD_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableFoepplvonKarmanElement< FVK_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableUnsteadyHeatElement< UNSTEADY_HEAT_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableHelmholtzElement< HELMHOLTZ_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableDisplacementBasedFoepplvonKarmanElement< FVK_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectablePoissonElement< POISSON_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectablePMLHelmholtzElement< HELMHOLTZ_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableAxisymmetricPoroelasticityElement< AXISYMMETRIC_POROELASTICITY_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableAxisymmetricTaylorHoodElement< TAYLOR_HOOD_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::GeneralisedNewtonianProjectableAxisymmetricTaylorHoodElement< TAYLOR_HOOD_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableAxisymmetricCrouzeixRaviartElement< CROUZEIX_RAVIART_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::GeneralisedNewtonianProjectableAxisymmetricCrouzeixRaviartElement< CROUZEIX_RAVIART_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectablePVDElementWithContinuousPressure< PVD_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableGeneralisedNewtonianTaylorHoodElement< TAYLOR_HOOD_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableGeneralisedNewtonianCrouzeixRaviartElement< CROUZEIX_RAVIART_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableTaylorHoodElement< TAYLOR_HOOD_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableCrouzeixRaviartElement< CROUZEIX_RAVIART_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::Node::Node(), oomph::PeriodicOrbitTimeDiscretisation::nprev_values(), oomph::Data::ntstorage(), oomph::ProjectableUnsteadyHeatElement< UNSTEADY_HEAT_ELEMENT >::output(), oomph::PRefineableQElement< 1, INITIAL_NNODE_1D >::p_refine(), oomph::PRefineableQElement< 2, INITIAL_NNODE_1D >::p_refine(), oomph::PRefineableQElement< 3, INITIAL_NNODE_1D >::p_refine(), oomph::Node::read(), oomph::Node::read_values_from_vector(), oomph::RefineableQSpectralElement< 3 >::rebuild_from_sons(), oomph::RefineableQSpectralElement< 2 >::rebuild_from_sons(), oomph::RefineableQSpectralElement< 1 >::rebuild_from_sons(), oomph::PRefineableQElement< 1, INITIAL_NNODE_1D >::rebuild_from_sons(), oomph::PRefineableQElement< 2, INITIAL_NNODE_1D >::rebuild_from_sons(), oomph::PRefineableQElement< 3, INITIAL_NNODE_1D >::rebuild_from_sons(), oomph::Node::set_position_time_stepper(), oomph::Data::set_time_stepper(), oomph::PeriodicOrbitEquations::set_timestepper_weights(), oomph::IMRBase::shift_time_positions(), oomph::Newmark< NSTEPS >::shift_time_positions(), oomph::NewmarkBDF< NSTEPS >::shift_time_positions(), and oomph::NewmarkBDF< NSTEPS >::shift_time_values().

◆ operator=()

void oomph::TimeStepper::operator= ( const TimeStepper )
inline

Broken assignment operator.

Definition at line 306 of file timesteppers.h.

References oomph::BrokenCopy::broken_assign().

◆ order()

virtual unsigned oomph::TimeStepper::order ( ) const
inlinevirtual

◆ predict_by_explicit_step()

bool oomph::TimeStepper::predict_by_explicit_step ( ) const
inline

Flag: is adaptivity done by taking a separate step using an ExplicitTimeStepper object?

Definition at line 382 of file timesteppers.h.

Referenced by oomph::Problem::calculate_predictions().

◆ predictor_storage_index()

unsigned oomph::TimeStepper::predictor_storage_index ( ) const
inline

Return the time-index in each Data where predicted values are stored if the timestepper is adaptive.

Definition at line 420 of file timesteppers.h.

References oomph::Global_string_for_annotation::string().

◆ set_error_weights()

virtual void oomph::TimeStepper::set_error_weights ( )
inlinevirtual

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

Reimplemented in oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::IMRBase, and oomph::TR.

Definition at line 598 of file timesteppers.h.

Referenced by oomph::Problem::adaptive_unsteady_newton_solve(), and oomph::Problem::initialise_dt().

◆ set_predictor_pt()

void oomph::TimeStepper::set_predictor_pt ( ExplicitTimeStepper _pred_pt)
inline

Set the pointer to the explicit timestepper to use as a predictor in adaptivity if Predict_by_explicit_step is set.

Definition at line 393 of file timesteppers.h.

◆ set_predictor_weights()

virtual void oomph::TimeStepper::set_predictor_weights ( )
inlinevirtual

Set the weights for the predictor previous timestep (currently empty – overwrite for specific scheme)

Reimplemented in oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::IMRBase, oomph::TR, and oomph::BDF< NSTEPS >.

Definition at line 586 of file timesteppers.h.

Referenced by oomph::Problem::adaptive_unsteady_newton_solve().

◆ set_weights()

virtual void oomph::TimeStepper::set_weights ( )
pure virtual

◆ shift_time_positions()

virtual void oomph::TimeStepper::shift_time_positions ( Node *const &  node_pt)
pure virtual

This function advances the time history of the positions at a node. The default should be OK, but would need to be overloaded.

Implemented in oomph::BDF< NSTEPS >, oomph::NewmarkBDF< NSTEPS >, oomph::Newmark< NSTEPS >, oomph::Steady< NSTEPS >, oomph::Steady< 0 >, oomph::TR, oomph::ContinuationStorageScheme, oomph::PeriodicOrbitTimeDiscretisation, and oomph::IMRBase.

◆ shift_time_values()

virtual void oomph::TimeStepper::shift_time_values ( Data *const &  data_pt)
pure virtual

◆ temporal_error_in_position()

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

Compute the error in the position i at a node zero here – overwrite for specific scheme.

Reimplemented in oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::TR, oomph::BDF< NSTEPS >, and oomph::IMRBase.

Definition at line 602 of file timesteppers.h.

◆ temporal_error_in_value()

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

Compute the error in the value i in a Data structure zero here – overwrite for specific scheme.

Reimplemented in oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::TR, oomph::BDF< NSTEPS >, and oomph::IMRBase.

Definition at line 608 of file timesteppers.h.

◆ time() [1/2]

double& oomph::TimeStepper::time ( )
inline

◆ time() [2/2]

double oomph::TimeStepper::time ( ) const
inline

Return current value of continous time.

Definition at line 327 of file timesteppers.h.

References oomph::Time::ndt(), oomph::Global_string_for_annotation::string(), and oomph::Time::time().

◆ time_derivative() [1/4]

void oomph::TimeStepper::time_derivative ( const unsigned &  i,
Data *const &  data_pt,
Vector< double > &  deriv 
)
inline

◆ time_derivative() [2/4]

double oomph::TimeStepper::time_derivative ( const unsigned &  i,
Data *const &  data_pt,
const unsigned &  j 
)
inline

Evaluate i-th derivative of j-th value in Data.

Definition at line 489 of file timesteppers.h.

References t, and oomph::Data::value().

◆ time_derivative() [3/4]

void oomph::TimeStepper::time_derivative ( const unsigned &  i,
Node *const &  node_pt,
Vector< double > &  deriv 
)
inline

Evaluate i-th derivative of all values in Node and return in Vector deriv[] (this can't be simply combined with time_derivative(.., Data, ...) because of differences with haning nodes).

Definition at line 506 of file timesteppers.h.

References oomph::Data::nvalue().

◆ time_derivative() [4/4]

double oomph::TimeStepper::time_derivative ( const unsigned &  i,
Node *const &  node_pt,
const unsigned &  j 
)
inline

Evaluate i-th derivative of j-th value in Node. Note the use of the node's value() function so that hanging nodes are taken into account (this is why the functions for Data and Node cannot be combined through simple polymorphism: value is not virtual).

Definition at line 524 of file timesteppers.h.

References t, and oomph::Node::value().

◆ time_pt() [1/2]

Time* const& oomph::TimeStepper::time_pt ( ) const
inline

Access function for the pointer to time (const version)

Definition at line 540 of file timesteppers.h.

References oomph::Global_string_for_annotation::string().

Referenced by oomph::PseudoBucklingRing::accel(), oomph::IMRByBDF::actions_before_timestep(), oomph::Problem::add_time_stepper_pt(), oomph::TimeHarmonicLinearElasticityEquationsBase< DIM >::body_force(), oomph::LinearElasticityEquationsBase< DIM >::body_force(), oomph::PVDEquationsBase< DIM >::body_force(), oomph::GeneralisedNewtonianNavierStokesEquations< DIM >::extrapolated_strain_rate(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::extrapolated_strain_rate(), oomph::AxisymmetricLinearElasticityEquations::fill_in_generic_contribution_to_residuals_axisymmetric_linear_elasticity(), oomph::PolarNavierStokesTractionElement< ELEMENT >::fill_in_generic_residual_contribution(), oomph::AxisymmetricPoroelasticityEquations::fill_in_generic_residual_contribution(), oomph::RefineableAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::RefineableGeneralisedNewtonianAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::AxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::NavierStokesTractionElement< ELEMENT >::fill_in_generic_residual_contribution_fluid_traction(), oomph::LinearWaveFluxElement< ELEMENT >::fill_in_generic_residual_contribution_lin_wave_flux(), oomph::SphericalNavierStokesEquations::fill_in_generic_residual_contribution_spherical_nst(), oomph::RefineableSphericalNavierStokesEquations::fill_in_generic_residual_contribution_spherical_nst(), oomph::UnsteadyHeatFluxElement< ELEMENT >::fill_in_generic_residual_contribution_ust_heat_flux(), oomph::RefineableAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::RefineableGeneralisedNewtonianAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::AxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::get_wind_adv_diff_react(), oomph::PseudoBucklingRing::position(), oomph::RefineableNavierStokesTractionElement< ELEMENT >::refineable_fill_in_generic_residual_contribution_fluid_traction(), oomph::SolidICProblem::set_static_initial_condition(), oomph::TR::shift_time_values(), and oomph::PseudoBucklingRing::veloc().

◆ time_pt() [2/2]

Time*& oomph::TimeStepper::time_pt ( )
inline

Access function for the pointer to time (can't have a paranoid test for null pointers because this is used as a set function).

Definition at line 554 of file timesteppers.h.

◆ type()

std::string oomph::TimeStepper::type ( ) const
inline

◆ undo_make_steady()

virtual void oomph::TimeStepper::undo_make_steady ( )
inlinevirtual

◆ update_predicted_time()

void oomph::TimeStepper::update_predicted_time ( const double &  new_time)
inline

Set the time that the current predictions apply for, only needed for paranoid checks when doing Predict_by_explicit_step.

Definition at line 400 of file timesteppers.h.

Referenced by oomph::Problem::calculate_predictions().

◆ weight()

virtual double oomph::TimeStepper::weight ( const unsigned &  i,
const unsigned &  j 
) const
inlinevirtual

Access function for j-th weight for the i-th derivative.

Reimplemented in oomph::Steady< NSTEPS >, and oomph::Steady< 0 >.

Definition at line 557 of file timesteppers.h.

Referenced by oomph::SurfactantTransportInterfaceElement::add_additional_residual_contributions_interface(), oomph::PoroelasticityEquations< 2 >::d2u_dt2(), oomph::AxisymmetricPoroelasticityEquations::d2u_dt2(), oomph::AxisymmetricLinearElasticityEquationsBase::d2u_dt2_axisymmetric_linear_elasticity(), oomph::LinearWaveEquations< DIM >::d2u_dt2_lin_wave(), oomph::LinearElasticityEquationsBase< DIM >::d2u_dt2_linear_elasticity(), oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::dc_dt_adv_diff_react(), oomph::SurfactantTransportInterfaceElement::dcdt_surface(), oomph::FSIAxisymFoepplvonKarmanElement< NNODE_1D, FLUID_ELEMENT >::dposition_dt(), oomph::Node::dposition_dt(), oomph::Node::dposition_gen_dt(), oomph::PoroelasticityEquations< 2 >::dq_edge_dt(), oomph::AxisymmetricPoroelasticityEquations::dq_edge_dt(), oomph::PoroelasticityEquations< 2 >::dq_internal_dt(), oomph::AxisymmetricPoroelasticityEquations::dq_internal_dt(), oomph::PoroelasticityEquations< 2 >::du_dt(), oomph::AxisymmetricPoroelasticityEquations::du_dt(), oomph::AdvectionDiffusionEquations< DIM >::du_dt_adv_diff(), oomph::AxisymAdvectionDiffusionEquations::du_dt_axi_adv_diff(), oomph::AxisymmetricNavierStokesEquations::du_dt_axi_nst(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::du_dt_axi_nst(), oomph::AxisymmetricLinearElasticityEquationsBase::du_dt_axisymmetric_linear_elasticity(), oomph::GeneralisedAdvectionDiffusionEquations< DIM >::du_dt_cons_adv_diff(), oomph::du_dt_cons_axisym_adv_diff(), oomph::FluxTransportEquations< DIM >::du_dt_flux_transport(), oomph::LinearWaveEquations< DIM >::du_dt_lin_wave(), oomph::LinearisedAxisymmetricNavierStokesEquations::du_dt_linearised_axi_nst(), oomph::LinearisedNavierStokesEquations::du_dt_linearised_nst(), oomph::GeneralisedNewtonianNavierStokesEquations< DIM >::du_dt_nst(), oomph::NavierStokesEquations< DIM >::du_dt_nst(), oomph::PolarNavierStokesEquations::du_dt_pnst(), oomph::SphericalAdvectionDiffusionEquations::du_dt_spherical_adv_diff(), oomph::SphericalNavierStokesEquations::du_dt_spherical_nst(), oomph::UnsteadyHeatEquations< DIM >::du_dt_ust_heat(), oomph::WomersleyEquations< DIM >::du_dt_womersley(), oomph::Node::dx_dt(), oomph::Node::dx_gen_dt(), oomph::ODEElement::fill_in_contribution_to_jacobian(), oomph::AxisymmetricLinearElasticityEquations::fill_in_generic_contribution_to_residuals_axisymmetric_linear_elasticity(), oomph::AxisymmetricNavierStokesEquations::fill_in_generic_dresidual_contribution_axi_nst(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::fill_in_generic_dresidual_contribution_axi_nst(), oomph::PoroelasticityEquations< 2 >::fill_in_generic_residual_contribution(), oomph::AxisymmetricPoroelasticityEquations::fill_in_generic_residual_contribution(), oomph::RefineableAxisymAdvectionDiffusionEquations::fill_in_generic_residual_contribution_axi_adv_diff(), oomph::AxisymAdvectionDiffusionEquations::fill_in_generic_residual_contribution_axi_adv_diff(), oomph::RefineableAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::RefineableGeneralisedNewtonianAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::AxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::RefineableGeneralisedAxisymAdvectionDiffusionEquations::fill_in_generic_residual_contribution_cons_axisym_adv_diff(), oomph::LinearisedAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_linearised_axi_nst(), oomph::RefineableLinearisedAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_linearised_axi_nst(), oomph::RefineableSphericalAdvectionDiffusionEquations::fill_in_generic_residual_contribution_spherical_adv_diff(), oomph::SphericalAdvectionDiffusionEquations::fill_in_generic_residual_contribution_spherical_adv_diff(), oomph::SphericalNavierStokesEquations::fill_in_generic_residual_contribution_spherical_nst(), oomph::RefineableSphericalNavierStokesEquations::fill_in_generic_residual_contribution_spherical_nst(), oomph::SolidFiniteElement::fill_in_jacobian_for_newmark_accel(), oomph::RefineableAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::RefineableGeneralisedNewtonianAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::AxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::ImmersedRigidBodyElement::get_residuals_rigid_body_generic(), oomph::Newmark< NSTEPS >::shift_time_positions(), oomph::NewmarkBDF< NSTEPS >::shift_time_positions(), and oomph::NewmarkBDF< NSTEPS >::shift_time_values().

◆ weights_pt()

const DenseMatrix<double>* oomph::TimeStepper::weights_pt ( ) const
inline

Get a (const) pointer to the weights.

Definition at line 453 of file timesteppers.h.

Member Data Documentation

◆ Adaptive_Flag

bool oomph::TimeStepper::Adaptive_Flag
protected

Boolean variable to indicate whether the timestepping scheme can be adaptive.

Definition at line 235 of file timesteppers.h.

Referenced by oomph::IMRBase::IMRBase(), and oomph::TR::TR().

◆ Explicit_predictor_pt

ExplicitTimeStepper* oomph::TimeStepper::Explicit_predictor_pt
protected

Pointer to explicit time stepper to use as predictor if Predict_by_explicit_step is set. ??ds merge the two? predict by explicit if pointer is set?

Definition at line 255 of file timesteppers.h.

Referenced by ~TimeStepper().

◆ Is_steady

bool oomph::TimeStepper::Is_steady
protected

Bool to indicate if the timestepper is steady, i.e. its time-derivatives evaluate to zero. This status may be achieved temporarily by calling make_steady(). It can be reset to the appropriate default by the function undo_make_steady().

Definition at line 241 of file timesteppers.h.

Referenced by oomph::ContinuationStorageScheme::ContinuationStorageScheme(), oomph::IMRBase::IMRBase(), and oomph::ContinuationStorageScheme::undo_make_steady().

◆ Predict_by_explicit_step

bool oomph::TimeStepper::Predict_by_explicit_step
protected

Flag: is adaptivity done by taking a separate step using an ExplicitTimeStepper object?

Definition at line 250 of file timesteppers.h.

Referenced by oomph::IMRBase::IMRBase().

◆ Predicted_time

double oomph::TimeStepper::Predicted_time
protected

Store the time that the predicted values currently stored are at, to compare for paranoid checks.

Definition at line 259 of file timesteppers.h.

◆ Predictor_storage_index

int oomph::TimeStepper::Predictor_storage_index
protected

◆ Shut_up_in_assign_initial_data_values

bool oomph::TimeStepper::Shut_up_in_assign_initial_data_values
protected

Boolean to indicate if the timestepper will output warnings when setting possibly an incorrect number of initial data values from function pointers.

Definition at line 246 of file timesteppers.h.

Referenced by oomph::Newmark< NSTEPS >::assign_initial_data_values().

◆ Time_pt

Time* oomph::TimeStepper::Time_pt
protected

◆ Type

std::string oomph::TimeStepper::Type
protected

String that indicates the type of the timestepper (e.g. "BDF", "Newmark", etc.)

Definition at line 231 of file timesteppers.h.

Referenced by oomph::ContinuationStorageScheme::ContinuationStorageScheme(), oomph::IMRBase::IMRBase(), and oomph::PeriodicOrbitTimeDiscretisation::PeriodicOrbitTimeDiscretisation().

◆ Weight

DenseMatrix<double> oomph::TimeStepper::Weight
protected

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