////////////////////////////////////////////////////////////////////// ////////////////////////////////////////////////////////////////////// More...
#include <timesteppers.h>
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 &)=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... | |
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... | |
ExplicitTimeStepper * | explicit_predictor_pt () |
Get the pointer to the explicit timestepper to use as a predictor in adaptivity if Predict_by_explicit_step is set. More... | |
void | set_predictor_pt (ExplicitTimeStepper *_pred_pt) |
Set the pointer to the explicit timestepper to use as a predictor in adaptivity if Predict_by_explicit_step is set. More... | |
void | update_predicted_time (const double &new_time) |
Set the time that the current predictions apply for, only needed for paranoid checks when doing Predict_by_explicit_step. More... | |
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 () |
Access function for the pointer to time (can't have a paranoid test for null pointers because this is used as a set function). More... | |
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) |
Compute the error in the position i at a node zero here – overwrite for specific scheme. More... | |
virtual double | temporal_error_in_value (Data *const &data_pt, const unsigned &i) |
Compute the error in the value i in a Data structure zero here – overwrite for specific scheme. More... | |
virtual void | actions_before_timestep (Problem *problem_pt) |
Interface for any actions that need to be performed before a time step. More... | |
virtual void | actions_after_timestep (Problem *problem_pt) |
Interface for any actions that need to be performed after a time step. More... | |
Protected Attributes | |
Time * | Time_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... | |
ExplicitTimeStepper * | Explicit_predictor_pt |
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? More... | |
double | Predicted_time |
Store the time that the predicted values currently stored are at, to compare for paranoid checks. More... | |
int | Predictor_storage_index |
The time-index in each Data object where predicted values are stored. -1 if not set. More... | |
////////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////////
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
where the time index t
is such that
t
=
0
refers to the current value of the Datat
>
0
refers to the ‘history’ values, i.e. the doubles that are stored in Data to represent the value's time derivatives. For BDF schemes these doubles are the values at previous times; for other timestepping schemes they can represent quantities such as previous first and second derivatives, etc.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 230 of file timesteppers.h.
|
inline |
Constructor. Pass the amount of storage required by timestepper (present value + history values) and the order of highest time-derivative.
Definition at line 279 of file timesteppers.h.
References oomph::DenseMatrix< T >::resize().
|
inline |
Broken empty constructor.
Definition at line 301 of file timesteppers.h.
|
delete |
Broken copy constructor.
|
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 37 of file timesteppers.cc.
References Explicit_predictor_pt.
|
inlinevirtual |
Interface for any actions that need to be performed after a time step.
Reimplemented in oomph::TR, and oomph::IMRByBDF.
Definition at line 666 of file timesteppers.h.
Referenced by oomph::Problem::adaptive_unsteady_newton_solve(), and oomph::Problem::unsteady_newton_solve().
|
inlinevirtual |
Interface for any actions that need to be performed before a time step.
Reimplemented in oomph::TR, and oomph::IMRByBDF.
Definition at line 662 of file timesteppers.h.
Referenced by oomph::Problem::adaptive_unsteady_newton_solve(), and oomph::Problem::unsteady_newton_solve().
|
inline |
Function to indicate whether the scheme is adaptive (false by default)
Definition at line 623 of file timesteppers.h.
Referenced by oomph::IMRBase::calculate_predicted_values(), oomph::IMRBase::shift_time_positions(), and oomph::IMRBase::temporal_error_in_value().
|
pure virtual |
Initialiset the positions for the node corresponding to an impulsive start.
Implemented in oomph::TR, oomph::BDF< NSTEPS >, oomph::Newmark< NSTEPS >, oomph::Steady< NSTEPS >, oomph::Steady< 0 >, oomph::PeriodicOrbitTimeDiscretisation, oomph::IMRBase, and oomph::ContinuationStorageScheme.
|
pure virtual |
Initialise the time-history for the Data values corresponding to an impulsive start.
Implemented in oomph::TR, oomph::BDF< NSTEPS >, oomph::Newmark< NSTEPS >, oomph::Steady< NSTEPS >, oomph::Steady< 0 >, oomph::PeriodicOrbitTimeDiscretisation, oomph::IMRBase, and oomph::ContinuationStorageScheme.
Referenced by oomph::Circle::Circle().
|
inlinevirtual |
Do the predictor step for the positions at a node (currently empty — overwrite for a specific scheme)
Reimplemented in oomph::TR, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, and oomph::IMRBase.
Definition at line 638 of file timesteppers.h.
|
inlinevirtual |
Do the predictor step for data stored in a Data object (currently empty – overwrite for specific scheme)
Reimplemented in oomph::TR, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, and oomph::IMRBase.
Definition at line 634 of file timesteppers.h.
|
inline |
Check that the predicted values are the ones we want.
Definition at line 423 of file timesteppers.h.
References e, and oomph::Time::time().
Referenced by oomph::IMRBase::calculate_predicted_values().
|
inline |
Disable the output of warnings due to possible fct pointer vector size mismatch in assign_initial_data_values.
Definition at line 469 of file timesteppers.h.
|
inline |
Enable the output of warnings due to possible fct pointer vector size mismatch in assign_initial_data_values (Default)
Definition at line 462 of file timesteppers.h.
|
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 403 of file timesteppers.h.
Referenced by oomph::Problem::calculate_predictions().
|
inline |
Highest order derivative that the scheme can compute.
Definition at line 318 of file timesteppers.h.
References oomph::DenseMatrix< T >::nrow().
|
inline |
Flag to indicate if a timestepper has been made steady (possibly temporarily to switch off time-dependence)
Definition at line 389 of file timesteppers.h.
Referenced by oomph::Problem::arc_length_step_solve_helper(), oomph::AxisymmetricPoroelasticityEquations::d2u_dt2(), oomph::PoroelasticityEquations< DIM >::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::Node::dposition_dt(), oomph::FSIAxisymFoepplvonKarmanElement< NNODE_1D, FLUID_ELEMENT >::dposition_dt(), oomph::Node::dposition_gen_dt(), oomph::AxisymmetricPoroelasticityEquations::dq_edge_dt(), oomph::PoroelasticityEquations< DIM >::dq_edge_dt(), oomph::AxisymmetricPoroelasticityEquations::dq_internal_dt(), oomph::PoroelasticityEquations< DIM >::dq_internal_dt(), oomph::AxisymmetricPoroelasticityEquations::du_dt(), oomph::PoroelasticityEquations< DIM >::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::GeneralisedAxisymAdvectionDiffusionEquations::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_adjoint_eigenproblem(), oomph::Problem::solve_adjoint_eigenproblem_legacy(), oomph::Problem::solve_eigenproblem(), oomph::Problem::solve_eigenproblem_legacy(), oomph::Problem::steady_newton_solve(), and oomph::SegregatableFSIProblem::steady_segregated_solve().
|
inline |
Function to make the time stepper temporarily steady. This is trivially achieved by setting all the weights to zero.
Definition at line 374 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_adjoint_eigenproblem(), oomph::Problem::solve_adjoint_eigenproblem_legacy(), oomph::Problem::solve_eigenproblem(), oomph::Problem::solve_eigenproblem_legacy(), oomph::Problem::steady_newton_solve(), and oomph::SegregatableFSIProblem::steady_segregated_solve().
|
pure virtual |
Number of timestep increments that are required by the scheme.
Implemented in oomph::TR, oomph::BDF< NSTEPS >, oomph::Newmark< NSTEPS >, oomph::Steady< NSTEPS >, oomph::Steady< 0 >, oomph::PeriodicOrbitTimeDiscretisation, oomph::IMRBase, and oomph::ContinuationStorageScheme.
Referenced by oomph::Problem::add_time_stepper_pt().
|
pure virtual |
Number of previous values available: 0 for static, 1 for BDF<1>,...
Implemented in oomph::TR, oomph::BDF< NSTEPS >, oomph::Newmark< NSTEPS >, oomph::Steady< NSTEPS >, oomph::Steady< 0 >, oomph::PeriodicOrbitTimeDiscretisation, oomph::IMRBase, and oomph::ContinuationStorageScheme.
Referenced by oomph::AlgebraicCollapsibleChannelMesh< ELEMENT >::algebraic_node_update(), oomph::AlgebraicFSIDrivenCavityMesh< ELEMENT >::algebraic_node_update(), oomph::RefineablePolarTaylorHoodElement::get_interpolated_values(), oomph::RefineablePolarCrouzeixRaviartElement::get_interpolated_values(), oomph::FiniteElement::get_x(), oomph::AlgebraicNode::node_update(), 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::SolidICProblem::set_newmark_initial_condition_consistently(), oomph::TetMeshBase::snap_nodes_onto_geometric_objects(), and oomph::UnstructuredTwoDMeshGeometryBase::snap_nodes_onto_geometric_objects().
|
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::IMRBase, oomph::IMRByBDF, and oomph::IMR.
Definition at line 359 of file timesteppers.h.
|
inline |
Return the number of doubles required to represent history (one for steady)
Definition at line 601 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< 3 >::build(), oomph::RefineableSolidQElement< 2 >::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::AxisymmetricPoroelasticityEquations::d2u_dt2(), oomph::PoroelasticityEquations< DIM >::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::Node::dposition_dt(), oomph::FSIAxisymFoepplvonKarmanElement< NNODE_1D, FLUID_ELEMENT >::dposition_dt(), oomph::Node::dposition_gen_dt(), oomph::AxisymmetricPoroelasticityEquations::dq_edge_dt(), oomph::PoroelasticityEquations< DIM >::dq_edge_dt(), oomph::AxisymmetricPoroelasticityEquations::dq_internal_dt(), oomph::PoroelasticityEquations< DIM >::dq_internal_dt(), oomph::AxisymmetricPoroelasticityEquations::du_dt(), oomph::PoroelasticityEquations< DIM >::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::GeneralisedAxisymAdvectionDiffusionEquations::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::ProjectableAdvectionDiffusionReactionElement< ADR_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableAxisymLinearElasticityElement< AXISYM_LINEAR_ELAST_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableAxisymmetricTaylorHoodElement< TAYLOR_HOOD_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableAxisymmetricCrouzeixRaviartElement< CROUZEIX_RAVIART_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableAxisymmetricPoroelasticityElement< AXISYMMETRIC_POROELASTICITY_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableDarcyElement< DARCY_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableDisplacementBasedFoepplvonKarmanElement< FVK_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableFoepplvonKarmanElement< FVK_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableFourierDecomposedHelmholtzElement< FOURIER_DECOMPOSED_HELMHOLTZ_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::GeneralisedNewtonianProjectableAxisymmetricTaylorHoodElement< TAYLOR_HOOD_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::GeneralisedNewtonianProjectableAxisymmetricCrouzeixRaviartElement< CROUZEIX_RAVIART_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::GenericLagrangeInterpolatedProjectableElement< ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableHelmholtzElement< HELMHOLTZ_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableLinearElasticityElement< LINEAR_ELAST_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::ProjectablePMLFourierDecomposedHelmholtzElement< FOURIER_DECOMPOSED_HELMHOLTZ_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectablePMLHelmholtzElement< HELMHOLTZ_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectablePMLTimeHarmonicLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectablePoissonElement< POISSON_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectablePVDElement< PVD_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectablePVDElementWithContinuousPressure< PVD_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableUnsteadyHeatSpaceTimeElement< UNSTEADY_HEAT_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableUnsteadyHeatMixedOrderSpaceTimeElement< UNSTEADY_HEAT_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableTimeHarmonicFourierDecomposedLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableTimeHarmonicLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::ProjectableUnsteadyHeatElement< UNSTEADY_HEAT_ELEMENT >::nhistory_values_for_coordinate_projection(), oomph::Node::Node(), oomph::PeriodicOrbitTimeDiscretisation::nprev_values(), oomph::Data::ntstorage(), oomph::ProjectableUnsteadyHeatSpaceTimeElement< UNSTEADY_HEAT_ELEMENT >::output(), oomph::ProjectableUnsteadyHeatMixedOrderSpaceTimeElement< UNSTEADY_HEAT_ELEMENT >::output(), 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::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::RefineableQSpectralElement< 3 >::rebuild_from_sons(), oomph::RefineableQSpectralElement< 1 >::rebuild_from_sons(), oomph::RefineableQSpectralElement< 2 >::rebuild_from_sons(), oomph::SolidICProblem::set_newmark_initial_condition_consistently(), oomph::Node::set_position_time_stepper(), oomph::Data::set_time_stepper(), oomph::PeriodicOrbitEquations::set_timestepper_weights(), oomph::IMRBase::shift_time_positions(), and oomph::Node::x_gen_range_check().
|
delete |
Broken assignment operator.
|
inlinevirtual |
Actual order (accuracy) of the scheme.
Reimplemented in oomph::TR, oomph::BDF< NSTEPS >, oomph::Newmark< NSTEPS >, oomph::Steady< NSTEPS >, oomph::Steady< 0 >, oomph::PeriodicOrbitTimeDiscretisation, oomph::IMRBase, and oomph::ContinuationStorageScheme.
Definition at line 324 of file timesteppers.h.
|
inline |
Flag: is adaptivity done by taking a separate step using an ExplicitTimeStepper object?
Definition at line 396 of file timesteppers.h.
|
inline |
Return the time-index in each Data where predicted values are stored if the timestepper is adaptive.
Definition at line 438 of file timesteppers.h.
References oomph::Global_string_for_annotation::string().
|
inlinevirtual |
Set the weights for the error computation, (currently empty – overwrite for specific scheme)
Reimplemented in oomph::TR, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, and oomph::IMRBase.
Definition at line 642 of file timesteppers.h.
Referenced by oomph::Problem::adaptive_unsteady_newton_solve(), and oomph::Problem::initialise_dt().
|
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 410 of file timesteppers.h.
|
inlinevirtual |
Set the weights for the predictor previous timestep (currently empty – overwrite for specific scheme)
Reimplemented in oomph::TR, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, and oomph::IMRBase.
Definition at line 630 of file timesteppers.h.
Referenced by oomph::Problem::adaptive_unsteady_newton_solve().
|
pure virtual |
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)
Implemented in oomph::IMRBase, oomph::TR, oomph::BDF< NSTEPS >, oomph::NewmarkBDF< NSTEPS >, oomph::Newmark< NSTEPS >, oomph::Steady< NSTEPS >, oomph::Steady< 0 >, oomph::NewmarkBDF< NSTEPS >, oomph::NewmarkBDF< NSTEPS >, oomph::NewmarkBDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::PeriodicOrbitTimeDiscretisation, oomph::IMRByBDF, oomph::IMR, and oomph::ContinuationStorageScheme.
Referenced by oomph::Problem::adaptive_unsteady_newton_solve(), oomph::Problem::copy(), oomph::Problem::initialise_dt(), oomph::Problem::unsteady_newton_solve(), and oomph::SegregatableFSIProblem::unsteady_segregated_solve().
|
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::TR, oomph::BDF< NSTEPS >, oomph::NewmarkBDF< NSTEPS >, oomph::Newmark< NSTEPS >, oomph::Steady< NSTEPS >, oomph::Steady< 0 >, oomph::PeriodicOrbitTimeDiscretisation, oomph::IMRBase, and oomph::ContinuationStorageScheme.
|
pure virtual |
This function advances the Data's time history so that we can move on to the next timestep.
Implemented in oomph::TR, oomph::BDF< NSTEPS >, oomph::NewmarkBDF< NSTEPS >, oomph::Newmark< NSTEPS >, oomph::Steady< NSTEPS >, oomph::Steady< 0 >, oomph::PeriodicOrbitTimeDiscretisation, oomph::IMRBase, and oomph::ContinuationStorageScheme.
|
inlinevirtual |
Compute the error in the position i at a node zero here – overwrite for specific scheme.
Reimplemented in oomph::TR, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, and oomph::IMRBase.
Definition at line 646 of file timesteppers.h.
|
inlinevirtual |
Compute the error in the value i in a Data structure zero here – overwrite for specific scheme.
Reimplemented in oomph::TR, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, oomph::BDF< NSTEPS >, and oomph::IMRBase.
Definition at line 654 of file timesteppers.h.
|
inline |
Return current value of continous time.
Definition at line 332 of file timesteppers.h.
References oomph::Time::time().
Referenced by oomph::PeriodicOrbitTimeDiscretisation::assign_initial_data_values(), oomph::ODEElement::fill_in_contribution_to_jacobian(), oomph::ODEElement::fill_in_contribution_to_residuals(), oomph::LinearisedAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_linearised_axi_nst(), oomph::RefineableLinearisedAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_linearised_axi_nst(), oomph::LinearisedNavierStokesEquations::fill_in_generic_residual_contribution_linearised_nst(), oomph::RefineableLinearisedNavierStokesEquations::fill_in_generic_residual_contribution_linearised_nst(), oomph::ImmersedRigidBodyElement::get_residuals_rigid_body_generic(), and oomph::ImmersedRigidBodyElement::output_centre_of_gravity().
|
inline |
Return current value of continous time.
Definition at line 338 of file timesteppers.h.
References oomph::Global_string_for_annotation::string(), and oomph::Time::time().
|
inline |
Evaluate i-th derivative of j-th value in Data.
Definition at line 518 of file timesteppers.h.
References i, t, and oomph::Data::value().
|
inline |
Evaluate i-th derivative of all values in Data and return in Vector deriv[].
Definition at line 502 of file timesteppers.h.
References i, and oomph::Data::nvalue().
Referenced by oomph::ImmersedRigidBodyElement::dposition_dt(), oomph::ODEElement::fill_in_contribution_to_residuals(), oomph::ImmersedRigidBodyElement::get_residuals_rigid_body_generic(), and oomph::ImmersedRigidBodyElement::output_centre_of_gravity().
|
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 555 of file timesteppers.h.
References i, t, and oomph::Node::value().
|
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 536 of file timesteppers.h.
References i, and oomph::Data::nvalue().
|
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 588 of file timesteppers.h.
|
inline |
Access function for the pointer to time (const version)
Definition at line 572 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::LinearElasticityEquationsBase< DIM >::body_force(), oomph::TimeHarmonicLinearElasticityEquationsBase< DIM >::body_force(), oomph::PVDEquationsBase< DIM >::body_force(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::extrapolated_strain_rate(), oomph::GeneralisedNewtonianNavierStokesEquations< DIM >::extrapolated_strain_rate(), oomph::AxisymmetricLinearElasticityEquations::fill_in_generic_contribution_to_residuals_axisymmetric_linear_elasticity(), oomph::AxisymmetricPoroelasticityEquations::fill_in_generic_residual_contribution(), oomph::AxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::RefineableAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::RefineableGeneralisedNewtonianAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::RefineableSphericalNavierStokesEquations::fill_in_generic_residual_contribution_spherical_nst(), oomph::SphericalNavierStokesEquations::fill_in_generic_residual_contribution_spherical_nst(), oomph::AxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::RefineableAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::RefineableGeneralisedNewtonianAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::get_wind_adv_diff_react(), oomph::AxisymmetricLinearElasticityTractionElement< ELEMENT >::output(), oomph::AxisymmetricNavierStokesTractionElement< ELEMENT >::output(), oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::output(), oomph::FSILinearisedAxisymPoroelasticTractionElement< POROELASTICITY_BULK_ELEMENT, NAVIER_STOKES_BULK_ELEMENT >::output(), oomph::PseudoBucklingRing::position(), oomph::AxisymmetricNavierStokesTractionElement< ELEMENT >::scalar_value_paraview(), oomph::SolidICProblem::set_newmark_initial_condition_consistently(), oomph::SolidICProblem::set_newmark_initial_condition_directly(), oomph::SolidICProblem::set_static_initial_condition(), oomph::TR::shift_time_values(), and oomph::PseudoBucklingRing::veloc().
|
inline |
Return string that indicates the type of the timestepper (e.g. "BDF", "Newmark", etc.)
Definition at line 490 of file timesteppers.h.
Referenced by oomph::SurfactantTransportInterfaceElement::dcdt_surface(), oomph::PolarNavierStokesEquations::du_dt_pnst(), oomph::SphericalNavierStokesEquations::du_dt_spherical_nst(), and oomph::SolidFiniteElement::fill_in_jacobian_for_newmark_accel().
|
inlinevirtual |
Reset the is_steady status of a specific TimeStepper to its default and re-assign the weights.
Reimplemented in oomph::ContinuationStorageScheme.
Definition at line 482 of file timesteppers.h.
Referenced by oomph::Problem::arc_length_step_solve_helper(), oomph::Problem::get_dvaluesdt(), oomph::Problem::solve_adjoint_eigenproblem(), oomph::Problem::solve_adjoint_eigenproblem_legacy(), oomph::Problem::solve_eigenproblem(), oomph::Problem::solve_eigenproblem_legacy(), oomph::Problem::steady_newton_solve(), and oomph::SegregatableFSIProblem::steady_segregated_solve().
|
inline |
Set the time that the current predictions apply for, only needed for paranoid checks when doing Predict_by_explicit_step.
Definition at line 417 of file timesteppers.h.
Referenced by oomph::Problem::calculate_predictions().
|
inlinevirtual |
Access function for j-th weight for the i-th derivative.
Reimplemented in oomph::Steady< NSTEPS >, and oomph::Steady< 0 >.
Definition at line 594 of file timesteppers.h.
References i.
Referenced by oomph::SurfactantTransportInterfaceElement::add_additional_residual_contributions_interface(), oomph::AxisymmetricPoroelasticityEquations::d2u_dt2(), oomph::PoroelasticityEquations< DIM >::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::Node::dposition_dt(), oomph::FSIAxisymFoepplvonKarmanElement< NNODE_1D, FLUID_ELEMENT >::dposition_dt(), oomph::Node::dposition_gen_dt(), oomph::AxisymmetricPoroelasticityEquations::dq_edge_dt(), oomph::PoroelasticityEquations< DIM >::dq_edge_dt(), oomph::AxisymmetricPoroelasticityEquations::dq_internal_dt(), oomph::PoroelasticityEquations< DIM >::dq_internal_dt(), oomph::AxisymmetricPoroelasticityEquations::du_dt(), oomph::PoroelasticityEquations< DIM >::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::GeneralisedAxisymAdvectionDiffusionEquations::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::AxisymmetricPoroelasticityEquations::fill_in_generic_residual_contribution(), oomph::AxisymAdvectionDiffusionEquations::fill_in_generic_residual_contribution_axi_adv_diff(), oomph::RefineableAxisymAdvectionDiffusionEquations::fill_in_generic_residual_contribution_axi_adv_diff(), oomph::AxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::RefineableAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::RefineableGeneralisedNewtonianAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::GeneralisedAxisymAdvectionDiffusionEquations::fill_in_generic_residual_contribution_cons_axisym_adv_diff(), 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::RefineableSphericalNavierStokesEquations::fill_in_generic_residual_contribution_spherical_nst(), oomph::SphericalNavierStokesEquations::fill_in_generic_residual_contribution_spherical_nst(), oomph::SolidFiniteElement::fill_in_jacobian_for_newmark_accel(), oomph::AxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::RefineableAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::RefineableGeneralisedNewtonianAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), and oomph::ImmersedRigidBodyElement::get_residuals_rigid_body_generic().
|
inline |
Get a (const) pointer to the weights.
Definition at line 475 of file timesteppers.h.
|
protected |
Boolean variable to indicate whether the timestepping scheme can be adaptive.
Definition at line 245 of file timesteppers.h.
Referenced by oomph::IMRBase::IMRBase(), and oomph::TR::TR().
|
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 265 of file timesteppers.h.
Referenced by ~TimeStepper().
|
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 251 of file timesteppers.h.
Referenced by oomph::ContinuationStorageScheme::ContinuationStorageScheme(), oomph::IMRBase::IMRBase(), and oomph::ContinuationStorageScheme::undo_make_steady().
|
protected |
Flag: is adaptivity done by taking a separate step using an ExplicitTimeStepper object?
Definition at line 260 of file timesteppers.h.
Referenced by oomph::IMRBase::IMRBase().
|
protected |
Store the time that the predicted values currently stored are at, to compare for paranoid checks.
Definition at line 269 of file timesteppers.h.
|
protected |
The time-index in each Data object where predicted values are stored. -1 if not set.
Definition at line 273 of file timesteppers.h.
Referenced by oomph::IMRBase::IMRBase(), and oomph::IMRBase::temporal_error_in_value().
|
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 256 of file timesteppers.h.
|
protected |
Pointer to discrete time storage scheme.
Definition at line 234 of file timesteppers.h.
Referenced by oomph::PeriodicOrbitTimeDiscretisation::assign_initial_data_values(), oomph::TR::set_error_weights(), oomph::TR::set_predictor_weights(), oomph::IMR::set_weights(), oomph::IMRByBDF::set_weights(), and oomph::TR::set_weights().
|
protected |
String that indicates the type of the timestepper (e.g. "BDF", "Newmark", etc.)
Definition at line 241 of file timesteppers.h.
Referenced by oomph::ContinuationStorageScheme::ContinuationStorageScheme(), oomph::IMRBase::IMRBase(), and oomph::PeriodicOrbitTimeDiscretisation::PeriodicOrbitTimeDiscretisation().
|
protected |
Storage for the weights associated with the timestepper.
Definition at line 237 of file timesteppers.h.
Referenced by oomph::IMRBase::IMRBase(), oomph::ContinuationStorageScheme::modify_storage(), oomph::PeriodicOrbitEquations::set_timestepper_weights(), oomph::IMR::set_weights(), oomph::IMRByBDF::set_weights(), oomph::TR::set_weights(), oomph::IMRBase::shift_time_positions(), and oomph::TR::TR().