![]() |
|
#include <nodes.h>
Public Member Functions | |
virtual void | clear_copied_pointers () |
Data () | |
Default constructor. More... | |
Data (const unsigned &initial_n_value) | |
Data (TimeStepper *const &time_stepper_pt, const unsigned &initial_n_value, const bool &allocate_storage=true) | |
Data (const Data &data)=delete | |
Broken copy constructor. More... | |
void | operator= (const Data &)=delete |
Broken assignment operator. More... | |
virtual | ~Data () |
Destructor, deallocates memory assigned for data. More... | |
void | set_time_stepper (TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data) |
TimeStepper *& | time_stepper_pt () |
Return the pointer to the timestepper. More... | |
TimeStepper *const & | time_stepper_pt () const |
Return the pointer to the timestepper (const version). More... | |
virtual bool | is_a_copy () const |
virtual bool | is_a_copy (const unsigned &i) const |
void | set_value (const unsigned &i, const double &value_) |
void | set_value (const unsigned &t, const unsigned &i, const double &value_) |
double | value (const unsigned &i) const |
double | value (const unsigned &t, const unsigned &i) const |
void | value (Vector< double > &values) const |
Compute Vector of values for the Data value. More... | |
void | value (const unsigned &t, Vector< double > &values) const |
double * | value_pt (const unsigned &i) const |
double * | value_pt (const unsigned &t, const unsigned &i) const |
bool | does_pointer_correspond_to_value (double *const ¶meter_pt) |
Check whether the pointer parameter_pt addresses internal data values. More... | |
void | copy (Data *orig_data_pt) |
Copy Data values from specified Data object. More... | |
void | dump (std::ostream &dump_file) const |
Dump the data object to a file. More... | |
void | read (std::ifstream &restart_file) |
Read data object from a file. More... | |
long * | eqn_number_pt (const unsigned &i) |
Return the pointer to the equation number of the i-th stored variable. More... | |
long & | eqn_number (const unsigned &i) |
Return the equation number of the i-th stored variable. More... | |
long | eqn_number (const unsigned &i) const |
Return the equation number of the i-th stored variable. More... | |
void | pin (const unsigned &i) |
Pin the i-th stored variable. More... | |
void | unpin (const unsigned &i) |
Unpin the i-th stored variable. More... | |
void | pin_all () |
Pin all the stored variables. More... | |
void | unpin_all () |
Unpin all the stored variables. More... | |
bool | is_pinned (const unsigned &i) const |
Test whether the i-th variable is pinned (1: true; 0: false). More... | |
bool | is_segregated_solve_pinned (const unsigned &i) |
void | constrain (const unsigned &i) |
void | unconstrain (const unsigned &i) |
void | constrain_all () |
Constrain all the stored variables when the data is made hanging. More... | |
void | unconstrain_all () |
Unconstrain all the stored variables when the data is made nonhanging. More... | |
bool | is_constrained (const unsigned &i) |
unsigned | self_test () |
unsigned | nvalue () const |
Return number of values stored in data object (incl pinned ones). More... | |
unsigned | ntstorage () const |
virtual void | assign_eqn_numbers (unsigned long &global_ndof, Vector< double * > &dof_pt) |
virtual void | describe_dofs (std::ostream &out, const std::string ¤t_string) const |
virtual void | resize (const unsigned &n_value) |
Change (increase) the number of values that may be stored. More... | |
virtual void | add_value_pt_to_map (std::map< unsigned, double * > &map_of_value_pt) |
Static Public Attributes | |
static long | Is_pinned = -1 |
Static "Magic number" to indicate pinned values. More... | |
static long | Is_segregated_solve_pinned = -3 |
static long | Is_unclassified = -10 |
static long | Is_constrained = -2 |
Protected Member Functions | |
virtual void | reset_copied_pointers () |
Protected Attributes | |
Data ** | Copy_of_data_pt |
unsigned | Ncopies |
Static Protected Attributes | |
static TimeStepper * | Default_static_time_stepper_pt = new Steady<0>() |
Default (static) timestepper used in steady problems. More... | |
Private Member Functions | |
void | range_check (const unsigned &t, const unsigned &i) const |
void | delete_value_storage () |
void | add_copy (Data *const &data_pt) |
void | remove_copy (Data *const &data_pt) |
Private Attributes | |
double ** | Value |
long * | Eqn_number |
TimeStepper * | Time_stepper_pt |
unsigned | Nvalue |
Number of values stored in the data object. More... | |
Friends | |
class | HijackedData |
class | CopiedData |
class | BoundaryNodeBase |
template<class NODE_TYPE > | |
class | BoundaryNode |
class | SolidNode |
std::ostream & | operator<< (std::ostream &out, const Data &d) |
A class that represents a collection of data; each Data object may contain many different individual values, as would be natural in non-scalar problems. Data provides storage for auxiliary ‘history’ values that are used by TimeStepper objects to calculate the time derivatives of the stored data and also stores a pointer to the appropriate TimeStepper object. In addition, an associated (global) equation number is stored for each value.
The Data class permits copies of the stored data values and equation numbers into another Data object using the copy() function. Shallow (pointer based) copies of the values can be obtained by using specific derived classes. In such cases pointers to the objects that contain the pointer-based copies should be stored in the original Data class (in the array Copy_of_data_pt) so that resize and destruction operations can be performed safely.
oomph::Data::Data | ( | ) |
Default constructor.
Default: Just set pointer to (steady) timestepper. No storage for values is allocated.
Referenced by oomph::SolidNode::SolidNode().
oomph::Data::Data | ( | const unsigned & | initial_n_value | ) |
Default constructor for steady problems: assign memory for initial_n_value values.
Default constructor for steady problems. Memory is assigned for a given number of values, which are assumed to be free (not pinned)
References Eqn_number, i, Is_unclassified, and Value.
oomph::Data::Data | ( | TimeStepper *const & | time_stepper_pt_, |
const unsigned & | initial_n_value, | ||
const bool & | allocate_storage = true |
||
) |
Constructor for unsteady problems: assign memory for initial_n_value values and any memory required by the Timestepper for the storage of history values.
Constructor for unsteady problems. Memory is assigned for a given number of values; and the additional storage required by the Timestepper. The values are assumed to be free (not pinned).
References Eqn_number, i, Is_unclassified, ntstorage(), plotPSD::t, and Value.
|
delete |
Broken copy constructor.
|
virtual |
Destructor, deallocates memory assigned for data.
Virtual destructor, deallocates memory assigned for data.
References clear_copied_pointers(), Copy_of_data_pt, delete_value_storage(), i, and Ncopies.
|
private |
Add the pointer data_pt to the array Copy_of_data_pt. This should be used whenever copies are made of the data.
Add the pointer data_pt to the internal storage used to keep track of copies of the Data object.
References Copy_of_data_pt, i, and Ncopies.
Referenced by oomph::CopiedData::CopiedData(), oomph::HijackedData::HijackedData(), oomph::BoundaryNodeBase::make_node_periodic(), and oomph::BoundaryNodeBase::make_nodes_periodic().
Add pointers to all unpinned and unconstrained data to a map indexed by (global) equation number
Reimplemented in oomph::SolidNode.
References eqn_number(), i, nvalue(), and value_pt().
Referenced by oomph::GeneralisedElement::add_internal_value_pt_to_map(), and oomph::SolidNode::add_value_pt_to_map().
|
virtual |
Assign global equation numbers; increment global number of unknowns, global_ndof; and add any new dofs to the dof_pt.
Assign (global) equation number. This function does NOT initialise the value because if we're using things like node position as variables in the problem they will have been set before the call to assign equation numbers and setting it to zero will wipe it out :(.
Pass:
Reimplemented in oomph::SolidNode, oomph::Node, oomph::CopiedData, and oomph::HijackedData.
References Eqn_number, eqn_number(), i, is_constrained(), Is_pinned, is_pinned(), is_segregated_solve_pinned(), Nvalue, and value_pt().
Referenced by oomph::Node::assign_eqn_numbers(), oomph::BoundaryNode< NODE_TYPE >::assign_eqn_numbers(), oomph::SolidNode::assign_eqn_numbers(), and oomph::GeneralisedElement::assign_internal_eqn_numbers().
|
virtual |
Helper function that should be overloaded derived classes that contain copies of data. The function must unset (NULL out) the internal pointers to the copied data. This is used when destructing data to ensure that all pointers remain valid. The default implementation throws an error because Data cannot be a copy.
Helper function that should be overloaded classes that contain copies of data. The function must unset (NULL out) the internal pointers to the copied data. This is used when destructing data to ensure that all pointers remain valid.
Reimplemented in oomph::CopiedData, and oomph::HijackedData.
References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.
Referenced by ~Data().
|
inline |
Constrain the i-th stored variable when making hanging data If the data is already pinned leave it along, otherwise mark as constrained (hanging)
References eqn_number(), i, Is_constrained, and Is_pinned.
Referenced by oomph::Node::assign_eqn_numbers(), constrain_all(), oomph::Node::resize(), and oomph::Node::set_hanging_pt().
|
inline |
Constrain all the stored variables when the data is made hanging.
References constrain(), i, and Nvalue.
Referenced by oomph::SolidNode::constrain_positions().
void oomph::Data::copy | ( | Data * | orig_data_pt | ) |
Copy Data values from specified Data object.
References j, ntstorage(), nvalue(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, set_value(), plotPSD::t, and value().
Referenced by oomph::Node::copy(), oomph::Problem::copy(), and oomph::SolidNode::copy().
|
private |
Delete all storage allocated by the Data object for values and equation numbers.
Delete all the storage allocated by the Data object and set its pointers to NULL
References Eqn_number, and Value.
Referenced by oomph::BoundaryNodeBase::make_node_periodic(), oomph::BoundaryNodeBase::make_nodes_periodic(), and ~Data().
|
virtual |
Function to describe the dofs of the Node. The ostream specifies the output stream to which the description is written; the string stores the currently assembled output that is ultimately written to the output stream by Data::describe_dofs(...); it is typically built up incrementally as we descend through the call hierarchy of this function when called from Problem::describe_dofs(...)
Function to describe the dofs of the Data. The ostream specifies the output stream to which the description is written; the string stores the currently assembled output that is ultimately written to the output stream by Data::describe_dofs(...); it is typically built up incrementally as we descend through the call hierarchy of this function when called from Problem::describe_dofs(...)
Reimplemented in oomph::SolidNode.
References Eqn_number, eqn_number(), i, Nvalue, and out().
Referenced by oomph::GeneralisedElement::describe_dofs(), oomph::SolidNode::describe_dofs(), oomph::ElementWithExternalElement::describe_local_dofs(), oomph::ElementWithMovingNodes::describe_local_dofs(), oomph::GeneralisedElement::describe_local_dofs(), oomph::FSIWallElement::describe_local_dofs(), oomph::SpectralElement::describe_local_dofs(), and oomph::FiniteElement::describe_nodal_local_dofs().
Check whether the pointer parameter_pt addresses internal data values.
If pointer parameter_pt addresses internal data values then return return true, otherwise return false
References i, ntstorage(), nvalue(), and Value.
Referenced by oomph::Mesh::does_pointer_correspond_to_mesh_data(), and oomph::SolidNode::does_pointer_correspond_to_position_data().
void oomph::Data::dump | ( | std::ostream & | dump_file | ) | const |
Dump the data object to a file.
Dump data object to a file.
References j, ntstorage(), nvalue(), plotPSD::t, and value().
Referenced by oomph::SpineMesh::dump(), oomph::PerturbedSpineMesh::dump(), oomph::Mesh::dump(), and oomph::Node::dump().
|
inline |
Return the equation number of the i-th stored variable.
References Eqn_number, i, and range_check().
Referenced by oomph::NavierStokesImpedanceTractionElement< BULK_NAVIER_STOKES_ELEMENT, WOMERSLEY_ELEMENT, DIM >::add_element_contribution_to_aux_integral(), add_value_pt_to_map(), oomph::ElementWithMovingNodes::assign_all_generic_local_eqn_numbers(), oomph::SpectralElement::assign_all_generic_local_eqn_numbers(), assign_eqn_numbers(), oomph::ElementWithExternalElement::assign_external_interaction_data_local_eqn_numbers(), oomph::RefineableElement::assign_hanging_local_eqn_numbers(), oomph::GeneralisedElement::assign_internal_and_external_local_eqn_numbers(), BrethertonElement< ELEMENT >::assign_local_eqn_numbers(), oomph::GeneralisedElement::assign_local_eqn_numbers(), oomph::FiniteElement::assign_nodal_local_eqn_numbers(), oomph::MyProblem::check_not_segregated(), oomph::BoundaryNode< NODE_TYPE >::clear_copied_pointers(), constrain(), oomph::AddedMainNumberingLookup::construct_added_to_main_mapping(), describe_dofs(), oomph::AdvectionDiffusionEquations< DIM >::dinterpolated_u_adv_diff_ddata(), oomph::RefineableAdvectionDiffusionEquations< DIM >::dinterpolated_u_adv_diff_ddata(), oomph::RefineableAxisymAdvectionDiffusionEquations::dinterpolated_u_adv_diff_ddata(), oomph::RefineableSphericalAdvectionDiffusionEquations::dinterpolated_u_adv_diff_ddata(), oomph::SphericalAdvectionDiffusionEquations::dinterpolated_u_adv_diff_ddata(), oomph::SteadyAxisymAdvectionDiffusionEquations::dinterpolated_u_adv_diff_ddata(), oomph::AxisymAdvectionDiffusionEquations::dinterpolated_u_axi_adv_diff_ddata(), oomph::AxisymmetricNavierStokesEquations::dinterpolated_u_axi_nst_ddata(), oomph::RefineableAxisymmetricNavierStokesEquations::dinterpolated_u_axi_nst_ddata(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::dinterpolated_u_axi_nst_ddata(), oomph::RefineableGeneralisedNewtonianAxisymmetricNavierStokesEquations::dinterpolated_u_axi_nst_ddata(), oomph::GeneralisedNewtonianNavierStokesEquations< DIM >::dinterpolated_u_nst_ddata(), oomph::RefineableGeneralisedNewtonianNavierStokesEquations< DIM >::dinterpolated_u_nst_ddata(), oomph::NavierStokesEquations< DIM >::dinterpolated_u_nst_ddata(), oomph::RefineableNavierStokesEquations< DIM >::dinterpolated_u_nst_ddata(), oomph::SpaceTimeNavierStokesEquations< DIM >::dinterpolated_u_nst_ddata(), oomph::RefineableSpaceTimeNavierStokesEquations< DIM >::dinterpolated_u_nst_ddata(), oomph::SpaceTimeNavierStokesMixedOrderEquations< DIM >::dinterpolated_u_nst_ddata(), oomph::RefineableSpaceTimeNavierStokesMixedOrderEquations< DIM >::dinterpolated_u_nst_ddata(), oomph::NavierStokesImpedanceTractionElement< BULK_NAVIER_STOKES_ELEMENT, WOMERSLEY_ELEMENT, DIM >::fill_in_generic_residual_contribution_fluid_traction(), oomph::Problem::get_dofs(), pin(), oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::pin_all_non_pressure_dofs(), oomph::GeneralisedNewtonianNavierStokesEquations< DIM >::pin_all_non_pressure_dofs(), oomph::NavierStokesEquations< DIM >::pin_all_non_pressure_dofs(), oomph::SpaceTimeNavierStokesEquations< DIM >::pin_all_non_pressure_dofs(), oomph::SpaceTimeNavierStokesMixedOrderEquations< DIM >::pin_all_non_pressure_dofs(), oomph::MGSolver< DIM >::plot(), oomph::UnsteadyHeatFluxPseudoMeltElement< ELEMENT >::plot_residual_landscape(), oomph::SolidNode::position_eqn_number(), oomph::post_midpoint_update(), print_connectivity_matrix(), BrethertonElement< ELEMENT >::reassign_inflow(), oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::reset_pin_status(), oomph::NavierStokesSchurComplementPreconditioner::reset_pin_status(), oomph::MyProblem::segregated_pin_indices(), oomph::Problem::set_dofs(), oomph::SolidICProblem::set_newmark_initial_condition_consistently(), oomph::MGSolver< DIM >::set_self_test_vector(), oomph::MGSolver< DIM >::setup_interpolation_matrices(), oomph::HelmholtzMGPreconditioner< DIM >::setup_interpolation_matrices(), oomph::MGSolver< DIM >::setup_interpolation_matrices_unstructured(), oomph::HelmholtzMGPreconditioner< DIM >::setup_interpolation_matrices_unstructured(), unconstrain(), oomph::MyProblem::undo_segregated_pinning(), and unpin().
|
inline |
Return the equation number of the i-th stored variable.
References Eqn_number, i, and range_check().
|
inline |
Return the pointer to the equation number of the i-th stored variable.
References Eqn_number, i, and range_check().
Referenced by oomph::Hijacked< ELEMENT >::hijack_nodal_position_value(), oomph::Hijacked< ELEMENT >::hijack_nodal_spine_value(), and oomph::VolumeConstraintBoundingElement::set_volume_constraint_element().
|
inlinevirtual |
Return a boolean to indicate whether the Data objact contains any copied values. A base Data object can never be a copy so the default implementation always returns false.
Reimplemented in oomph::CopiedData, and oomph::HijackedData.
Referenced by oomph::Newmark< 2 >::assign_initial_data_values(), oomph::SelfStartingBDF2::assign_initial_values_impulsive(), oomph::Steady< NSTEPS >::assign_initial_values_impulsive(), oomph::Newmark< NSTEPS >::assign_initial_values_impulsive(), oomph::BDF< NSTEPS >::assign_initial_values_impulsive(), oomph::BDF< NSTEPS >::calculate_predicted_values(), oomph::TR::calculate_predicted_values(), oomph::SelfStartingBDF2::calculate_predicted_values_bdf2(), oomph::Mesh::check_for_repeated_nodes(), oomph::FourierDecomposedHelmholtzDtNBoundaryElement< ELEMENT >::complete_setup_of_dependencies(), oomph::HelmholtzDtNBoundaryElement< ELEMENT >::complete_setup_of_dependencies(), oomph::CopiedData::CopiedData(), oomph::HijackedData::HijackedData(), oomph::BoundaryNodeBase::make_node_periodic(), oomph::BoundaryNodeBase::make_nodes_periodic(), oomph::SolidNode::position_is_a_copy(), oomph::post_midpoint_update(), oomph::ContinuationStorageScheme::set_consistent_pinned_values(), set_time_stepper(), oomph::SelfStartingBDF2::shift_time_values(), oomph::IMRBase::shift_time_values(), oomph::Steady< NSTEPS >::shift_time_values(), oomph::NewmarkBDF< NSTEPS >::shift_time_values(), oomph::BDF< NSTEPS >::shift_time_values(), and oomph::TR::shift_time_values().
Return flag to indicate whether the i-th value is a copy. A base Data object can never be a copy so the default implementation always returns false.
Reimplemented in oomph::CopiedData, and oomph::HijackedData.
Test whether the i-th variable is constrained (1: true; 0: false).
References Eqn_number, i, and Is_constrained.
Referenced by SurfactantProblem< ELEMENT, INTERFACE_ELEMENT >::actions_after_adapt(), assign_eqn_numbers(), oomph::Node::assign_eqn_numbers(), FlowAroundCylinderProblem< ELEMENT >::pin_all_base_flow_dofs(), oomph::LinearisedQCrouzeixRaviartElement::pin_real_or_imag(), FlowAroundCylinderProblem< ELEMENT >::unpin_all_base_flow_dofs(), and oomph::LinearisedQCrouzeixRaviartElement::unpin_real_or_imag().
Test whether the i-th variable is pinned (1: true; 0: false).
References Eqn_number, i, and Is_pinned.
Referenced by oomph::TreeBasedRefineableMeshBase::adapt_mesh(), oomph::Problem::arc_length_step_solve(), assign_eqn_numbers(), oomph::Node::assign_eqn_numbers(), oomph::SolidICProblem::backup_original_state(), RefineableDrivenCavityProblem< ELEMENT >::doc_solution(), GlobalParameters::find_node_on_centerline(), TwoDDGProblem< ELEMENT >::limit(), main(), oomph::TreeBasedRefineableMeshBase::p_adapt_mesh(), oomph::SolidNode::position_is_pinned(), print_connectivity_matrix(), oomph::MyProblem::segregated_pin_indices(), oomph::ContinuationStorageScheme::set_consistent_pinned_values(), oomph::Problem::set_pinned_values_to_zero(), and oomph::WomersleyMesh< WOMERSLEY_ELEMENT >::WomersleyMesh().
Test whether the i-th variable is temporaily pinned for a segregated solve.
References Eqn_number, i, and Is_segregated_solve_pinned.
Referenced by assign_eqn_numbers().
unsigned oomph::Data::ntstorage | ( | ) | const |
Return total number of doubles stored per value to record time history of each value (one for steady problems).
Return the total number of doubles stored per value to record the time history of ecah value. The information is read from the time stepper
References oomph::TimeStepper::ntstorage(), and Time_stepper_pt.
Referenced by oomph::TreeBasedRefineableMeshBase::adapt_mesh(), oomph::BackupMeshForProjection< GEOMETRIC_ELEMENT >::BackupMeshForProjection(), copy(), oomph::BackupMeshForProjection< GEOMETRIC_ELEMENT >::copy_onto_original_mesh(), Data(), oomph::Mesh::delete_all_external_storage(), does_pointer_correspond_to_value(), dump(), oomph::SpectralPeriodicOrbitElement< NNODE_1D >::get_Z2_flux(), oomph::ProjectableAdvectionDiffusionReactionElement< ADR_ELEMENT >::nhistory_values_for_projection(), oomph::ProjectableAxisymLinearElasticityElement< AXISYM_LINEAR_ELAST_ELEMENT >::nhistory_values_for_projection(), oomph::ProjectableAxisymmetricTaylorHoodElement< TAYLOR_HOOD_ELEMENT >::nhistory_values_for_projection(), oomph::ProjectableAxisymmetricCrouzeixRaviartElement< CROUZEIX_RAVIART_ELEMENT >::nhistory_values_for_projection(), oomph::ProjectableAxisymmetricPoroelasticityElement< AXISYMMETRIC_POROELASTICITY_ELEMENT >::nhistory_values_for_projection(), oomph::ProjectableDarcyElement< DARCY_ELEMENT >::nhistory_values_for_projection(), oomph::ProjectableDisplacementBasedFoepplvonKarmanElement< FVK_ELEMENT >::nhistory_values_for_projection(), oomph::ProjectableFoepplvonKarmanElement< FVK_ELEMENT >::nhistory_values_for_projection(), oomph::ProjectableFourierDecomposedHelmholtzElement< FOURIER_DECOMPOSED_HELMHOLTZ_ELEMENT >::nhistory_values_for_projection(), oomph::GeneralisedNewtonianProjectableAxisymmetricTaylorHoodElement< TAYLOR_HOOD_ELEMENT >::nhistory_values_for_projection(), oomph::GeneralisedNewtonianProjectableAxisymmetricCrouzeixRaviartElement< CROUZEIX_RAVIART_ELEMENT >::nhistory_values_for_projection(), oomph::ProjectableGeneralisedNewtonianTaylorHoodElement< TAYLOR_HOOD_ELEMENT >::nhistory_values_for_projection(), oomph::ProjectableGeneralisedNewtonianCrouzeixRaviartElement< CROUZEIX_RAVIART_ELEMENT >::nhistory_values_for_projection(), oomph::GenericLagrangeInterpolatedProjectableElement< ELEMENT >::nhistory_values_for_projection(), oomph::ProjectableHelmholtzElement< HELMHOLTZ_ELEMENT >::nhistory_values_for_projection(), oomph::ProjectableLinearElasticityElement< LINEAR_ELAST_ELEMENT >::nhistory_values_for_projection(), oomph::ProjectablePMLFourierDecomposedHelmholtzElement< FOURIER_DECOMPOSED_HELMHOLTZ_ELEMENT >::nhistory_values_for_projection(), oomph::ProjectablePMLHelmholtzElement< HELMHOLTZ_ELEMENT >::nhistory_values_for_projection(), oomph::ProjectablePMLTimeHarmonicLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::nhistory_values_for_projection(), oomph::ProjectablePoissonElement< POISSON_ELEMENT >::nhistory_values_for_projection(), oomph::ProjectableUnsteadyHeatSpaceTimeElement< UNSTEADY_HEAT_ELEMENT >::nhistory_values_for_projection(), oomph::ProjectableUnsteadyHeatMixedOrderSpaceTimeElement< UNSTEADY_HEAT_ELEMENT >::nhistory_values_for_projection(), oomph::ProjectableTimeHarmonicFourierDecomposedLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::nhistory_values_for_projection(), oomph::ProjectableTimeHarmonicLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::nhistory_values_for_projection(), oomph::ProjectableUnsteadyHeatElement< UNSTEADY_HEAT_ELEMENT >::nhistory_values_for_projection(), oomph::SpectralPeriodicOrbitElement< NNODE_1D >::num_Z2_flux_terms(), oomph::TreeBasedRefineableMeshBase::p_adapt_mesh(), range_check(), read(), resize(), and set_time_stepper().
|
inline |
Return number of values stored in data object (incl pinned ones).
References Nvalue.
Referenced by oomph::TreeBasedRefineableMeshBase::adapt_mesh(), oomph::PeriodicOrbitAssemblyHandler< NNODE_1D >::adapt_temporal_mesh(), add_value_pt_to_map(), MovingBlockProblem< ELEMENT >::apply_boundary_conditions(), TwoDDGProblem< ELEMENT >::apply_boundary_conditions(), oomph::ElementWithMovingNodes::assign_all_generic_local_eqn_numbers(), oomph::SpectralElement::assign_all_generic_local_eqn_numbers(), oomph::Node::assign_eqn_numbers(), oomph::ElementWithExternalElement::assign_external_interaction_data_local_eqn_numbers(), oomph::RefineableElement::assign_hanging_local_eqn_numbers(), oomph::SelfStartingBDF2::assign_initial_data_values(), oomph::PeriodicOrbitTimeDiscretisation::assign_initial_data_values(), oomph::Steady< NSTEPS >::assign_initial_data_values(), oomph::BDF< NSTEPS >::assign_initial_data_values(), oomph::Newmark< NSTEPS >::assign_initial_data_values(), oomph::Newmark< NSTEPS >::assign_initial_data_values_stage1(), oomph::Newmark< NSTEPS >::assign_initial_data_values_stage2(), oomph::Newmark< 2 >::assign_initial_positions_impulsive(), oomph::SelfStartingBDF2::assign_initial_values_impulsive(), oomph::Steady< NSTEPS >::assign_initial_values_impulsive(), oomph::Newmark< NSTEPS >::assign_initial_values_impulsive(), oomph::BDF< NSTEPS >::assign_initial_values_impulsive(), oomph::GeneralisedElement::assign_internal_and_external_local_eqn_numbers(), BrethertonElement< ELEMENT >::assign_local_eqn_numbers(), oomph::GeneralisedElement::assign_local_eqn_numbers(), oomph::FiniteElement::assign_nodal_local_eqn_numbers(), oomph::SolidICProblem::backup_original_state(), oomph::BackupMeshForProjection< GEOMETRIC_ELEMENT >::BackupMeshForProjection(), oomph::RefineableQElement< 3 >::build(), oomph::RefineableQElement< 1 >::build(), oomph::RefineableQElement< 2 >::build(), oomph::BDF< NSTEPS >::calculate_predicted_values(), oomph::TR::calculate_predicted_values(), oomph::SelfStartingBDF2::calculate_predicted_values_bdf2(), oomph::MyProblem::check_not_segregated(), oomph::Circle::Circle(), AxiPoroProblem< ELEMENT, TIMESTEPPER >::complete_problem_setup(), AxisymmetricVibratingShellProblem< ELEMENT >::complete_problem_setup(), VibratingShellProblem< ELEMENT >::complete_problem_setup(), UnsteadyHeatMeltProblem< ELEMENT >::complete_problem_setup(), oomph::NodeNodeMortaringElement::construct_lagrange_multipliers(), copy(), oomph::BackupMeshForProjection< GEOMETRIC_ELEMENT >::copy_onto_original_mesh(), ContactProblem< ELEMENT >::create_displ_imposition_elements(), UnstructuredTorusProblem< ELEMENT >::create_lagrange_multiplier_elements(), FSICollapsibleChannelProblem< ELEMENT >::create_lagrange_multiplier_elements(), PseudoElasticCollapsibleChannelProblem< FLUID_ELEMENT, SOLID_ELEMENT >::create_lagrange_multiplier_elements(), FSIChannelWithLeafletProblem< ELEMENT >::create_lagrange_multiplier_elements(), UnstructuredFSIProblem< FLUID_ELEMENT, SOLID_ELEMENT >::create_lagrange_multiplier_elements(), PrescribedBoundaryDisplacementProblem< ELEMENT >::create_lagrange_multiplier_elements(), UnstructuredImmersedEllipseProblem< ELEMENT >::create_lagrange_multiplier_elements(), UnstructuredFluidProblem< ELEMENT >::create_lagrange_multiplier_elements(), CollapsibleChannelProblem< ELEMENT >::create_lagrange_multiplier_elements(), TiltedCavityProblem< ELEMENT >::create_parall_outflow_lagrange_elements(), oomph::ProjectableAxisymmetricPoroelasticityElement< AXISYMMETRIC_POROELASTICITY_ELEMENT >::data_values_of_field(), oomph::ProjectableDarcyElement< DARCY_ELEMENT >::data_values_of_field(), oomph::Mesh::delete_all_external_storage(), oomph::DisplacementControlElement::DisplacementControlElement(), oomph::FSI_functions::doc_fsi(), RefineableDrivenCavityProblem< ELEMENT >::doc_solution(), PRefineableDrivenCavityProblem< ELEMENT >::doc_solution(), does_pointer_correspond_to_value(), dump(), oomph::Ellipse::Ellipse(), oomph::GeneralisedElement::external_local_eqn(), oomph::ConstraintElement::fill_in_contribution_to_residuals(), oomph::GeneralisedElement::fill_in_jacobian_from_external_by_fd(), oomph::GeneralisedElement::fill_in_jacobian_from_internal_by_fd(), oomph::FiniteElement::fill_in_jacobian_from_nodal_by_fd(), oomph::RefineableElement::fill_in_jacobian_from_nodal_by_fd(), oomph::GeneralCircle::GeneralCircle(), oomph::AxisymmetricQCrouzeixRaviartElement::get_dof_numbers_for_unknowns(), oomph::AxisymmetricQTaylorHoodElement::get_dof_numbers_for_unknowns(), oomph::AxisymmetricTCrouzeixRaviartElement::get_dof_numbers_for_unknowns(), oomph::BiharmonicEquations< DIM >::get_dof_numbers_for_unknowns(), oomph::GeneralisedNewtonianAxisymmetricQCrouzeixRaviartElement::get_dof_numbers_for_unknowns(), oomph::GeneralisedNewtonianAxisymmetricQTaylorHoodElement::get_dof_numbers_for_unknowns(), oomph::GeneralisedNewtonianAxisymmetricTCrouzeixRaviartElement::get_dof_numbers_for_unknowns(), oomph::RefineableAdvectionDiffusionBoussinesqElement< AD_ELEMENT, NST_ELEMENT >::get_dof_numbers_for_unknowns(), oomph::ClampedHermiteShellBoundaryConditionElement::get_dof_numbers_for_unknowns(), oomph::BlockPrecQTaylorHoodSpaceTimeElement::get_dof_numbers_for_unknowns(), oomph::BlockPrecQTaylorHoodMixedOrderSpaceTimeElement::get_dof_numbers_for_unknowns(), oomph::QSphericalCrouzeixRaviartElement::get_dof_numbers_for_unknowns(), oomph::QSphericalTaylorHoodElement::get_dof_numbers_for_unknowns(), oomph::Problem::get_dofs(), oomph::MyProblem::get_error_norm(), oomph::ProjectablePMLTimeHarmonicLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::get_field(), oomph::ProjectableTimeHarmonicFourierDecomposedLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::get_field(), oomph::ProjectableTimeHarmonicLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::get_field(), BrethertonElement< ELEMENT >::get_jacobian(), oomph::PolarCrouzeixRaviartElement::get_load_data(), oomph::Node::hang_code(), oomph::Node::hanging_pt(), oomph::FiniteElement::identify_field_data_for_interactions(), oomph::RefineableElement::identify_field_data_for_interactions(), oomph::AxisymmetricTCrouzeixRaviartElement::identify_pressure_data(), oomph::GeneralisedNewtonianAxisymmetricTCrouzeixRaviartElement::identify_pressure_data(), oomph::QSphericalCrouzeixRaviartElement::identify_pressure_data(), oomph::GeneralisedElement::internal_local_eqn(), oomph::BiharmonicEquations< DIM >::interpolated_dudx(), oomph::BiharmonicEquations< DIM >::interpolated_u_biharmonic(), oomph::Node::is_hanging(), oomph::ProjectablePMLTimeHarmonicLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::local_equation(), oomph::ProjectableTimeHarmonicFourierDecomposedLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::local_equation(), oomph::ProjectableTimeHarmonicLinearElasticityElement< TIME_HARMONIC_LINEAR_ELAST_ELEMENT >::local_equation(), oomph::FixedVolumeSpineLineMarangoniFluidInterfaceElement< ELEMENT >::make_bounding_element(), MortaringValidationProblem< ELEMENT, NON_MORTAR_ELEMENT >::MortaringValidationProblem(), oomph::MyProblem::my_set_initial_condition(), oomph::RefineableImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::ncont_interpolated_values(), oomph::RefineableFSIImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::ncont_interpolated_values(), oomph::GenericLagrangeInterpolatedProjectableElement< ELEMENT >::nfields_for_projection(), oomph::FiniteElement::nodal_local_eqn(), oomph::ODEElement::nvalue(), oomph::MyProblem::output_ltes(), oomph::TreeBasedRefineableMeshBase::p_adapt_mesh(), 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::ProjectionProblem< PROJECTABLE_ELEMENT >::pin_all(), FlowAroundCylinderProblem< ELEMENT >::pin_all_base_flow_dofs(), oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::pin_all_non_pressure_dofs(), oomph::GeneralisedNewtonianNavierStokesEquations< DIM >::pin_all_non_pressure_dofs(), oomph::NavierStokesEquations< DIM >::pin_all_non_pressure_dofs(), oomph::SpaceTimeNavierStokesEquations< DIM >::pin_all_non_pressure_dofs(), oomph::SpaceTimeNavierStokesMixedOrderEquations< DIM >::pin_all_non_pressure_dofs(), NavierStokesProblem< ELEMENT >::pin_redundant_temporal_nodes(), oomph::MGSolver< DIM >::plot(), oomph::post_midpoint_update(), oomph::PRefineableQElement< 2, INITIAL_NNODE_1D >::pre_build(), oomph::PRefineableQElement< 3, INITIAL_NNODE_1D >::pre_build(), print_connectivity_matrix(), oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::project(), oomph::PseudoBucklingRing::PseudoBucklingRing(), read(), oomph::PRefineableQElement< 2, INITIAL_NNODE_1D >::rebuild_from_sons(), oomph::PRefineableQElement< 3, INITIAL_NNODE_1D >::rebuild_from_sons(), oomph::CopiedData::reset_copied_pointers(), oomph::BoundaryNode< NODE_TYPE >::reset_copied_pointers(), oomph::SolidICProblem::reset_original_state(), oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::reset_pin_status(), oomph::NavierStokesSchurComplementPreconditioner::reset_pin_status(), resize(), oomph::Node::resize(), oomph::BoundaryNode< NODE_TYPE >::resize(), oomph::FaceElement::resize_nodes(), oomph::ContinuationStorageScheme::set_consistent_pinned_values(), oomph::Problem::set_dofs(), oomph::FluidInterfaceElement::set_external_pressure_data(), oomph::Node::set_hanging_pt(), oomph::YoungLaplaceEquations::set_kappa(), oomph::ElasticallySupportedRingElement::set_load_pt(), ElasticFishBackElement::set_load_pt(), oomph::SolidICProblem::set_newmark_initial_condition_consistently(), oomph::Node::set_nonhanging(), oomph::Problem::set_pinned_values_to_zero(), oomph::WomersleyEquations< DIM >::set_pressure_gradient_pt(), oomph::MGSolver< DIM >::set_self_test_vector(), set_time_stepper(), FixSpineHeightElement::set_traded_pressure_data(), oomph::SpineVolumeConstraintPointElement< ELEMENT >::set_traded_pressure_data(), oomph::FixedVolumeSpineLineMarangoniFluidInterfaceElement< ELEMENT >::set_traded_pressure_data(), oomph::MyProblem::set_up_impulsive_initial_condition(), oomph::VolumeConstraintBoundingElement::set_volume_constraint_element(), oomph::FoepplvonKarmanEquations::set_volume_constraint_pressure_data_as_external_data(), oomph::SolidICProblem::setup_problem(), oomph::SelfStartingBDF2::shift_time_values(), oomph::IMRBase::shift_time_values(), oomph::Steady< NSTEPS >::shift_time_values(), oomph::NewmarkBDF< NSTEPS >::shift_time_values(), oomph::BDF< NSTEPS >::shift_time_values(), oomph::TR::shift_time_values(), oomph::MGSolver< DIM >::solve(), oomph::StraightLine::StraightLine(), oomph::TimeStepper::time_derivative(), oomph::MyProblem::trace_values(), TwoDDGProblem< ELEMENT >::TwoDDGProblem(), oomph::MyProblem::undo_segregated_pinning(), FlowAroundCylinderProblem< ELEMENT >::unpin_all_base_flow_dofs(), oomph::PRefineableGeneralisedNewtonianQCrouzeixRaviartElement< DIM >::unpin_elemental_pressure_dofs(), oomph::PRefineableQCrouzeixRaviartElement< DIM >::unpin_elemental_pressure_dofs(), value(), oomph::Node::value(), and oomph::Node::~Node().
|
delete |
Broken assignment operator.
|
inline |
Pin the i-th stored variable.
References eqn_number(), i, and Is_pinned.
Referenced by UnstructuredTorusProblem< ELEMENT >::actions_after_adapt(), SurfactantProblem< ELEMENT, INTERFACE_ELEMENT >::actions_after_adapt(), RefineableTwoDPoissonProblem< ELEMENT >::actions_before_newton_solve(), TiltedCavityProblem< ELEMENT >::actions_before_newton_solve(), EighthSpherePoissonProblem< ELEMENT >::actions_before_newton_solve(), AnnularDiskProblem< ELASTICITY_ELEMENT >::AnnularDiskProblem(), MovingBlockProblem< ELEMENT >::apply_boundary_conditions(), ExtrudedMovingCylinderProblem< TWO_D_ELEMENT, THREE_D_ELEMENT >::apply_boundary_conditions(), PMLStructuredCubicHelmholtz< ELEMENT >::apply_boundary_conditions(), PMLHelmholtzMGProblem< ELEMENT >::apply_boundary_conditions(), FlowAroundCylinderProblem< ELEMENT >::apply_boundary_conditions(), NavierStokesProblem< ELEMENT >::apply_boundary_conditions(), VorticityRecoveryProblem< ELEMENT >::apply_boundary_conditions(), PMLProblem< ELEMENT >::apply_boundary_conditions(), UnsteadyHeatProblem< ELEMENT >::apply_boundary_conditions(), PMLFourierDecomposedHelmholtzProblem< ELEMENT >::apply_zero_dirichlet_boundary_conditions(), AxisymFreeSurfaceNozzleAdvDiffRobinProblem< ELEMENT >::AxisymFreeSurfaceNozzleAdvDiffRobinProblem(), oomph::RefineableQElement< 3 >::build(), oomph::RefineableQElement< 2 >::build(), oomph::ChannelSpineMesh< ELEMENT >::build_channel_spine_mesh(), RectangularWomersleyImpedanceTube< ELEMENT >::build_mesh_and_apply_boundary_conditions(), CapProblem< ELEMENT >::CapProblem(), CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::CoatedDiskProblem(), CoatedSphereProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::CoatedSphereProblem(), InclinedPlaneProblem< ELEMENT, INTERFACE_ELEMENT >::complete_build(), UnstructuredFvKProblem< ELEMENT >::complete_problem_setup(), AxiPoroProblem< ELEMENT, TIMESTEPPER >::complete_problem_setup(), DarcyProblem< ELEMENT >::complete_problem_setup(), AxisymmetricVibratingShellProblem< ELEMENT >::complete_problem_setup(), VibratingShellProblem< ELEMENT >::complete_problem_setup(), ContactProblem< ELEMENT >::complete_problem_setup(), UnsteadyHeatMeltProblem< ELEMENT >::complete_problem_setup(), CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::complete_problem_setup(), CoatedSphereProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::complete_problem_setup(), UnstructuredPoissonProblem< ELEMENT >::complete_problem_setup(), UnstructuredImmersedEllipseProblem< ELEMENT >::complete_problem_setup(), TwoLayerInterfaceProblem< ELEMENT >::complete_problem_setup(), UnstructuredFluidProblem< ELEMENT >::complete_problem_setup(), BubbleInChannelProblem< ELEMENT >::complete_problem_setup(), DropInChannelProblem< ELEMENT >::complete_problem_setup(), ElasticAnnulusProblem< ELASTICITY_ELEMENT >::complete_problem_setup(), FourierDecomposedTimeHarmonicLinearElasticityProblem< ELEMENT >::complete_problem_setup(), RingWithTRibProblem< ELASTICITY_ELEMENT >::complete_problem_setup(), oomph::NodeNodeMortaringElement::construct_lagrange_multipliers(), oomph::TwoDimensionalPMLHelper::create_bottom_left_pml_mesh(), oomph::TwoDimensionalPMLHelper::create_bottom_pml_mesh(), oomph::TwoDimensionalPMLHelper::create_bottom_right_pml_mesh(), ContactProblem< ELEMENT >::create_displ_imposition_elements(), UnstructuredTorusProblem< ELEMENT >::create_lagrange_multiplier_elements(), FSICollapsibleChannelProblem< ELEMENT >::create_lagrange_multiplier_elements(), PseudoElasticCollapsibleChannelProblem< FLUID_ELEMENT, SOLID_ELEMENT >::create_lagrange_multiplier_elements(), FSIChannelWithLeafletProblem< ELEMENT >::create_lagrange_multiplier_elements(), UnstructuredFSIProblem< FLUID_ELEMENT, SOLID_ELEMENT >::create_lagrange_multiplier_elements(), PrescribedBoundaryDisplacementProblem< ELEMENT >::create_lagrange_multiplier_elements(), UnstructuredImmersedEllipseProblem< ELEMENT >::create_lagrange_multiplier_elements(), UnstructuredFluidProblem< ELEMENT >::create_lagrange_multiplier_elements(), CollapsibleChannelProblem< ELEMENT >::create_lagrange_multiplier_elements(), oomph::TwoDimensionalPMLHelper::create_left_pml_mesh(), TiltedCavityProblem< ELEMENT >::create_parall_outflow_lagrange_elements(), oomph::TwoDimensionalPMLHelper::create_right_pml_mesh(), oomph::TwoDimensionalPMLHelper::create_top_left_pml_mesh(), oomph::TwoDimensionalPMLHelper::create_top_pml_mesh(), oomph::TwoDimensionalPMLHelper::create_top_right_pml_mesh(), oomph::SurfaceMeltElement< ELEMENT >::disable_melting(), oomph::UnsteadyHeatFluxPseudoMeltElement< ELEMENT >::disable_melting(), ElasticInterfaceProblem< ELEMENT, TIMESTEPPER >::ElasticInterfaceProblem(), ElasticRingProblem< ELEMENT >::ElasticRingProblem(), oomph::LinearisedAxisymmetricQTaylorHoodElement::fix_cosine_component_of_pressure(), oomph::LinearisedAxisymmetricQCrouzeixRaviartElement::fix_cosine_component_of_pressure(), oomph::LinearisedAxisymmetricQTaylorHoodElement::fix_pressure(), oomph::AxisymmetricQTaylorHoodElement::fix_pressure(), oomph::GeneralisedNewtonianAxisymmetricQTaylorHoodElement::fix_pressure(), oomph::LinearisedQTaylorHoodElement::fix_pressure(), oomph::AxisymmetricTCrouzeixRaviartElement::fix_pressure(), oomph::AxisymmetricTTaylorHoodElement::fix_pressure(), oomph::GeneralisedNewtonianAxisymmetricTCrouzeixRaviartElement::fix_pressure(), oomph::GeneralisedNewtonianAxisymmetricTTaylorHoodElement::fix_pressure(), oomph::GeneralisedNewtonianQCrouzeixRaviartElement< DIM >::fix_pressure(), oomph::GeneralisedNewtonianQTaylorHoodElement< DIM >::fix_pressure(), oomph::PRefineableGeneralisedNewtonianQCrouzeixRaviartElement< DIM >::fix_pressure(), oomph::GeneralisedNewtonianTCrouzeixRaviartElement< DIM >::fix_pressure(), oomph::GeneralisedNewtonianTTaylorHoodElement< DIM >::fix_pressure(), oomph::QCrouzeixRaviartElement< DIM >::fix_pressure(), oomph::QTaylorHoodElement< DIM >::fix_pressure(), oomph::PRefineableQCrouzeixRaviartElement< DIM >::fix_pressure(), oomph::TCrouzeixRaviartElement< DIM >::fix_pressure(), oomph::TTaylorHoodElement< DIM >::fix_pressure(), oomph::PolarCrouzeixRaviartElement::fix_pressure(), oomph::PolarTaylorHoodElement::fix_pressure(), oomph::QTaylorHoodSpaceTimeElement< DIM >::fix_pressure(), oomph::QTaylorHoodMixedOrderSpaceTimeElement< DIM >::fix_pressure(), oomph::QSphericalCrouzeixRaviartElement::fix_pressure(), oomph::QSphericalTaylorHoodElement::fix_pressure(), oomph::LinearisedAxisymmetricQCrouzeixRaviartElement::fix_pressure(), oomph::AxisymmetricQCrouzeixRaviartElement::fix_pressure(), oomph::GeneralisedNewtonianAxisymmetricQCrouzeixRaviartElement::fix_pressure(), oomph::LinearisedQCrouzeixRaviartElement::fix_pressure(), oomph::LinearisedAxisymmetricQTaylorHoodElement::fix_sine_component_of_pressure(), oomph::LinearisedAxisymmetricQCrouzeixRaviartElement::fix_sine_component_of_pressure(), oomph::AxisymQPVDElementWithPressure::fix_solid_pressure(), oomph::QPVDElementWithPressure< DIM >::fix_solid_pressure(), oomph::QPVDElementWithContinuousPressure< DIM >::fix_solid_pressure(), oomph::TPVDElementWithContinuousPressure< DIM >::fix_solid_pressure(), FluxPoissonProblem< ELEMENT >::FluxPoissonProblem(), FourierDecomposedTimeHarmonicLinearElasticityProblem< ELEMENT >::FourierDecomposedTimeHarmonicLinearElasticityProblem(), FreeSurfaceRotationProblem< ELEMENT >::FreeSurfaceRotationProblem(), FSICollapsibleChannelProblem< ELEMENT >::FSICollapsibleChannelProblem(), FSIRingProblem::FSIRingProblem(), InterfaceProblem< ELEMENT, TIMESTEPPER >::InterfaceProblem(), MeltSpinningProblem< ELEMENT >::MeltSpinningProblem(), OneDPoissonProblem< ELEMENT >::OneDPoissonProblem(), oomph::PoissonSmoothMesh< POISSON_ELEMENT >::operator()(), OrrSommerfeldProblem< ELEMENT >::OrrSommerfeldProblem(), OscRingNStProblem< ELEMENT >::OscRingNStProblem(), oomph::PRefineableQElement< 2, INITIAL_NNODE_1D >::p_refine(), oomph::PRefineableQElement< 3, INITIAL_NNODE_1D >::p_refine(), PeriodicLoadProblem< ELEMENT >::PeriodicLoadProblem(), PerturbedStateProblem< BASE_ELEMENT, PERTURBED_ELEMENT >::PerturbedStateProblem(), oomph::ProjectionProblem< PROJECTABLE_ELEMENT >::pin_all(), FlowAroundCylinderProblem< ELEMENT >::pin_all_base_flow_dofs(), oomph::AxisymmetricTTaylorHoodElement::pin_all_nodal_pressure_dofs(), oomph::GeneralisedNewtonianAxisymmetricTTaylorHoodElement::pin_all_nodal_pressure_dofs(), oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::pin_all_non_pressure_dofs(), oomph::GeneralisedNewtonianNavierStokesEquations< DIM >::pin_all_non_pressure_dofs(), oomph::NavierStokesEquations< DIM >::pin_all_non_pressure_dofs(), oomph::SpaceTimeNavierStokesEquations< DIM >::pin_all_non_pressure_dofs(), oomph::SpaceTimeNavierStokesMixedOrderEquations< DIM >::pin_all_non_pressure_dofs(), oomph::BrethertonSpineMesh< ELEMENT, INTERFACE_ELEMENT >::pin_all_spines(), oomph::StreamfunctionProblem::pin_boundaries(), oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::pin_c(), oomph::ImmersedRigidBodyElement::pin_centre_of_mass_coordinate(), oomph::RefineableLinearisedAxisymmetricQTaylorHoodElement::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableAxisymmetricQTaylorHoodElement::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableGeneralisedNewtonianAxisymmetricQTaylorHoodElement::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableGeneralisedNewtonianQTaylorHoodElement< DIM >::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableLinearisedQTaylorHoodElement::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableQTaylorHoodElement< DIM >::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineablePolarTaylorHoodElement::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableQTaylorHoodSpaceTimeElement< DIM >::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableQTaylorHoodMixedOrderSpaceTimeElement< DIM >::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableQSphericalTaylorHoodElement::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableQPVDElementWithContinuousPressure< DIM >::pin_elemental_redundant_nodal_solid_pressures(), oomph::TAxisymmetricPoroelasticityElement< ORDER >::pin_p_value(), oomph::TRaviartThomasDarcyElement< ORDER >::pin_p_value(), oomph::TPoroelasticityElement< ORDER >::pin_p_value(), oomph::SolidNode::pin_position(), oomph::TAxisymmetricPoroelasticityElement< ORDER >::pin_q_edge_value(), oomph::TAxisymmetricPoroelasticityElement< ORDER >::pin_q_internal_value(), oomph::TRaviartThomasDarcyElement< ORDER >::pin_q_internal_value(), oomph::TPoroelasticityElement< ORDER >::pin_q_internal_value(), oomph::LinearisedQCrouzeixRaviartElement::pin_real_or_imag(), NavierStokesProblem< ELEMENT >::pin_redundant_temporal_nodes(), UnsteadyHeatProblem< ELEMENT >::pin_redundant_temporal_nodes(), oomph::ImmersedRigidBodyElement::pin_rotation_angle(), oomph::VorticitySmootherElement< ELEMENT >::pin_smoothed_vorticity(), oomph::ProjectableDarcyElement< DARCY_ELEMENT >::pin_superfluous_darcy_dofs(), oomph::StreamfunctionProblem::pin_velocities(), oomph::ElasticallySupportedRingElement::pin_yc(), PressureWaveFSIProblem< FLUID_ELEMENT, SOLID_ELEMENT >::PressureWaveFSIProblem(), PseudoSolidCapProblem< ELEMENT >::PseudoSolidCapProblem(), QuarterCircleDrivenCavityProblem< ELEMENT >::QuarterCircleDrivenCavityProblem(), QuarterCircleDrivenCavityProblem2< ELEMENT >::QuarterCircleDrivenCavityProblem2(), oomph::PRefineableQElement< 2, INITIAL_NNODE_1D >::rebuild_from_sons(), oomph::PRefineableQElement< 3, INITIAL_NNODE_1D >::rebuild_from_sons(), RefineableFishPoissonProblem< ELEMENT >::RefineableFishPoissonProblem(), RefineablePeriodicLoadProblem< ELEMENT >::RefineablePeriodicLoadProblem(), RisingBubbleProblem< ELEMENT >::RisingBubbleProblem(), run_navier_stokes_outflow(), run_prescribed_flux(), run_prescribed_pressure_gradient(), AxisymmetricLinearElasticityProblem< ELEMENT, TIMESTEPPER >::set_boundary_conditions(), UnstructuredFluidProblem< ELEMENT >::set_boundary_conditions(), oomph::SolidICProblem::setup_problem(), SolidFreeSurfaceRotationProblem< ELEMENT >::SolidFreeSurfaceRotationProblem(), ShellProblem< ELEMENT >::solve(), SphericalSteadyRotationProblem< ELEMENT >::SphericalSteadyRotationProblem(), StefanBoltzmannProblem< ELEMENT >::StefanBoltzmannProblem(), SurfactantProblem< ELEMENT, INTERFACE_ELEMENT >::SurfactantProblem(), ThermalProblem< ELEMENT >::ThermalProblem(), TiltedCavityProblem< ELEMENT >::TiltedCavityProblem(), TurekProblem< FLUID_ELEMENT, SOLID_ELEMENT >::TurekProblem(), TwoDDGProblem< ELEMENT >::TwoDDGProblem(), UnstructuredFSIProblem< FLUID_ELEMENT, SOLID_ELEMENT >::UnstructuredFSIProblem(), oomph::AxisymFoepplvonKarmanEquations::use_linear_bending_model(), oomph::DisplacementBasedFoepplvonKarmanEquations::use_linear_bending_model(), oomph::FoepplvonKarmanEquations::use_linear_bending_model(), oomph::WomersleyMesh< WOMERSLEY_ELEMENT >::WomersleyMesh(), oomph::WomersleyProblem< ELEMENT, DIM >::WomersleyProblem(), and YoungLaplaceProblem< ELEMENT >::YoungLaplaceProblem().
|
inline |
Pin all the stored variables.
References Eqn_number, i, Is_pinned, and Nvalue.
Referenced by oomph::Node::pin_all(), oomph::SolidNode::pin_all(), oomph::LinearisedQCrouzeixRaviartElement::pin_pressure_normalisation_dofs(), and oomph::LinearisedQCrouzeixRaviartElement::pin_real_or_imag().
Check that the arguments are within the range of the stored data values and timesteps.
//////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////// Private function to check that the arguments are within the range of the stored data values and timesteps.
References i, ntstorage(), Nvalue, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and plotPSD::t.
Referenced by eqn_number(), eqn_number_pt(), set_value(), value(), and value_pt().
void oomph::Data::read | ( | std::ifstream & | restart_file | ) |
Read data object from a file.
Read data object from file.
References i, j, oomph::Node::ndim(), ntstorage(), nvalue(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, set_value(), oomph::Global_string_for_annotation::string(), plotPSD::t, and oomph::Node::x().
Referenced by oomph::PerturbedSpineMesh::read(), oomph::Mesh::read(), oomph::Node::read(), and oomph::SpineMesh::read().
|
private |
Remove the pointer data_pt from the array Copy_of_data_pt. This should be used whenever copies of the data are deleted.
Remove the pointer data_pt from the internal storage used to keep track of copies
References Copy_of_data_pt, i, Ncopies, OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.
Referenced by oomph::BoundaryNode< NODE_TYPE >::~BoundaryNode(), oomph::CopiedData::~CopiedData(), and oomph::HijackedData::~HijackedData().
|
protectedvirtual |
Helper function that should be overloaded in derived classes that can contain copies of Data. The function must reset the internal pointers to the copied data. This is used when resizing data to ensure that all the pointers remain valid. The default implementation throws an error beacause Data cannot be a copy.
Helper function that should be overloaded in classes that contain copies of Data. The function must reset the internal pointers to the copied data. This is used when resizing data to ensure that all the pointers remain valid. The base Data class cannot be a copy, so throw an error
Reimplemented in oomph::CopiedData, and oomph::HijackedData.
References OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.
Referenced by resize(), and set_time_stepper().
|
virtual |
Change (increase) the number of values that may be stored.
Increase the number of data values stored, useful when adding additional data at a node, almost always Lagrange multipliers. Note if any of the unresized data is copied, then we assume all the resized data is copied from the same node as the unresized data.
Reimplemented in oomph::Node, oomph::CopiedData, and oomph::HijackedData.
References Copy_of_data_pt, Eqn_number, i, Is_unclassified, Ncopies, ntstorage(), Nvalue, nvalue(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, reset_copied_pointers(), plotPSD::t, and Value.
Referenced by oomph::PRefineableGeneralisedNewtonianQCrouzeixRaviartElement< DIM >::PRefineableGeneralisedNewtonianQCrouzeixRaviartElement(), oomph::PRefineableQCrouzeixRaviartElement< DIM >::PRefineableQCrouzeixRaviartElement(), and oomph::Node::resize().
unsigned oomph::Data::self_test | ( | ) |
Self-test: Have all values been classified as pinned/unpinned? Return 0 if OK.
References Eqn_number, i, Is_unclassified, Nvalue, and oomph::oomph_info.
Referenced by oomph::AlgebraicNode::self_test().
void oomph::Data::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 will not be altered
Set a new TimeStepper be resizing the appropriate storage. Equation numbering (if already performed) will be unaffected. The current (zero) values will be unaffected, but all other entries will be set to zero.
References Copy_of_data_pt, i, is_a_copy(), Ncopies, ntstorage(), oomph::TimeStepper::ntstorage(), nvalue(), reset_copied_pointers(), plotPSD::t, Time_stepper_pt, time_stepper_pt(), and Value.
Referenced by oomph::PeriodicOrbitAssemblyHandler< NNODE_1D >::adapt_temporal_mesh(), oomph::GeneralisedElement::set_internal_data_time_stepper(), oomph::SolidNode::set_position_time_stepper(), oomph::AxisymmetricPoroelasticityEquations::set_q_internal_timestepper(), and oomph::PoroelasticityEquations< DIM >::set_q_internal_timestepper().
Set the i-th stored data value to specified value. The only reason that we require an explicit set function is because we redefine value() in the Node class to interpolate the values for nodes that are hanging and so we cannot return a reference to the value in this case.
References i, range_check(), and Value.
Referenced by SurfactantProblem< ELEMENT, INTERFACE_ELEMENT >::actions_after_adapt(), ContactProblem< ELEMENT >::actions_before_implicit_timestep(), StefanBoltzmannProblem< ELEMENT >::actions_before_implicit_timestep(), UnsteadyHeatProblem< ELEMENT >::actions_before_implicit_timestep(), FSIDrivenCavityProblem< ELEMENT >::actions_before_implicit_timestep(), LinearWaveProblem< ELEMENT, TIMESTEPPER >::actions_before_implicit_timestep(), RefineableUnsteadyHeatProblem< ELEMENT >::actions_before_implicit_timestep(), SurfactantProblem< ELEMENT, INTERFACE_ELEMENT >::actions_before_implicit_timestep(), InclinedPlaneProblem< ELEMENT, INTERFACE_ELEMENT >::actions_before_implicit_timestep(), SphericalSpinUpProblem< ELEMENT >::actions_before_implicit_timestep(), oomph::WomersleyProblem< ELEMENT, DIM >::actions_before_implicit_timestep(), RefineableAdvectionDiffusionProblem< ELEMENT >::actions_before_newton_solve(), TwoMeshFluxAdvectionDiffusionProblem< ELEMENT >::actions_before_newton_solve(), SUPGAdvectionDiffusionProblem< ELEMENT >::actions_before_newton_solve(), OneDPoissonProblem< ELEMENT >::actions_before_newton_solve(), PoissonProblem< ELEMENT >::actions_before_newton_solve(), MultiPoissonProblem< ELEMENT >::actions_before_newton_solve(), RefineableTwoDPoissonProblem< ELEMENT >::actions_before_newton_solve(), RefineableTwoMeshFluxPoissonProblem< ELEMENT >::actions_before_newton_solve(), RefineableConvectionProblem< NST_ELEMENT, AD_ELEMENT >::actions_before_newton_solve(), FluxPoissonMGProblem< ELEMENT, MESH >::actions_before_newton_solve(), UnitCubePoissonMGProblem< ELEMENT, MESH >::actions_before_newton_solve(), TiltedCavityProblem< ELEMENT >::actions_before_newton_solve(), FpTestProblem::actions_before_newton_solve(), FallingBlockProblem< ELEMENT >::actions_before_newton_solve(), RisingBubbleProblem< ELEMENT >::actions_before_newton_solve(), RefineablePoissonProblem< ELEMENT >::actions_before_newton_solve(), EighthSpherePoissonProblem< ELEMENT >::actions_before_newton_solve(), RefineableOneDPoissonProblem< ELEMENT >::actions_before_newton_solve(), PRefineableOneDPoissonProblem< ELEMENT >::actions_before_newton_solve(), FluxPoissonProblem< ELEMENT >::actions_before_newton_solve(), TwoMeshFluxPoissonProblem< ELEMENT >::actions_before_newton_solve(), YoungLaplaceProblem< ELEMENT >::actions_before_newton_solve(), ThreeDPoissonProblem< ELEMENT >::actions_before_newton_solve(), TestPoissonProblem< ELEMENT >::actions_before_newton_solve(), TestRefineablePoissonProblem< ELEMENT >::actions_before_newton_solve(), oomph::StreamfunctionProblem::actions_before_solve(), oomph::TreeBasedRefineableMeshBase::adapt_mesh(), oomph::PeriodicOrbitAssemblyHandler< NNODE_1D >::adapt_temporal_mesh(), AnnularDiskProblem< ELASTICITY_ELEMENT >::AnnularDiskProblem(), MovingBlockProblem< ELEMENT >::apply_boundary_conditions(), ExtrudedMovingCylinderProblem< TWO_D_ELEMENT, THREE_D_ELEMENT >::apply_boundary_conditions(), UnstructuredPoissonProblem< ELEMENT >::apply_boundary_conditions(), PoissonProblem< ELEMENT >::apply_boundary_conditions(), PMLStructuredCubicHelmholtz< ELEMENT >::apply_boundary_conditions(), PMLHelmholtzMGProblem< ELEMENT >::apply_boundary_conditions(), FlowAroundCylinderProblem< ELEMENT >::apply_boundary_conditions(), NavierStokesProblem< ELEMENT >::apply_boundary_conditions(), VorticityRecoveryProblem< ELEMENT >::apply_boundary_conditions(), PMLProblem< ELEMENT >::apply_boundary_conditions(), UnsteadyHeatProblem< ELEMENT >::apply_boundary_conditions(), TwoDDGProblem< ELEMENT >::apply_boundary_conditions(), TwoDDGProblem< ELEMENT >::apply_initial_conditions(), oomph::FSI_functions::apply_no_slip_on_moving_wall(), PMLFourierDecomposedHelmholtzProblem< ELEMENT >::apply_zero_dirichlet_boundary_conditions(), oomph::SelfStartingBDF2::assign_initial_data_values(), oomph::PeriodicOrbitTimeDiscretisation::assign_initial_data_values(), oomph::Steady< NSTEPS >::assign_initial_data_values(), oomph::BDF< NSTEPS >::assign_initial_data_values(), oomph::Newmark< NSTEPS >::assign_initial_data_values(), oomph::Newmark< 2 >::assign_initial_data_values(), oomph::Newmark< NSTEPS >::assign_initial_data_values_stage1(), oomph::Newmark< NSTEPS >::assign_initial_data_values_stage2(), oomph::SelfStartingBDF2::assign_initial_values_impulsive(), oomph::Steady< NSTEPS >::assign_initial_values_impulsive(), oomph::Newmark< NSTEPS >::assign_initial_values_impulsive(), oomph::BDF< NSTEPS >::assign_initial_values_impulsive(), VorticityRecoveryProblem< ELEMENT >::assign_synthetic_veloc_field(), oomph::StreamfunctionProblem::assign_velocities(), AxisymFreeSurfaceNozzleAdvDiffRobinProblem< ELEMENT >::AxisymFreeSurfaceNozzleAdvDiffRobinProblem(), oomph::BackupMeshForProjection< GEOMETRIC_ELEMENT >::BackupMeshForProjection(), oomph::RefineableQElement< 3 >::build(), oomph::RefineableQElement< 1 >::build(), oomph::RefineableQElement< 2 >::build(), oomph::BrickFromTetMesh< ELEMENT >::build_mesh(), oomph::BDF< NSTEPS >::calculate_predicted_values(), oomph::TR::calculate_predicted_values(), oomph::SelfStartingBDF2::calculate_predicted_values_bdf2(), CapProblem< ELEMENT >::CapProblem(), CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::CoatedDiskProblem(), CoatedSphereProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::CoatedSphereProblem(), DarcyProblem< ELEMENT >::complete_problem_setup(), AxisymmetricVibratingShellProblem< ELEMENT >::complete_problem_setup(), VibratingShellProblem< ELEMENT >::complete_problem_setup(), CoatedDiskProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::complete_problem_setup(), CoatedSphereProblem< ELASTICITY_ELEMENT, HELMHOLTZ_ELEMENT >::complete_problem_setup(), TwoLayerInterfaceProblem< ELEMENT >::complete_problem_setup(), UnstructuredFluidProblem< ELEMENT >::complete_problem_setup(), BubbleInChannelProblem< ELEMENT >::complete_problem_setup(), DropInChannelProblem< ELEMENT >::complete_problem_setup(), ElasticAnnulusProblem< ELASTICITY_ELEMENT >::complete_problem_setup(), FourierDecomposedTimeHarmonicLinearElasticityProblem< ELEMENT >::complete_problem_setup(), RingWithTRibProblem< ELASTICITY_ELEMENT >::complete_problem_setup(), copy(), oomph::LinearisedQCrouzeixRaviartElement::copy_efunction_to_normalisation(), oomph::BackupMeshForProjection< GEOMETRIC_ELEMENT >::copy_onto_original_mesh(), oomph::TwoDimensionalPMLHelper::create_bottom_left_pml_mesh(), oomph::TwoDimensionalPMLHelper::create_bottom_pml_mesh(), oomph::TwoDimensionalPMLHelper::create_bottom_right_pml_mesh(), oomph::TwoDimensionalPMLHelper::create_left_pml_mesh(), oomph::TwoDimensionalPMLHelper::create_right_pml_mesh(), oomph::TwoDimensionalPMLHelper::create_top_left_pml_mesh(), oomph::TwoDimensionalPMLHelper::create_top_pml_mesh(), oomph::TwoDimensionalPMLHelper::create_top_right_pml_mesh(), oomph::Mesh::delete_all_external_storage(), ElasticInterfaceProblem< ELEMENT, TIMESTEPPER >::ElasticInterfaceProblem(), EntryFlowProblem< ELEMENT >::EntryFlowProblem(), oomph::LinearisedAxisymmetricQTaylorHoodElement::fix_cosine_component_of_pressure(), oomph::LinearisedAxisymmetricQCrouzeixRaviartElement::fix_cosine_component_of_pressure(), oomph::LinearisedAxisymmetricQTaylorHoodElement::fix_pressure(), oomph::AxisymmetricQTaylorHoodElement::fix_pressure(), oomph::GeneralisedNewtonianAxisymmetricQTaylorHoodElement::fix_pressure(), oomph::LinearisedQTaylorHoodElement::fix_pressure(), oomph::AxisymmetricTCrouzeixRaviartElement::fix_pressure(), oomph::AxisymmetricTTaylorHoodElement::fix_pressure(), oomph::GeneralisedNewtonianAxisymmetricTCrouzeixRaviartElement::fix_pressure(), oomph::GeneralisedNewtonianAxisymmetricTTaylorHoodElement::fix_pressure(), oomph::GeneralisedNewtonianQCrouzeixRaviartElement< DIM >::fix_pressure(), oomph::GeneralisedNewtonianQTaylorHoodElement< DIM >::fix_pressure(), oomph::PRefineableGeneralisedNewtonianQCrouzeixRaviartElement< DIM >::fix_pressure(), oomph::GeneralisedNewtonianTCrouzeixRaviartElement< DIM >::fix_pressure(), oomph::GeneralisedNewtonianTTaylorHoodElement< DIM >::fix_pressure(), oomph::QCrouzeixRaviartElement< DIM >::fix_pressure(), oomph::QTaylorHoodElement< DIM >::fix_pressure(), oomph::PRefineableQCrouzeixRaviartElement< DIM >::fix_pressure(), oomph::TCrouzeixRaviartElement< DIM >::fix_pressure(), oomph::TTaylorHoodElement< DIM >::fix_pressure(), oomph::QTaylorHoodSpaceTimeElement< DIM >::fix_pressure(), oomph::QTaylorHoodMixedOrderSpaceTimeElement< DIM >::fix_pressure(), oomph::QSphericalCrouzeixRaviartElement::fix_pressure(), oomph::QSphericalTaylorHoodElement::fix_pressure(), oomph::LinearisedAxisymmetricQCrouzeixRaviartElement::fix_pressure(), oomph::AxisymmetricQCrouzeixRaviartElement::fix_pressure(), oomph::GeneralisedNewtonianAxisymmetricQCrouzeixRaviartElement::fix_pressure(), oomph::LinearisedQCrouzeixRaviartElement::fix_pressure(), oomph::LinearisedAxisymmetricQTaylorHoodElement::fix_sine_component_of_pressure(), oomph::LinearisedAxisymmetricQCrouzeixRaviartElement::fix_sine_component_of_pressure(), oomph::AxisymQPVDElementWithPressure::fix_solid_pressure(), oomph::QPVDElementWithPressure< DIM >::fix_solid_pressure(), oomph::QPVDElementWithContinuousPressure< DIM >::fix_solid_pressure(), oomph::TPVDElementWithContinuousPressure< DIM >::fix_solid_pressure(), oomph::FoepplvonKarmanVolumeConstraintElement< ELEMENT, MESH >::FoepplvonKarmanVolumeConstraintElement(), for(), FourierDecomposedTimeHarmonicLinearElasticityProblem< ELEMENT >::FourierDecomposedTimeHarmonicLinearElasticityProblem(), FreeSurfaceRotationProblem< ELEMENT >::FreeSurfaceRotationProblem(), oomph::RefineableAxisymmetricQCrouzeixRaviartElement::further_build(), oomph::RefineableGeneralisedNewtonianAxisymmetricQCrouzeixRaviartElement::further_build(), oomph::RefineableLinearisedAxisymmetricQCrouzeixRaviartElement::further_build(), oomph::RefineableLinearisedQCrouzeixRaviartElement::further_build(), oomph::RefineablePolarCrouzeixRaviartElement::further_build(), oomph::RefineableQSphericalCrouzeixRaviartElement::further_build(), HeatedCircularPenetratorElement::HeatedCircularPenetratorElement(), oomph::ImposeFluxForWomersleyElement< DIM >::ImposeFluxForWomersleyElement(), RefineableYoungLaplaceProblem< ELEMENT >::increment_parameters(), InterfaceProblem< ELEMENT, TIMESTEPPER >::InterfaceProblem(), TwoDDGProblem< ELEMENT >::limit(), main(), MeltSpinningProblem< ELEMENT >::MeltSpinningProblem(), MortaringValidationProblem< ELEMENT, NON_MORTAR_ELEMENT >::MortaringValidationProblem(), oomph::MyProblem::my_set_initial_condition(), oomph::ODEProblem::my_set_initial_condition(), oomph::PoissonSmoothMesh< POISSON_ELEMENT >::operator()(), PolarNSProblem< ELEMENT >::output_streamfunction(), oomph::TreeBasedRefineableMeshBase::p_adapt_mesh(), oomph::PRefineableQElement< 1, INITIAL_NNODE_1D >::p_refine(), oomph::PRefineableQElement< 2, INITIAL_NNODE_1D >::p_refine(), oomph::PRefineableQElement< 3, INITIAL_NNODE_1D >::p_refine(), ElasticRingProblem< ELEMENT >::parameter_study(), PeriodicLoadProblem< ELEMENT >::PeriodicLoadProblem(), oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::pin_all_non_pressure_dofs(), oomph::TAxisymmetricPoroelasticityElement< ORDER >::pin_p_value(), oomph::TPoroelasticityElement< ORDER >::pin_p_value(), oomph::TAxisymmetricPoroelasticityElement< ORDER >::pin_q_edge_value(), oomph::TAxisymmetricPoroelasticityElement< ORDER >::pin_q_internal_value(), oomph::UnsteadyHeatFluxPseudoMeltElement< ELEMENT >::plot_residual_landscape(), oomph::post_midpoint_update(), oomph::PRefineableQElement< 2, INITIAL_NNODE_1D >::pre_build(), oomph::PRefineableQElement< 3, INITIAL_NNODE_1D >::pre_build(), PseudoSolidCapProblem< ELEMENT >::PseudoSolidCapProblem(), read(), BrethertonElement< ELEMENT >::reassign_inflow(), oomph::RefineableLinearisedAxisymmetricQCrouzeixRaviartElement::rebuild_from_sons(), oomph::RefineableAxisymmetricQCrouzeixRaviartElement::rebuild_from_sons(), oomph::RefineableGeneralisedNewtonianAxisymmetricQCrouzeixRaviartElement::rebuild_from_sons(), oomph::PRefineableGeneralisedNewtonianQCrouzeixRaviartElement< DIM >::rebuild_from_sons(), oomph::PRefineableQElement< 2, INITIAL_NNODE_1D >::rebuild_from_sons(), oomph::PRefineableQElement< 3, INITIAL_NNODE_1D >::rebuild_from_sons(), oomph::RefineableLinearisedQCrouzeixRaviartElement::rebuild_from_sons(), oomph::PRefineableQCrouzeixRaviartElement< DIM >::rebuild_from_sons(), oomph::RefineablePolarCrouzeixRaviartElement::rebuild_from_sons(), oomph::RefineableQSphericalCrouzeixRaviartElement::rebuild_from_sons(), oomph::VorticitySmoother< ELEMENT >::recover_vorticity(), RefineablePeriodicLoadProblem< ELEMENT >::RefineablePeriodicLoadProblem(), oomph::ImmersedRigidBodyTriangleMeshPolygon::reset_reference_configuration(), run_navier_stokes_outflow(), HeatedCircularPenetratorElement::set_angle(), AxisymmetricLinearElasticityProblem< ELEMENT, TIMESTEPPER >::set_boundary_conditions(), BaseStateProblem< BASE_ELEMENT >::set_boundary_conditions(), PerturbedStateProblem< BASE_ELEMENT, PERTURBED_ELEMENT >::set_boundary_conditions(), EultingSphereProblem< ELEMENT >::set_boundary_conditions(), RefineableElutingSphereProblem< ELEMENT >::set_boundary_conditions(), RefineableSphericalCouetteProblem< ELEMENT >::set_boundary_conditions(), RefineableSphericalSpinUpProblem< ELEMENT >::set_boundary_conditions(), SphericalSpinUpProblem< ELEMENT >::set_boundary_conditions(), ConvectionProblem< NST_ELEMENT, AD_ELEMENT >::set_boundary_conditions(), DDConvectionProblem< NST_ELEMENT, AD_ELEMENT >::set_boundary_conditions(), RefineableDDConvectionProblem< NST_ELEMENT, AD_ELEMENT >::set_boundary_conditions(), RefineableSphereConvectionProblem< ELEMENT >::set_boundary_conditions(), SurfactantProblem< ELEMENT, INTERFACE_ELEMENT >::set_boundary_conditions(), AxiPoroProblem< ELEMENT, TIMESTEPPER >::set_boundary_values(), UnstructuredImmersedEllipseProblem< ELEMENT >::set_boundary_velocity(), oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::set_c(), oomph::ContinuationStorageScheme::set_consistent_pinned_values(), oomph::Problem::set_dofs(), FlowAroundCylinderProblem< ELEMENT >::set_eigenvalue(), RefineableAdvectionDiffusionPipeProblem< ELEMENT >::set_initial_condition(), LinearWaveProblem< ELEMENT, TIMESTEPPER >::set_initial_condition(), RefineableOneDAdvectionDiffusionReactionProblem< ELEMENT >::set_initial_condition(), RefineableActivatorInhibitorProblem< ELEMENT >::set_initial_condition(), AxisymmetricLinearElasticityProblem< ELEMENT, TIMESTEPPER >::set_initial_conditions(), AdvectionProblem::set_initial_conditions(), EulerProblem::set_initial_conditions(), RefineableAdvectionDiffusionPipeProblem< ELEMENT >::set_inlet_concentration(), oomph::SurfaceMeltElement< ELEMENT >::set_lagrange_multiplier_pressure_to_zero(), oomph::SolidICProblem::set_newmark_initial_condition_consistently(), oomph::TAxisymmetricPoroelasticityElement< ORDER >::set_p_value(), oomph::TRaviartThomasDarcyElement< ORDER >::set_p_value(), oomph::Problem::set_pinned_values_to_zero(), oomph::TAxisymmetricPoroelasticityElement< ORDER >::set_q_edge(), oomph::TRaviartThomasDarcyElement< ORDER >::set_q_edge(), oomph::TAxisymmetricPoroelasticityElement< ORDER >::set_q_internal(), oomph::TRaviartThomasDarcyElement< ORDER >::set_q_internal(), oomph::QPVDElementWithPressure< DIM >::set_solid_p(), oomph::QPVDElementWithContinuousPressure< DIM >::set_solid_p(), oomph::TPVDElementWithContinuousPressure< DIM >::set_solid_p(), oomph::MyProblem::set_up_impulsive_initial_condition(), PolarNSProblem< ELEMENT >::setup_external_pressure(), oomph::SelfStartingBDF2::shift_time_values(), oomph::IMRBase::shift_time_values(), oomph::Steady< NSTEPS >::shift_time_values(), oomph::NewmarkBDF< NSTEPS >::shift_time_values(), oomph::BDF< NSTEPS >::shift_time_values(), oomph::TR::shift_time_values(), SolidFreeSurfaceRotationProblem< ELEMENT >::SolidFreeSurfaceRotationProblem(), SphericalSteadyRotationProblem< ELEMENT >::SphericalSteadyRotationProblem(), oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::spine_node_update_lower(), oomph::TwoLayerPerturbedSpineMesh< ELEMENT >::spine_node_update_upper(), SegregatedFSICollapsibleChannelProblem< ELEMENT >::steady_run(), StefanBoltzmannProblem< ELEMENT >::StefanBoltzmannProblem(), SurfactantProblem< ELEMENT, INTERFACE_ELEMENT >::SurfactantProblem(), RefineablePorousChannelProblem< ELEMENT >::symmetrise_eigenfunction_for_adaptive_pitchfork_tracking(), ThermalProblem< ELEMENT >::ThermalProblem(), TiltedCavityProblem< ELEMENT >::TiltedCavityProblem(), ActivatorInhibitorProblem< ELEMENT >::timestep(), TurekProblem< FLUID_ELEMENT, SOLID_ELEMENT >::TurekProblem(), TwoDDGMesh< ELEMENT >::TwoDDGMesh(), and TwoDDGProblem< ELEMENT >::TwoDDGProblem().
|
inline |
Set the t-th history value of the i-th stored data value to specified value.
References i, range_check(), plotPSD::t, and Value.
|
inline |
Return the pointer to the timestepper.
References Time_stepper_pt.
Referenced by SurfactantProblem< ELEMENT, INTERFACE_ELEMENT >::actions_after_adapt(), oomph::SurfactantTransportInterfaceElement::add_additional_residual_contributions_interface(), oomph::BackupMeshForProjection< GEOMETRIC_ELEMENT >::BackupMeshForProjection(), oomph::LinearElasticityEquationsBase< DIM >::body_force(), oomph::TimeHarmonicLinearElasticityEquationsBase< DIM >::body_force(), oomph::PVDEquationsBase< DIM >::body_force(), oomph::RefineableQElement< 3 >::build(), oomph::RefineableSolidQElement< 3 >::build(), oomph::RefineableQElement< 1 >::build(), oomph::RefineableQElement< 2 >::build(), oomph::RefineableSolidQElement< 2 >::build(), AxisymmetricVibratingShellProblem< ELEMENT >::complete_problem_setup(), VibratingShellProblem< ELEMENT >::complete_problem_setup(), TwoLayerInterfaceProblem< ELEMENT >::complete_problem_setup(), UnstructuredFluidProblem< ELEMENT >::complete_problem_setup(), BubbleInChannelProblem< ELEMENT >::complete_problem_setup(), DropInChannelProblem< ELEMENT >::complete_problem_setup(), oomph::Mesh::convert_to_boundary_node(), 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::PerturbedSpineLinearisedAxisymmetricFluidInterfaceElement< ELEMENT >::dH_dt(), oomph::FSIAxisymFoepplvonKarmanElement< NNODE_1D, FLUID_ELEMENT >::dposition_dt(), oomph::ImmersedRigidBodyElement::dposition_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::LinearisedAxisymmetricNavierStokesEquations::du_dt_lin_axi_nst(), 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::GeneralisedNewtonianAxisymmetricNavierStokesEquations::extrapolated_strain_rate(), oomph::IMRODEElement::fill_in_contribution_to_jacobian(), oomph::ODEElement::fill_in_contribution_to_jacobian(), oomph::IMRODEElement::fill_in_contribution_to_residuals(), oomph::ODEElement::fill_in_contribution_to_residuals(), 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::RefineableGeneralisedAxisymAdvectionDiffusionEquations::fill_in_generic_residual_contribution_cons_axisym_adv_diff(), oomph::PerturbedSpineLinearisedAxisymmetricFluidInterfaceElement< ELEMENT >::fill_in_generic_residual_contribution_interface(), oomph::LinearisedAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_lin_axi_nst(), 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::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::AxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::RefineableAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::GeneralisedNewtonianAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::RefineableGeneralisedNewtonianAxisymmetricNavierStokesEquations::get_dresidual_dnodal_coordinates(), oomph::RefineablePolarTaylorHoodElement::get_interpolated_values(), oomph::RefineablePolarCrouzeixRaviartElement::get_interpolated_values(), oomph::ImmersedRigidBodyElement::get_residuals_rigid_body_generic(), oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::get_wind_adv_diff_react(), oomph::FiniteElement::get_x(), RefineableUnsteadyHeatProblem< ELEMENT >::global_temporal_error_norm(), oomph::UnsteadyHeatBaseFaceElement< ELEMENT >::output(), oomph::UnsteadyHeatFluxPseudoMeltElement< ELEMENT >::output(), oomph::AxisymmetricLinearElasticityTractionElement< ELEMENT >::output(), oomph::AxisymmetricNavierStokesTractionElement< ELEMENT >::output(), oomph::AxisymmetricPoroelasticityTractionElement< ELEMENT >::output(), oomph::FSILinearisedAxisymPoroelasticTractionElement< POROELASTICITY_BULK_ELEMENT, NAVIER_STOKES_BULK_ELEMENT >::output(), oomph::MyTaylorHoodElement< DIM >::output(), oomph::ProjectableUnsteadyHeatSpaceTimeElement< UNSTEADY_HEAT_ELEMENT >::output(), oomph::ProjectableUnsteadyHeatMixedOrderSpaceTimeElement< UNSTEADY_HEAT_ELEMENT >::output(), oomph::ProjectableUnsteadyHeatElement< UNSTEADY_HEAT_ELEMENT >::output(), oomph::SolarRadiationBase::output_atmospheric_radiation(), oomph::ImmersedRigidBodyElement::output_centre_of_gravity(), oomph::MyProblem::output_ltes(), oomph::SurfaceMeltElement< ELEMENT >::output_melt(), oomph::StefanBoltzmannUnsteadyHeatFluxElement< ELEMENT >::output_stefan_boltzmann_radiation(), oomph::PRefineableQElement< 2, INITIAL_NNODE_1D >::pre_build(), oomph::PRefineableQElement< 3, INITIAL_NNODE_1D >::pre_build(), oomph::ImmersedRigidBodyTriangleMeshPolygon::reset_reference_configuration(), oomph::AxisymmetricNavierStokesTractionElement< ELEMENT >::scalar_value_paraview(), UnstructuredImmersedEllipseProblem< ELEMENT >::set_boundary_velocity(), oomph::ContinuationStorageScheme::set_consistent_pinned_values(), AdvectionProblem::set_initial_conditions(), EulerProblem::set_initial_conditions(), set_time_stepper(), oomph::IMRODEElement::time_interpolate_time(), oomph::IMRODEElement::time_interpolate_u(), oomph::TemplateFreeContactElementBase::traction_fct(), and oomph::ODEProblem::ts_pt().
|
inline |
Return the pointer to the timestepper (const version).
References Time_stepper_pt.
|
inline |
Unconstrain the i-th stored variable when make the data nonhanging. Only unconstrain if it was actually constrained (hanging)
References eqn_number(), i, Is_constrained, and Is_unclassified.
Referenced by oomph::Node::resize(), oomph::Node::set_hanging_pt(), oomph::Node::set_nonhanging(), and unconstrain_all().
|
inline |
Unconstrain all the stored variables when the data is made nonhanging.
References i, Nvalue, and unconstrain().
Referenced by oomph::SolidNode::unconstrain_positions().
|
inline |
Unpin the i-th stored variable.
References eqn_number(), i, and Is_unclassified.
Referenced by SurfactantProblem< ELEMENT, INTERFACE_ELEMENT >::actions_after_adapt(), oomph::ElasticallySupportedRingElement::ElasticallySupportedRingElement(), FpTestProblem::FpTestProblem(), oomph::RefineableLinearisedAxisymmetricQTaylorHoodElement::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableAxisymmetricQTaylorHoodElement::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableGeneralisedNewtonianAxisymmetricQTaylorHoodElement::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableGeneralisedNewtonianQTaylorHoodElement< DIM >::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableLinearisedQTaylorHoodElement::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableQTaylorHoodElement< DIM >::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineablePolarTaylorHoodElement::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableQTaylorHoodSpaceTimeElement< DIM >::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableQTaylorHoodMixedOrderSpaceTimeElement< DIM >::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableQSphericalTaylorHoodElement::pin_elemental_redundant_nodal_pressure_dofs(), oomph::RefineableQPVDElementWithContinuousPressure< DIM >::pin_elemental_redundant_nodal_solid_pressures(), oomph::SolidICProblem::reset_original_state(), oomph::DoubleBuoyantQCrouzeixRaviartElement< DIM >::unfix_pressure(), oomph::RefineableDoubleBuoyantQCrouzeixRaviartElement< DIM >::unfix_pressure(), oomph::BuoyantQCrouzeixRaviartElement< DIM >::unfix_pressure(), FlowAroundCylinderProblem< ELEMENT >::unpin_all_base_flow_dofs(), oomph::AxisymmetricTCrouzeixRaviartElement::unpin_all_internal_pressure_dofs(), oomph::GeneralisedNewtonianAxisymmetricTCrouzeixRaviartElement::unpin_all_internal_pressure_dofs(), oomph::AxisymmetricTTaylorHoodElement::unpin_all_nodal_pressure_dofs(), oomph::GeneralisedNewtonianAxisymmetricTTaylorHoodElement::unpin_all_nodal_pressure_dofs(), oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::unpin_c(), oomph::ImmersedRigidBodyElement::unpin_centre_of_mass_coordinate(), oomph::RefineableLinearisedAxisymmetricQTaylorHoodElement::unpin_elemental_pressure_dofs(), oomph::RefineableLinearisedAxisymmetricQCrouzeixRaviartElement::unpin_elemental_pressure_dofs(), oomph::RefineableAxisymmetricQTaylorHoodElement::unpin_elemental_pressure_dofs(), oomph::RefineableAxisymmetricQCrouzeixRaviartElement::unpin_elemental_pressure_dofs(), oomph::RefineableGeneralisedNewtonianAxisymmetricQTaylorHoodElement::unpin_elemental_pressure_dofs(), oomph::RefineableGeneralisedNewtonianAxisymmetricQCrouzeixRaviartElement::unpin_elemental_pressure_dofs(), oomph::RefineableGeneralisedNewtonianQTaylorHoodElement< DIM >::unpin_elemental_pressure_dofs(), oomph::RefineableGeneralisedNewtonianQCrouzeixRaviartElement< DIM >::unpin_elemental_pressure_dofs(), oomph::PRefineableGeneralisedNewtonianQCrouzeixRaviartElement< DIM >::unpin_elemental_pressure_dofs(), oomph::RefineableLinearisedQCrouzeixRaviartElement::unpin_elemental_pressure_dofs(), oomph::RefineableLinearisedQTaylorHoodElement::unpin_elemental_pressure_dofs(), oomph::RefineableQTaylorHoodElement< DIM >::unpin_elemental_pressure_dofs(), oomph::RefineableQCrouzeixRaviartElement< DIM >::unpin_elemental_pressure_dofs(), oomph::PRefineableQCrouzeixRaviartElement< DIM >::unpin_elemental_pressure_dofs(), oomph::RefineablePolarTaylorHoodElement::unpin_elemental_pressure_dofs(), oomph::RefineablePolarCrouzeixRaviartElement::unpin_elemental_pressure_dofs(), oomph::RefineableQTaylorHoodSpaceTimeElement< DIM >::unpin_elemental_pressure_dofs(), oomph::RefineableQTaylorHoodMixedOrderSpaceTimeElement< DIM >::unpin_elemental_pressure_dofs(), oomph::RefineableQSphericalTaylorHoodElement::unpin_elemental_pressure_dofs(), oomph::RefineableQSphericalCrouzeixRaviartElement::unpin_elemental_pressure_dofs(), oomph::RefineableQPVDElementWithPressure< DIM >::unpin_elemental_solid_pressure_dofs(), oomph::RefineableQPVDElementWithContinuousPressure< DIM >::unpin_elemental_solid_pressure_dofs(), oomph::QPVDElementWithPressure< DIM >::unpin_elemental_solid_pressure_dofs(), oomph::QPVDElementWithContinuousPressure< DIM >::unpin_elemental_solid_pressure_dofs(), oomph::TPVDElementWithContinuousPressure< DIM >::unpin_elemental_solid_pressure_dofs(), oomph::SolidNode::unpin_position(), oomph::AxisymmetricTTaylorHoodElement::unpin_proper_nodal_pressure_dofs(), oomph::GeneralisedNewtonianAxisymmetricTTaylorHoodElement::unpin_proper_nodal_pressure_dofs(), oomph::GeneralisedNewtonianTTaylorHoodElement< DIM >::unpin_proper_nodal_pressure_dofs(), oomph::TTaylorHoodElement< DIM >::unpin_proper_nodal_pressure_dofs(), oomph::LinearisedQCrouzeixRaviartElement::unpin_real_or_imag(), oomph::ImmersedRigidBodyElement::unpin_rotation_angle(), SurfactantProblem< ELEMENT, INTERFACE_ELEMENT >::unpin_surface(), oomph::ElasticallySupportedRingElement::unpin_yc(), and oomph::WomersleyMesh< WOMERSLEY_ELEMENT >::WomersleyMesh().
|
inline |
Unpin all the stored variables.
References Eqn_number, i, Is_unclassified, and Nvalue.
Referenced by oomph::Node::unpin_all(), oomph::SolidNode::unpin_all(), and oomph::LinearisedQCrouzeixRaviartElement::unpin_real_or_imag().
Return i-th stored value. This function is not virtual so that it can be inlined. This means that if we have an explicit pointer to a Data object Data* data_pt->value() always returns the "raw" stored value.
References i, range_check(), and Value.
Referenced by YoungLaplaceProblem< ELEMENT >::actions_before_newton_solve(), oomph::PeriodicOrbitAssemblyHandler< NNODE_1D >::adapt_temporal_mesh(), HeatedCircularPenetratorElement::angle(), oomph::ImmersedRigidBodyElement::apply_rigid_body_motion(), oomph::Newmark< 2 >::assign_initial_data_values(), oomph::Newmark< NSTEPS >::assign_initial_data_values_stage1(), oomph::Newmark< NSTEPS >::assign_initial_data_values_stage2(), oomph::SelfStartingBDF2::assign_initial_values_impulsive(), oomph::Steady< NSTEPS >::assign_initial_values_impulsive(), oomph::Newmark< NSTEPS >::assign_initial_values_impulsive(), oomph::BDF< NSTEPS >::assign_initial_values_impulsive(), oomph::GeneralisedElement::assign_local_eqn_numbers(), oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::c(), oomph::BDF< NSTEPS >::calculate_predicted_values(), oomph::TR::calculate_predicted_values(), oomph::SelfStartingBDF2::calculate_predicted_values_bdf2(), oomph::ImmersedRigidBodyElement::centre_of_gravity(), copy(), oomph::LinearisedQCrouzeixRaviartElement::copy_efunction_to_normalisation(), SpineGravityTractionElement< ELEMENT >::delta_p(), UnstructuredFvKProblem< ELEMENT >::doc_solution(), ElasticRingProblem< ELEMENT >::doc_solution(), YoungLaplaceProblem< ELEMENT >::doc_solution(), RefineableYoungLaplaceProblem< ELEMENT >::doc_solution(), SegregatedFSICollapsibleChannelProblem< ELEMENT >::doc_solution_steady(), oomph::ImmersedRigidBodyElement::dposition_dt(), dump(), oomph::LinearisedNavierStokesEigenfunctionNormalisationElement::eigenvalue(), oomph::ODEElement::fill_in_contribution_to_jacobian(), oomph::NodeElementSolidOnlyMortaringElement::fill_in_contribution_to_jacobian_mortared_nodes(), oomph::NodeNodeMortaringElement::fill_in_contribution_to_jacobian_mortared_nodes(), oomph::ConstraintElement::fill_in_contribution_to_residuals(), HeatedCircularPenetratorElement::fill_in_contribution_to_residuals(), oomph::ODEElement::fill_in_contribution_to_residuals(), oomph::PolarNavierStokesTractionElement< ELEMENT >::fill_in_generic_residual_contribution(), oomph::NavierStokesFluxControlElement< ELEMENT >::fill_in_generic_residual_contribution_fluid_traction(), oomph::NavierStokesWomersleyPressureControlElement::fill_in_generic_residual_contribution_pressure_control(), oomph::RefineableAxisymmetricQCrouzeixRaviartElement::further_build(), oomph::RefineableGeneralisedNewtonianAxisymmetricQCrouzeixRaviartElement::further_build(), oomph::RefineableGeneralisedNewtonianQCrouzeixRaviartElement< DIM >::further_build(), oomph::RefineableQCrouzeixRaviartElement< DIM >::further_build(), oomph::RefineablePolarCrouzeixRaviartElement::further_build(), oomph::Problem::get_dofs(), oomph::NavierStokesSurfaceDragTorqueElement< ELEMENT >::get_drag_and_torque(), oomph::YoungLaplaceEquations::get_kappa(), TestSoln::get_pressure(), oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::grad_p_bar(), oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::grad_u_bar(), RefineableYoungLaplaceProblem< ELEMENT >::increment_parameters(), oomph::ElasticallySupportedRingElement::load(), oomph::NavierStokesSurfaceDragTorqueElement< ELEMENT >::output(), oomph::ImmersedRigidBodyElement::output_centre_of_gravity(), oomph::AxisymmetricQCrouzeixRaviartElement::p_axi_nst(), oomph::AxisymmetricTCrouzeixRaviartElement::p_axi_nst(), oomph::GeneralisedNewtonianAxisymmetricQCrouzeixRaviartElement::p_axi_nst(), oomph::GeneralisedNewtonianAxisymmetricTCrouzeixRaviartElement::p_axi_nst(), oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::p_bar(), oomph::LinearisedAxisymmetricQCrouzeixRaviartElement::p_lin_axi_nst(), oomph::LinearisedAxisymmetricQCrouzeixRaviartElement::p_linearised_axi_nst(), oomph::LinearisedQCrouzeixRaviartElement::p_linearised_nst(), oomph::GeneralisedNewtonianQCrouzeixRaviartElement< DIM >::p_nst(), oomph::PRefineableGeneralisedNewtonianQCrouzeixRaviartElement< DIM >::p_nst(), oomph::GeneralisedNewtonianTCrouzeixRaviartElement< DIM >::p_nst(), oomph::QCrouzeixRaviartElement< DIM >::p_nst(), oomph::PRefineableQCrouzeixRaviartElement< DIM >::p_nst(), oomph::TCrouzeixRaviartElement< DIM >::p_nst(), oomph::QSphericalCrouzeixRaviartElement::p_spherical_nst(), oomph::VolumeConstraintElement::p_traded(), oomph::TAxisymmetricPoroelasticityElement< ORDER >::p_value(), oomph::TPoroelasticityElement< ORDER >::p_value(), oomph::FluidInterfaceElement::pext(), oomph::NodeElementSolidOnlyMortaringElement::position_in_element(), oomph::post_midpoint_update(), Global_Physical_Variables::press_load(), oomph::TAxisymmetricPoroelasticityElement< ORDER >::q_internal(), oomph::TRaviartThomasDarcyElement< ORDER >::q_internal(), oomph::TPoroelasticityElement< ORDER >::q_internal(), oomph::Node::raw_value(), oomph::RefineableLinearisedAxisymmetricQCrouzeixRaviartElement::rebuild_from_sons(), oomph::RefineableAxisymmetricQCrouzeixRaviartElement::rebuild_from_sons(), oomph::RefineableGeneralisedNewtonianAxisymmetricQCrouzeixRaviartElement::rebuild_from_sons(), oomph::RefineableLinearisedQCrouzeixRaviartElement::rebuild_from_sons(), oomph::RefineablePolarCrouzeixRaviartElement::rebuild_from_sons(), oomph::RefineableQSphericalCrouzeixRaviartElement::rebuild_from_sons(), oomph::PseudoBucklingRingElement::reference_pressure(), oomph::RefineableNavierStokesFluxControlElement< ELEMENT >::refineable_fill_in_generic_residual_contribution_fluid_traction(), oomph::ImmersedRigidBodyTriangleMeshPolygon::reset_reference_configuration(), oomph::ContinuationStorageScheme::set_consistent_pinned_values(), oomph::SolidICProblem::set_newmark_initial_condition_consistently(), oomph::SelfStartingBDF2::shift_time_values(), oomph::IMRBase::shift_time_values(), oomph::Steady< NSTEPS >::shift_time_values(), oomph::NewmarkBDF< NSTEPS >::shift_time_values(), oomph::BDF< NSTEPS >::shift_time_values(), oomph::TR::shift_time_values(), oomph::AxisymQPVDElementWithPressure::solid_p(), oomph::QPVDElementWithPressure< DIM >::solid_p(), oomph::ODEProblem::solution(), SegregatedFSICollapsibleChannelProblem< ELEMENT >::steady_run(), oomph::IMRBase::temporal_error_in_value(), oomph::BDF< NSTEPS >::temporal_error_in_value(), oomph::TR::temporal_error_in_value(), oomph::SelfStartingBDF2::temporal_error_in_value_bdf2(), oomph::TimeStepper::time_derivative(), oomph::IMRODEElement::time_interpolate_u(), oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::u_bar(), SegregatedFSICollapsibleChannelProblem< ELEMENT >::unsteady_run(), and value().
Return i-th value at time level t (t=0: present, t>0: previous) This function is not virtual so that it can be inlined. This means that if we have an explicit pointer to a Data object Data* data_pt->value() always returns to the "raw" stored value.
References i, range_check(), plotPSD::t, and Value.
Compute Vector of values (dofs or pinned) in this data at time level t (t=0: present; t>0: previous).
Compute Vector of values (dofs or pinned) at this node at time level t (t=0: present; t>0: previous)
References i, nvalue(), plotPSD::t, and value().
Compute Vector of values for the Data value.
Compute Vector of values (dofs or pinned) at this Data object.
Return the pointer to the i-the stored value. Typically this is required when direct access to the stored value is required, e.g. when writing functions that return a reference to a variable that is stored in a Data object.
References i, range_check(), and Value.
Referenced by add_value_pt_to_map(), oomph::Problem::arc_length_step_solve(), oomph::ElementWithMovingNodes::assign_all_generic_local_eqn_numbers(), oomph::SpectralElement::assign_all_generic_local_eqn_numbers(), assign_eqn_numbers(), oomph::RefineableElement::assign_hanging_local_eqn_numbers(), oomph::GeneralisedElement::assign_internal_and_external_local_eqn_numbers(), oomph::FiniteElement::assign_nodal_local_eqn_numbers(), oomph::ImmersedRigidBodyElement::centre_rotation_angle(), oomph::ImmersedRigidBodyElement::centre_x_displacement(), oomph::ImmersedRigidBodyElement::centre_y_displacement(), oomph::FpPressureAdvectionDiffusionProblem< ELEMENT >::doc_solution(), ZeroResidualGenElement::eval(), ZeroResidualGenElement::eval_pt(), Global_Physical_Variables::external_pressure(), oomph::SpinePointMarangoniSurfactantFluidInterfaceBoundingElement< ELEMENT >::fill_in_contribution_to_jacobian(), oomph::ConstraintElement::fill_in_contribution_to_residuals(), oomph::GeneralisedElement::fill_in_jacobian_from_external_by_fd(), oomph::ElementWithExternalElement::fill_in_jacobian_from_external_interaction_field_by_fd(), oomph::ElementWithExternalElement::fill_in_jacobian_from_external_interaction_geometric_by_fd(), oomph::GeneralisedElement::fill_in_jacobian_from_internal_by_fd(), oomph::FiniteElement::fill_in_jacobian_from_nodal_by_fd(), oomph::RefineableElement::fill_in_jacobian_from_nodal_by_fd(), oomph::DoubleBuoyantQCrouzeixRaviartElement< DIM >::fill_in_off_diagonal_jacobian_blocks_by_fd(), oomph::AxisymmetricQAdvectionCrouzeixRaviartElement::fill_in_off_diagonal_jacobian_blocks_by_fd(), oomph::RefineableDoubleBuoyantQCrouzeixRaviartElement< DIM >::fill_in_off_diagonal_jacobian_blocks_by_fd(), oomph::PolarCrouzeixRaviartElement::fix_pressure(), oomph::PolarTaylorHoodElement::fix_pressure(), 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::GeneralisedNewtonianNavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::RefineableGeneralisedNewtonianNavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::NavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::RefineableNavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::SpaceTimeNavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::RefineableSpaceTimeNavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), BrethertonElement< ELEMENT >::get_jacobian(), No_Slip::no_slip_condition_first(), No_Slip::no_slip_condition_second(), oomph::PolarCrouzeixRaviartElement::p_pnst(), oomph::SpaceTimeNavierStokesEquations< DIM >::re_st(), oomph::SpaceTimeNavierStokesMixedOrderEquations< DIM >::re_st(), oomph::SpaceTimeNavierStokesEquations< DIM >::re_st_pt(), oomph::SpaceTimeNavierStokesMixedOrderEquations< DIM >::re_st_pt(), oomph::LinearisedNavierStokesEquations::set_eigenfunction_normalisation_element(), oomph::SolidICProblem::set_newmark_initial_condition_consistently(), oomph::SpaceTimeNavierStokesEquations< DIM >::st(), oomph::SpaceTimeNavierStokesEquations< DIM >::st_pt(), oomph::YoungLaplaceContactAngleElement< ELEMENT >::u(), and ZeroResidualGenElement::ZeroResidualGenElement().
Return the pointer to the i-th stored value, or any of its history values (const version). Typically this is required when direct access to the stored value is required, e.g. when writing functions that return a reference to a variable that is stored in a Data object.
References i, range_check(), plotPSD::t, and Value.
|
friend |
|
friend |
|
friend |
Referenced by oomph::SolidNode::set_external_variable_position_pt().
|
friend |
|
friend |
Output operator: output all values at all times, along with any extra information stored for the timestepper.
Data output operator: output equation numbers and values at all times, along with any extra information stored for the timestepper.
|
friend |
|
protected |
C-style array of any Data objects that contain copies of the current Data object's data values.
Referenced by add_copy(), remove_copy(), resize(), set_time_stepper(), and ~Data().
|
staticprotected |
Default (static) timestepper used in steady problems.
Default (steady) timestepper for steady Data.
|
private |
C-style array of pointers to the (global) equation numbers of the values.
Referenced by assign_eqn_numbers(), oomph::HijackedData::clear_copied_pointers(), oomph::CopiedData::clear_copied_pointers(), oomph::CopiedData::CopiedData(), Data(), delete_value_storage(), describe_dofs(), eqn_number(), eqn_number_pt(), oomph::HijackedData::HijackedData(), is_constrained(), is_pinned(), is_segregated_solve_pinned(), oomph::BoundaryNodeBase::make_node_periodic(), oomph::BoundaryNodeBase::make_nodes_periodic(), pin_all(), oomph::HijackedData::reset_copied_pointers(), oomph::CopiedData::reset_copied_pointers(), oomph::BoundaryNode< NODE_TYPE >::reset_copied_pointers(), resize(), self_test(), unpin_all(), oomph::CopiedData::~CopiedData(), and oomph::HijackedData::~HijackedData().
|
static |
Static "Magic number" used in place of the equation number to indicate that the value is constrained because it is associated with non-conforming element boundaries — a hanging node — (and is therefore pinned)
Static "Magic number" to indicate that the value is constrained, usually because is it associated with non-conforming data, otherwise known as hanging nodes
Referenced by constrain(), is_constrained(), and unconstrain().
|
static |
Static "Magic number" to indicate pinned values.
Static "Magic number" used in place of the equation number to indicate that the value is pinned.
Referenced by oomph::ElementWithMovingNodes::assign_all_generic_local_eqn_numbers(), oomph::SpectralElement::assign_all_generic_local_eqn_numbers(), assign_eqn_numbers(), oomph::ElementWithExternalElement::assign_external_interaction_data_local_eqn_numbers(), oomph::RefineableElement::assign_hanging_local_eqn_numbers(), oomph::GeneralisedElement::assign_internal_and_external_local_eqn_numbers(), oomph::FiniteElement::assign_nodal_local_eqn_numbers(), oomph::RefineableSolidElement::assign_solid_hanging_local_eqn_numbers(), oomph::SolidFiniteElement::assign_solid_local_eqn_numbers(), constrain(), is_pinned(), pin(), and pin_all().
|
static |
Static "Magic number" used in place of the equation number to indicate that the value is pinned, but only for the duration of a segregated solve.
Referenced by oomph::MyProblem::check_not_segregated(), is_segregated_solve_pinned(), oomph::MyProblem::segregated_pin_indices(), and oomph::MyProblem::undo_segregated_pinning().
|
static |
Static "Magic number" used in place of the equation number to denote a value that hasn't been classified as pinned or free.
Static "Magic number" to indicate values that haven't been classified as pinned or free
Referenced by oomph::ElementWithMovingNodes::assign_all_generic_local_eqn_numbers(), oomph::SpectralElement::assign_all_generic_local_eqn_numbers(), oomph::GeneralisedElement::assign_internal_and_external_local_eqn_numbers(), oomph::FiniteElement::assign_nodal_local_eqn_numbers(), Data(), resize(), self_test(), unconstrain(), oomph::MyProblem::undo_segregated_pinning(), unpin(), and unpin_all().
|
protected |
Number of Data that contain copies of this Data object's values
Referenced by add_copy(), remove_copy(), resize(), set_time_stepper(), and ~Data().
|
private |
Number of values stored in the data object.
Referenced by assign_eqn_numbers(), constrain_all(), describe_dofs(), nvalue(), pin_all(), range_check(), oomph::CopiedData::reset_copied_pointers(), resize(), self_test(), unconstrain_all(), and unpin_all().
|
private |
Pointer to a Timestepper. The inclusion of a Timestepper pointer in the Data class, ensures that time-derivatives can be calculated and storage can be managed at the low (Data) level.
Referenced by ntstorage(), set_time_stepper(), and time_stepper_pt().
|
private |
C-style array of pointers to data values and possible history values. The data must be ordered in such a way that Value[i][t] gives the i-th data value at the time value t. The ordering is chosen so that all the time levels of a particular value can be access from a single pointer to a double. This is required for copying/hijacking functionality. The data should be accessed by using the member functions value(time,ival) and set_value(time,ival,value), where time=0: present.
Referenced by oomph::HijackedData::clear_copied_pointers(), oomph::CopiedData::clear_copied_pointers(), oomph::CopiedData::CopiedData(), Data(), delete_value_storage(), does_pointer_correspond_to_value(), oomph::HijackedData::HijackedData(), oomph::BoundaryNodeBase::make_node_periodic(), oomph::BoundaryNodeBase::make_nodes_periodic(), oomph::HijackedData::reset_copied_pointers(), oomph::CopiedData::reset_copied_pointers(), oomph::BoundaryNode< NODE_TYPE >::reset_copied_pointers(), resize(), oomph::SolidNode::set_external_variable_position_pt(), oomph::SolidNode::set_position_time_stepper(), set_time_stepper(), set_value(), oomph::SolidNode::SolidNode(), value(), value_pt(), oomph::CopiedData::~CopiedData(), and oomph::HijackedData::~HijackedData().