26 #ifndef OOMPH_PERIODIC_ORBIT_HANDLER_CLASS_HEADER
27 #define OOMPH_PERIODIC_ORBIT_HANDLER_CLASS_HEADER
31 #include <oomph-lib-config.h>
41 #include "../meshes/one_d_mesh.template.h"
42 #include "../meshes/one_d_mesh.template.cc"
46 class PeriodicOrbitEquations;
48 class PeriodicOrbitAssemblyHandlerBase;
65 Type =
"PeriodicOrbitTimeDiscretisation";
87 "Cannot perform impulsive start for PeriodicOrbitTimeDiscretisation",
88 OOMPH_CURRENT_FUNCTION,
89 OOMPH_EXCEPTION_LOCATION);
97 "Cannot perform impulsive start for PeriodicOrbitTimeDiscretisation",
98 OOMPH_CURRENT_FUNCTION,
99 OOMPH_EXCEPTION_LOCATION);
117 unsigned n_value = data_pt->
nvalue();
120 for (
unsigned t = 0;
t < n_time_value;
t++)
126 for (
unsigned j = 0; j < n_value; j++)
137 "Cannot shift time values for PeriodicOrbitTimeDiscretisation",
138 OOMPH_CURRENT_FUNCTION,
139 OOMPH_EXCEPTION_LOCATION);
146 "Cannot shift time positions for PeriodicOrbitTimeDiscretisation",
147 OOMPH_CURRENT_FUNCTION,
148 OOMPH_EXCEPTION_LOCATION);
265 assembly_handler_pt, elem_pt, residuals, jacobian, 1);
270 std::ostream& outfile,
271 const unsigned& n_plot);
281 const unsigned& flag);
291 unsigned n_time_dof = cast_time_stepper_pt->
ntstorage();
292 for (
unsigned i = 0;
i < n_time_dof;
i++)
294 cast_time_stepper_pt->
Weight(0,
i) = 0.0;
295 cast_time_stepper_pt->
Weight(1,
i) = 0.0;
299 const double inverse_timescale = this->
omega();
301 const unsigned n_node = this->
nnode();
306 for (
unsigned l = 0; l < n_node; l++)
309 cast_time_stepper_pt->
Weight(0, global_eqn) = psi(l);
310 cast_time_stepper_pt->
Weight(1, global_eqn) =
311 dpsidt(l, 0) * inverse_timescale;
321 DShape& dtestdt)
const = 0;
331 DShape& dtestdt)
const = 0;
347 template<
unsigned NNODE_1D>
401 const unsigned n_node = this->
nnode();
405 for (
unsigned n = 0; n < n_node; n++)
428 const unsigned n_node = this->
nnode();
441 for (
unsigned t = 0;
t < n_tstorage;
t++)
447 for (
unsigned n = 0; n < n_node; n++)
449 const double dpsidx_ = dpsidx(n, 0);
450 for (
unsigned t = 0;
t < n_tstorage;
t++)
483 void output(std::ostream& outfile,
const unsigned& n_plot)
499 void output(FILE* file_pt,
const unsigned& n_plot)
535 template<
unsigned NNODE_1D>
536 double SpectralPeriodicOrbitElement<
544 const double J = this->dshape_eulerian(
s, psi, dpsidt);
561 template<
unsigned NNODE_1D>
563 NNODE_1D>::dshape_and_dtest_eulerian_at_knot_orbit(
const unsigned& ipt,
570 const double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidt);
581 template<
unsigned NNODE_1D>
588 template<
class ELEMENT>
592 template<
unsigned NNODE_1D>
598 :
OneDMesh<ELEMENT>(n_element, 1.0),
610 std::ostream& outfile,
611 const unsigned& n_plot)
614 const unsigned n_element = this->
nelement();
615 for (
unsigned e = 0;
e < n_element;
e++)
631 const unsigned n_element = this->
nelement();
632 for (
unsigned e = 0;
e < n_element;
e++)
635 ->fill_in_contribution_to_integrated_residuals(
636 assembly_handler_pt, elem_pt, residuals);
653 const unsigned n_element = this->
nelement();
654 for (
unsigned e = 0;
e < n_element;
e++)
657 ->fill_in_contribution_to_integrated_jacobian(
658 assembly_handler_pt, elem_pt, residuals, jacobian);
689 template<
unsigned NNODE_1D>
725 const unsigned& n_element_in_period,
730 Omega(omega / (2.0 * MathematicalConstants::
Pi))
740 n_element_in_period);
766 for (
unsigned e = 0;
e < n_element_in_period;
e++)
793 "PeriodicOrbitHandler can't cope with submeshes yet",
794 OOMPH_CURRENT_FUNCTION,
795 OOMPH_EXCEPTION_LOCATION);
800 for (
unsigned n = 0; n < n_node; n++)
805 unsigned n_value = nod_pt->
nvalue();
806 for (
unsigned i = 0;
i < n_value;
i++)
810 const unsigned n_tstorage = nod_pt->
ntstorage();
811 const double value = nod_pt->
value(
i);
812 for (
unsigned t = 1;
t < n_tstorage;
t++)
821 for (
unsigned e = 0;
e < n_element;
e++)
823 unsigned n_internal =
825 for (
unsigned i = 0;
i < n_internal;
i++)
828 Data*
const data_pt =
834 unsigned n_value = data_pt->
nvalue();
835 for (
unsigned j = 0; j < n_value; j++)
839 const unsigned n_tstorage = data_pt->
ntstorage();
840 const double value = data_pt->
value(j);
841 for (
unsigned t = 1;
t < n_tstorage;
t++)
853 <<
" equation numbers\n";
862 unsigned offset =
Ndof *
i;
863 for (
unsigned n = 0; n <
Ndof; n++)
880 unsigned offset =
Ndof *
i;
881 for (
unsigned n = 0; n <
Ndof; n++)
883 problem_pt->
dof(offset + n) =
919 const unsigned& ieqn_local)
922 unsigned raw_ndof = elem_pt->
ndof();
932 unsigned t = ieqn_local / raw_ndof;
934 unsigned raw_ieqn = ieqn_local % raw_ndof;
944 Time_mesh_pt->assemble_residuals(
this, elem_pt, residuals);
953 const unsigned n_elem_dof = this->
ndof(elem_pt);
954 dofs.resize(n_elem_dof);
957 for (
unsigned i = 0;
i < n_elem_dof;
i++)
967 const unsigned n_elem_dof = this->
ndof(elem_pt);
968 dofs.resize(n_elem_dof);
971 for (
unsigned i = 0;
i < n_elem_dof;
i++)
982 const unsigned n_elem_dof = this->
ndof(elem_pt);
985 for (
unsigned i = 0;
i < n_elem_dof;
i++)
999 this, elem_pt, residuals, jacobian);
1025 for (
unsigned e = 0;
e < n_element;
e++)
1038 for (
unsigned n = 0; n < n_node; n++)
1058 unsigned total_n_value = 0;
1062 for (
unsigned i = 0;
i < n_global_data;
i++)
1069 for (
unsigned n = 0; n < n_node; n++)
1074 if (solid_nod_pt != 0)
1082 for (
unsigned e = 0;
e < n_space_element;
e++)
1084 const unsigned n_internal_data =
1086 for (
unsigned i = 0;
i < n_internal_data;
i++)
1100 for (
unsigned t = 0;
t < n_time_node;
t++)
1102 Time_mesh_pt->node_pt(
t)->set_time_stepper(fake_space_time_stepper_pt,
1108 for (
unsigned i = 0;
i < n_global_data;
i++)
1111 const unsigned n_value = glob_data_pt->
nvalue();
1112 for (
unsigned j = 0; j < n_value; j++)
1119 count, 0, glob_data_pt->
value(
t, j));
1126 for (
unsigned n = 0; n < n_node; n++)
1129 const unsigned n_value = nod_pt->
nvalue();
1130 for (
unsigned i = 0;
i < n_value;
i++)
1141 if (solid_nod_pt != 0)
1143 const unsigned n_solid_value =
1145 for (
unsigned i = 0;
i < n_solid_value;
i++)
1158 for (
unsigned e = 0;
e < n_space_element;
e++)
1160 const unsigned n_internal =
1162 for (
unsigned i = 0;
i < n_internal;
i++)
1164 Data*
const internal_dat_pt =
1167 const unsigned n_value = internal_dat_pt->
nvalue();
1168 for (
unsigned j = 0; j < n_value; j++)
1173 count, 0, internal_dat_pt->
value(
t, j));
1204 const unsigned n_time_element =
Time_mesh_pt->nelement();
1206 Time_mesh_pt->spatial_error_estimator_pt()->get_element_errors(
1210 for (
unsigned e = 0;
e < n_time_element;
e++)
1299 "PeriodicOrbitHandler can't cope with submeshes yet",
1300 OOMPH_CURRENT_FUNCTION,
1301 OOMPH_EXCEPTION_LOCATION);
1305 for (
unsigned n = 0; n < n_node; n++)
1311 for (
unsigned e = 0;
e < n_space_element;
e++)
1313 unsigned n_internal =
1315 for (
unsigned i = 0;
i < n_internal;
i++)
1318 Data*
const data_pt =
1332 <<
" equation numbers\n";
1341 unsigned offset =
Ndof *
i;
1342 for (
unsigned n = 0; n <
Ndof; n++)
1360 for (
unsigned i = 0;
i < n_global_data;
i++)
1363 const unsigned n_value = glob_data_pt->
nvalue();
1364 for (
unsigned j = 0; j < n_value; j++)
1378 for (
unsigned n = 0; n < n_node; n++)
1381 const unsigned n_value = nod_pt->
nvalue();
1382 for (
unsigned i = 0;
i < n_value;
i++)
1393 if (solid_nod_pt != 0)
1395 const unsigned n_solid_value =
1397 for (
unsigned i = 0;
i < n_solid_value;
i++)
1410 for (
unsigned e = 0;
e < n_space_element;
e++)
1412 const unsigned n_internal =
1414 for (
unsigned i = 0;
i < n_internal;
i++)
1416 Data*
const internal_dat_pt =
1419 const unsigned n_value = internal_dat_pt->
nvalue();
1420 for (
unsigned j = 0; j < n_value; j++)
1435 for (
unsigned t = 0;
t < n_time_node;
t++)
1441 delete fake_space_time_stepper_pt;
1474 const unsigned n_dof = this->
ndof();
1476 for (
unsigned i = 0;
i < n_dof;
i++)
1478 inner_product(
i,
i) = 1.0;
1483 const unsigned& Nplot,
1484 const double& time = 0.0)
A class that is used to define the functions used to assemble the elemental contributions to the resi...
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
A class that represents a collection of data; each Data object may contain many different individual ...
void set_value(const unsigned &i, const double &value_)
Set the i-th stored data value to specified value. The only reason that we require an explicit set fu...
unsigned nvalue() const
Return number of values stored in data object (incl pinned ones).
double value(const unsigned &i) const
Return i-th stored value. This function is not virtual so that it can be inlined. This means that if ...
unsigned ntstorage() const
Return total number of doubles stored per value to record time history of each value (one for steady ...
void set_time_stepper(TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
Set a new timestepper by resizing the appropriate storage. If already assigned the equation numbering...
bool is_pinned(const unsigned &i) const
Test whether the i-th variable is pinned (1: true; 0: false).
void initialise(const T &val)
Initialize all values in the matrix to val.
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
A general Finite Element class.
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
double nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n. Produces suitably interpolated values for hanging nodes...
virtual void output(std::ostream &outfile)
Output the element data — typically the values at the nodes in a format suitable for post-processing.
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
unsigned nnode() const
Return the number of nodes.
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Compute the geometric shape functions and also first derivatives w.r.t. global coordinates at local c...
A Generalised Element class.
unsigned ndof() const
Return the number of equations/dofs in the element.
unsigned long eqn_number(const unsigned &ieqn_local) const
Return the global equation number corresponding to the ieqn_local-th local equation number.
Data *& internal_data_pt(const unsigned &i)
Return a pointer to i-th internal data object.
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
unsigned ninternal_data() const
Return the number of internal data objects.
TimeStepper *& time_stepper_pt()
Access function for pointer to time stepper: Null if object is not time-dependent.
void build(const OomphCommunicator *const comm_pt, const unsigned &first_row, const unsigned &nrow_local, const unsigned &nrow=0)
Sets the distribution. Takes first_row, nrow_local and nrow as arguments. If nrow is not provided or ...
static Steady< 0 > Default_TimeStepper
Default Steady Timestepper, to be used in default arguments to Mesh constructors.
GeneralisedElement *& element_pt(const unsigned long &e)
Return pointer to element e.
unsigned long nnode() const
Return number of nodes in the mesh.
const Vector< GeneralisedElement * > & element_pt() const
Return reference to the Vector of elements.
Node *& node_pt(const unsigned long &n)
Return pointer to global node n.
Node *& boundary_node_pt(const unsigned &b, const unsigned &n)
Return pointer to node n on boundary b.
unsigned long nelement() const
Return number of elements in the mesh.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
virtual void make_periodic(Node *const &node_pt)
Make the node periodic by copying the values from node_pt. Note that the coordinates will always rema...
double value(const unsigned &i) const
Return i-th value (dofs or pinned) at this node either directly or via hanging node representation....
1D mesh consisting of N one-dimensional elements from the QElement family.
An OomphLibError object which should be thrown when an run-time error is encountered....
=============================================================== Base class to avoid template complica...
virtual void set_dofs_for_element(GeneralisedElement *const elem_pt, Vector< double > const &dofs)=0
virtual void get_dofs_for_element(GeneralisedElement *const elem_pt, Vector< double > &dofs)=0
virtual void get_previous_dofs_for_element(GeneralisedElement *const elem_pt, Vector< double > &dofs)=0
PeriodicOrbitAssemblyHandlerBase()
A class that is used to assemble and solve the augmented system of equations associated with calculat...
void adapt_temporal_mesh()
Adapt the time mesh.
void set_dofs_for_element(GeneralisedElement *const elem_pt, Vector< double > const &dofs)
void orbit_output(std::ostream &outfile, const unsigned &n_plot)
Calculate all desired vectors and matrices provided by the element elem_pt.
unsigned N_element_in_period
Storage for number of elements in the period.
void discrete_times(Vector< double > &t)
Tell me the times at which you want the solution.
void get_previous_dofs_for_element(GeneralisedElement *const elem_pt, Vector< double > &dofs)
void get_jacobian(GeneralisedElement *const &elem_pt, Vector< double > &residuals, DenseMatrix< double > &jacobian)
Calculate the elemental Jacobian matrix "d equation / d variable" for elem_pt.
~PeriodicOrbitAssemblyHandler()
Destructor, destroy the time mesh.
void get_dofs_for_element(GeneralisedElement *const elem_pt, Vector< double > &dofs)
unsigned Ndof
Store number of degrees of freedom in the original problem.
PeriodicOrbitTemporalMesh< SpectralPeriodicOrbitElement< NNODE_1D > > * Time_mesh_pt
Storage for mesh of temporal elements.
Vector< double > Previous_dofs
Storage for the previous solution.
unsigned ndof(GeneralisedElement *const &elem_pt)
Return the number of degrees of freedom in the element elem_pt.
void set_previous_dofs_to_current_dofs()
Update the previous dofs.
Problem * Problem_pt
Pointer to the problem.
unsigned long eqn_number(GeneralisedElement *const &elem_pt, const unsigned &ieqn_local)
Return the global equation number of the local unknown ieqn_local in elem_pt.
unsigned N_tstorage
Storage for the number of unknown time values.
PeriodicOrbitTimeDiscretisation * Time_stepper_pt
Pointer to the timestepper.
Mesh * Basic_time_mesh_pt
Storage for the mesh of temporal elements with a simple mesh pointer.
double Omega
Storage for the frequency of the orbit (scaled by 2pi)
void get_residuals(GeneralisedElement *const &elem_pt, Vector< double > &residuals)
Return the contribution to the residuals of the element elem_pt.
virtual void get_non_external_ddofs_dt(Vector< double > &du_dt)
Interface to get the current value of the time derivative of all (internal and shared) unknowns.
PeriodicOrbitBaseElement()
virtual void spacetime_output(std::ostream &outilfe, const unsigned &Nplot, const double &time=0.0)
virtual void get_non_external_dofs(Vector< double > &u)
Interface to get the current value of all (internal and shared) unknowns.
virtual void get_inner_product_matrix(DenseMatrix< double > &inner_product)
Get the inner product matrix.
Time * Time_pt
Pointer to global time.
void orbit_output(GeneralisedElement *const &elem_pt, std::ostream &outfile, const unsigned &n_plot)
void fill_in_contribution_to_integrated_jacobian(PeriodicOrbitAssemblyHandlerBase *const &assembly_handler_pt, GeneralisedElement *const &elem_pt, Vector< double > &residuals, DenseMatrix< double > &jacobian)
Add the element's contribution to its residual vector and element Jacobian matrix (wrapper)
void set_timestepper_weights(const Shape &psi, const DShape &dpsidt)
Set the timestepper weights.
virtual double dshape_and_dtest_eulerian_at_knot_orbit(const unsigned &ipt, Shape &psi, DShape &dpsidt, Shape &test, DShape &dtestdt) const =0
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of ...
void fill_in_contribution_to_integrated_residuals(PeriodicOrbitAssemblyHandlerBase *const &assembly_handler_pt, GeneralisedElement *const &elem_pt, Vector< double > &residuals)
Add the element's contribution to its residual vector (wrapper)
double *& omega_pt()
Broken assignment operator.
Time *& time_pt()
Retun the pointer to the global time.
virtual double dshape_and_dtest_eulerian_orbit(const Vector< double > &s, Shape &psi, DShape &dpsidt, Shape &test, DShape &dtestdt) const =0
Shape/test functions and derivs w.r.t. to global coords at local coord. s; return Jacobian of mapping...
double time() const
Return the global time, accessed via the time pointer.
void fill_in_generic_residual_contribution_orbit(PeriodicOrbitAssemblyHandlerBase *const &assembly_handler_pt, GeneralisedElement *const &elem_pt, Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
The routine that actually does all the work!
Time *const & time_pt() const
Return the pointer to the global time (const version)
void set_ntstorage(const unsigned &n_tstorage)
Set the total number of time storage values.
double omega()
Return the frequency.
PeriodicOrbitEquations(const PeriodicOrbitEquations &dummy)=delete
Broken copy constructor.
A special temporal mesh class.
PeriodicOrbitTemporalMesh(const unsigned &n_element)
Constructor, create a 1D mesh from 0 to 1 that is periodic.
void assemble_residuals_and_jacobian(PeriodicOrbitAssemblyHandlerBase *const &assembly_handler_pt, GeneralisedElement *const &elem_pt, Vector< double > &residuals, DenseMatrix< double > &jacobian)
void orbit_output(GeneralisedElement *const &elem_pt, std::ostream &outfile, const unsigned &n_plot)
void assemble_residuals(PeriodicOrbitAssemblyHandlerBase *const &assembly_handler_pt, GeneralisedElement *const &elem_pt, Vector< double > &residuals)
Timestepper used to calculate periodic orbits directly. It's not really a "timestepper" per se,...
double(* InitialConditionFctPt)(const double &t)
Typedef for function that returns the (scalar) initial value at a given value of the continuous time ...
void shift_time_positions(Node *const &node_pt)
Broken shifting of time positions.
void set_weights()
Set the weights.
unsigned order() const
Return the actual order of the scheme.
void assign_initial_values_impulsive(Data *const &data_pt)
Broken initialisation the time-history for the Data values corresponding to an impulsive start.
void shift_time_values(Data *const &data_pt)
Broken shifting of time values.
PeriodicOrbitTimeDiscretisation(const unsigned &n_tstorage)
Constructor for the case when we allow adaptive timestepping.
void assign_initial_positions_impulsive(Node *const &node_pt)
Broken initialisation of the positions for the node corresponding to an impulsive start.
void assign_initial_data_values(Data *const &data_pt, Vector< InitialConditionFctPt > initial_value_fct)
Initialise the time-history for the Data values, corresponding to given time history,...
unsigned ndt() const
Number of timestep increments that need to be stored by the scheme.
PeriodicOrbitTimeDiscretisation(const PeriodicOrbitTimeDiscretisation &)=delete
Broken copy constructor.
void operator=(const PeriodicOrbitTimeDiscretisation &)=delete
Broken assignment operator.
unsigned nprev_values() const
Number of previous values available.
////////////////////////////////////////////////////////////////// //////////////////////////////////...
unsigned long assign_eqn_numbers(const bool &assign_local_eqn_numbers=true)
Assign all equation numbers for problem: Deals with global data (= data that isn't attached to any el...
unsigned long ndof() const
Return the number of dofs.
Time *& time_pt()
Return a pointer to the global time object.
OomphCommunicator * communicator_pt()
access function to the oomph-lib communicator
unsigned nglobal_data() const
Return the number of global data values.
double & dof(const unsigned &i)
i-th dof in the problem
Vector< double * > Dof_pt
Vector of pointers to dofs.
Mesh *& mesh_pt()
Return a pointer to the global mesh.
unsigned nsub_mesh() const
Return number of submeshes.
LinearAlgebraDistribution * Dof_distribution_pt
The distribution of the DOFs in this problem. This object is created in the Problem constructor and s...
Data *& global_data_pt(const unsigned &i)
Return a pointer to the the i-th global data object.
void shape(const Vector< double > &s, Shape &psi) const
Calculate the geometric shape functions at local coordinate s.
General QLegendreElement class.
Refineable version of the OneDMesh.
A class that is used to template the refineable Q spectral elements by dimension. It's really nothing...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
A Class for nodes that deform elastically (i.e. position is an unknown in the problem)....
Data *const & variable_position_pt() const
Pointer to variable_position data (const version)
QPoissonElement elements are linear/quadrilateral/brick-shaped Poisson elements with isoparametric in...
double dshape_and_dtest_eulerian_at_knot_orbit(const unsigned &ipt, Shape &psi, DShape &dpsidt, Shape &test, DShape &dtestdt) const
Shape, test functions & derivs. w.r.t. to global coords. at integration point ipt....
double dshape_and_dtest_eulerian_orbit(const Vector< double > &s, Shape &psi, DShape &dpsidt, Shape &test, DShape &dtestdt) const
Shape, test functions & derivs. w.r.t. to global coords. Return Jacobian.
void get_interpolated_values(const Vector< double > &s, Vector< double > &value)
Return the dummy values.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
void output(FILE *file_pt)
C-style output function: x,y,u or x,y,z,u.
unsigned nrecovery_order()
Order of recovery shape functions.
unsigned nvertex_node() const
Return the number of vertex nodes in this element. Broken virtual function in "pure" finite elements.
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &value)
Return the temporal dummy values.
unsigned num_Z2_flux_terms()
Number of flux terms for Z2 error estimation This will be used to represent all spatial values,...
unsigned ncont_interpolated_values() const
Number of continuously interpolated values (1)
unsigned required_nvalue(const unsigned &n) const
Broken assignment operator.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function: x,y,u or x,y,z,u at n_plot^DIM plot points.
void output(std::ostream &outfile)
Function to return the number of values.
SpectralPeriodicOrbitElement(const SpectralPeriodicOrbitElement< NNODE_1D > &dummy)=delete
Broken copy constructor.
SpectralPeriodicOrbitElement()
Constructor: Call constructors for QElement and Poisson equations.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function: x,y,u or x,y,z,u at n_plot^DIM plot points.
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get the fluxes for the recovert.
////////////////////////////////////////////////////////////////////// //////////////////////////////...
DenseMatrix< double > Weight
Storage for the weights associated with the timestepper.
Time * Time_pt
Pointer to discrete time storage scheme.
unsigned ntstorage() const
Return the number of doubles required to represent history (one for steady)
std::string Type
String that indicates the type of the timestepper (e.g. "BDF", "Newmark", etc.)
double & time()
Return current value of continous time.
Class to keep track of discrete/continous time. It is essential to have a single Time object when usi...
double & time()
Return the current value of the continuous time.
A slight extension to the standard template vector class so that we can include "graceful" array rang...
void initialise(const _Tp &__value)
Iterate over all values and set to the desired value.
Z2-error-estimator: Elements that can be used with Z2 error estimation should be derived from the bas...
const double Pi
50 digits from maple
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
OomphInfo oomph_info
Single (global) instantiation of the OomphInfo object – this is used throughout the library as a "rep...