Example codes
and tutorials
The (Not-so) Quick Guide
List of tutorials/demo codes
Single-Physics Problems
Poisson
Adaptivity illustrated for Poisson
Advection-Diffusion
Unsteady heat equation
Linear wave equation
The Young-Laplace equation
Navier-Stokes
Free-surface Navier-Stokes
Axisymmetric Navier-Stokes
Solid mechanics
Beam structures
Shell structures
Multi-physics Problems
Fluid-structure interaction
Boussinesq convection
Steady thermoelasticity
Methods-based example codes and tutorials
Mesh generation
Linear solvers and preconditioners
Visualisation of the results
Parallel processing
How to write a new element
How to write a new refineable element
Default nonlinear solvers -- the sequence of action functions
...
Documentation
FE theory and top-down discussion of the data structure
The (Not-so) Quick Guide
Comprehensive bottom-up discussion of the data structure
List of available structured and unstructured meshes
Linear solvers and preconditioners
Visualisation of the results
Parallel processing
Coding conventions and C++ style
Creating documentation
Optimisation - robustness vs. "raw speed"
Linear vs. nonlinear problems
Storing shape functions
Changing the default "full" integration scheme
Disabling the ALE formulation of unsteady equations
C vs. C++ output
Different sparse assembly techniques and the STL memory pool
Publications
Publications
Talks
Journal publications
Theses
Picture show
Download
Copyright
Download/installation instructions
Download page
FAQ & Contact
FAQ
Change log
Bugs and other known problems
Completeness of the library & our "To-Do List"
Contact the developers
Get involved

 


Beta release!

Please note that the library has not been "officially" released. While we continue to work on the documentation, these web pages are likely to contain broken links and documents in draft form. Please send an email to

oomph-lib AT maths DOT man DOT ac DOT uk

if you wish to be informed of the library's "official" release.

oomph::GeneralisedElement Class Reference

#include <elements.h>

Inheritance diagram for oomph::GeneralisedElement:

oomph::DisplacementControlElement oomph::FiniteElement oomph::ImposeFluxForWomersleyElement< DIM > oomph::NavierStokesWomersleyPressureControlElement oomph::NetFluxControlElement oomph::PseudoBucklingRingElement oomph::TemplateFreeNavierStokesFluxControlElementBase oomph::AdvectionDiffusionEquations< DIM > oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM > oomph::AxisymmetricNavierStokesEquations oomph::DGElement oomph::ElementWithExternalElement oomph::ElementWithMovingNodes oomph::ElementWithZ2ErrorEstimator oomph::FaceElement oomph::FSIFluidElement oomph::GeneralisedAdvectionDiffusionEquations< DIM > oomph::LinearElasticityEquationsBase< DIM > oomph::LinearWaveEquations< DIM > oomph::PointElement oomph::PoissonEquations< DIM > oomph::QElementBase oomph::QHermiteElement< DIM > oomph::RefineableElement oomph::SolidFiniteElement oomph::SpectralElement oomph::StorableShapeElementBase oomph::TElementBase oomph::UnsteadyHeatEquations< DIM > oomph::WomersleyEquations< DIM > oomph::NetFluxControlElementForWomersleyPressureControl oomph::NavierStokesFluxControlElement< ELEMENT > List of all members.

Detailed Description

A Generalised Element class.

The main components of a GeneralisedElement are:

We also provide interfaces for functions that compute the element's Jacobian matrix and/or the Vector of residuals. In addition, an interface that returns a mass matrix --- the matrix of terms that multiply any time derivatives in the problem --- is also provided to permit explicit time-stepping and the solution of the generalised eigenproblems. Object size: 28 bytes (32-bit), 56 bytes (64-bit)

Definition at line 73 of file elements.h.


Public Member Functions

 GeneralisedElement ()
 Constructor: Initialise all pointers and all values to zero.
virtual ~GeneralisedElement ()
 Virtual destructor to clean up any memory allocated by the object.
 GeneralisedElement (const GeneralisedElement &)
 Broken copy constructor.
void operator= (const GeneralisedElement &)
 Broken assignment operator.
Time *& time_pt ()
 Retun the pointer to the global time.
Time *const & time_pt () const
 Return the pointer to the global time (const version).
double time () const
 Return the global time, accessed via the time pointer.
Data *& internal_data_pt (const unsigned &i)
 Return a pointer to i-th internal data object.
Data *const & internal_data_pt (const unsigned &i) const
 Return a pointer to i-th internal data object (const version).
Data *& external_data_pt (const unsigned &i)
 Return a pointer to i-th external data object.
Data *const & external_data_pt (const unsigned &i) const
 Return a pointer to i-th external data object (const version).
unsigned long eqn_number (const unsigned &ieqn_local)
 Return the global equation number corresponding to the ieqn_local-th local equation number.
int local_eqn_number (const unsigned long &ieqn_global)
 Return the local equation number corresponding to the ieqn_global-th global equation number. Returns minus one (-1) if there is no local degree of freedom corresponding to the chosen global equation number.
unsigned add_external_data (Data *const &data_pt, const bool &fd=true)
bool external_data_fd (const unsigned &i) const
 Return the status of the boolean flag indicating whether the external data is included in the finite difference loop.
void exclude_external_data_fd (const unsigned &i)
 Set the boolean flag to exclude the external datum from the the finite difference loop when computing the jacobian matrix.
void include_external_data_fd (const unsigned &i)
 Set the boolean flag to include the external datum in the the finite difference loop when computing the jacobian matrix.
void flush_external_data ()
 Flush all external data.
void flush_external_data (Data *const &data_pt)
 Flush the object addressed by data_pt from the external data array.
unsigned ninternal_data () const
 Return the number of internal data objects.
unsigned nexternal_data () const
 Return the number of external data objects.
unsigned ndof () const
 Return the number of equations/dofs in the element.
void assign_internal_eqn_numbers (unsigned long &global_number, Vector< double * > &Dof_pt)
 Assign the global equation numbers to the internal Data. The arguments are the current highest global equation number (which will be incremented) and a Vector of pointers to the global variables (to which any unpinned values in the internal Data are added).
virtual void assign_local_eqn_numbers ()
 Setup the arrays of local equation numbers for the element.
virtual void complete_setup_of_dependencies ()
 Complete the setup of any additional dependencies that the element may have. Empty virtual function that may be overloaded for specific derived elements. Used, e.g., for elements with algebraic node update functions to determine the "geometric Data", i.e. the Data that affects the element's shape. This function is called (for all elements) at the very beginning of the equation numbering procedure to ensure that all dependencies are accounted for.
virtual void get_residuals (Vector< double > &residuals)
 Calculate the vector of residuals of the equations in the element. By default initialise the vector to zero and then call the fill_in_contribution_to_residuals() function. Note that this entire function can be overloaded if desired.
virtual void get_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 Calculate the elemental Jacobian matrix "d equation / d variable".
virtual void get_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix)
 Calculate the residuals and the elemental "mass" matrix, the matrix that multiplies the time derivative terms in a problem.
virtual void get_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 Calculate the residuals and jacobian and elemental "mass" matrix, the matrix that multiplies the time derivative terms.
virtual unsigned self_test ()
 Self-test: Have all internal values been classified as pinned/unpinned? Return 0 if OK.
virtual unsigned ndof_types ()
 The number of types of degrees of freedom in this element are sub-divided into.
virtual void get_dof_numbers_for_unknowns (std::list< std::pair< unsigned long, unsigned > > &block_lookup_list)
 Create a list of pairs for the unknowns that this element is "in charge of" -- ignore any unknowns associated with external Data. The first entry in each pair must contain the global equation number of the unknown, while the second one contains the number of the DOF that this unknown is associated with. (The function can obviously only be called if the equation numbering scheme has been set up.).

Static Public Attributes

static bool Suppress_warning_about_repeated_internal_data = false
 Static boolean to suppress warnings about repeated internal data. Defaults to false.
static bool Suppress_warning_about_repeated_external_data = true
 Static boolean to suppress warnings about repeated external data. Defaults to true.
static double Default_fd_jacobian_step = 1.0e-8
 Double used for the default finite difference step in elemental jacobian calculations.

Protected Member Functions

unsigned add_internal_data (Data *const &data_pt, const bool &fd=true)
 Add a (pointer to an) internal data object to the element and return the index required to obtain it from the access function internal_data_pt(). The boolean indicates whether the datum should be included in the general finite-difference loop when calculating the jacobian. The default value is true, i.e. the data will be included in the finite differencing.
bool internal_data_fd (const unsigned &i) const
 Return the status of the boolean flag indicating whether the internal data is included in the finite difference loop.
void exclude_internal_data_fd (const unsigned &i)
 Set the boolean flag to exclude the internal datum from the finite difference loop when computing the jacobian matrix.
void include_internal_data_fd (const unsigned &i)
 Set the boolean flag to include the internal datum in the finite difference loop when computing the jacobian matrix.
void clear_global_eqn_numbers ()
 Clear the storage for the global equation numbers.
void add_global_eqn_numbers (std::deque< unsigned long > const &global_eqn_numbers)
 Add the contents of the queue global_eqn_numbers to the local storage for the local-to-global translation scheme. It is essential that the entries in the queue are added IN ORDER i.e. from the front.
virtual void assign_internal_and_external_local_eqn_numbers ()
 Assign the local equation numbers for the internal and external Data This must be called after the global equation numbers have all been assigned. It is virtual so that it can be overloaded by ElementWithExternalElements so that any external data from the external elements in included in the numbering scheme.
virtual void assign_all_generic_local_eqn_numbers ()
 Assign all the local equation numbering schemes that can be applied generically for the element. In most cases, this is the function that will be overloaded by inherited classes. It is required to ensure that assign_additional_local_eqn_numbers() can always be called after ALL other local equation numbering has been performed. The default for the GeneralisedElement is simply to call internal and external local equation numbering.
virtual void assign_additional_local_eqn_numbers ()
 Setup any additional look-up schemes for local equation numbers. Examples of use include using local storage to refer to explicit degrees of freedom. The additional memory cost of such storage may or may not be offset by fast local access.
int internal_local_eqn (const unsigned &i, const unsigned &j)
 Return the local equation number corresponding to the j-th value stored at the i-th internal data.
int external_local_eqn (const unsigned &i, const unsigned &j)
 Return the local equation number corresponding to the j-th value stored at the i-th external data.
virtual void fill_in_contribution_to_residuals (Vector< double > &residuals)
 Add the elemental contribution to the residuals vector. Note that this function will NOT initialise the residuals vector. It must be called after the residuals vector has been initialised to zero.
void fill_in_jacobian_from_internal_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 Calculate the contributions to the jacobian from the internal degrees of freedom using finite differences. This version of the function assumes that the residuals vector has already been calculated. If the boolean argument is true, the finite differencing will be performed for all internal data, irrespective of the information in Data_fd. The default value (false) uses the information in Data_fd to selectively difference only certain data.
void fill_in_jacobian_from_internal_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 Calculate the contributions to the jacobian from the internal degrees of freedom using finite differences. This version computes the residuals vector before calculating the jacobian terms. If the boolean argument is true, the finite differencing will be performed for all internal data, irrespective of the information in Data_fd. The default value (false) uses the information in Data_fd to selectively difference only certain data.
void fill_in_jacobian_from_external_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 Calculate the contributions to the jacobian from the external degrees of freedom using finite differences. This version of the function assumes that the residuals vector has already been calculated. If the boolean argument is true, the finite differencing will be performed for all external data, irrespective of the information in Data_fd. The default value (false) uses the information in Data_fd to selectively difference only certain data.
void fill_in_jacobian_from_external_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 Calculate the contributions to the jacobian from the external degrees of freedom using finite differences. This version computes the residuals vector before calculating the jacobian terms. If the boolean argument is true, the finite differencing will be performed for all internal data, irrespective of the information in Data_fd. The default value (false) uses the information in Data_fd to selectively difference only certain data.
virtual void update_before_internal_fd ()
 Function that is called before the finite differencing of any internal data. This may be overloaded to update any slaved data before finite differencing takes place.
virtual void reset_after_internal_fd ()
 Function that is call after the finite differencing of the internal data. This may be overloaded to reset any slaved variables that may have changed during the finite differencing.
virtual void update_in_internal_fd (const unsigned &i)
 Function called within the finite difference loop for internal data after a change in any values in the i-th internal data object.
virtual void reset_in_internal_fd (const unsigned &i)
 Function called within the finite difference loop for internal data after the values in the i-th external data object are reset. The default behaviour is to call the update function.
virtual void update_before_external_fd ()
 Function that is called before the finite differencing of any external data. This may be overloaded to update any slaved data before finite differencing takes place.
virtual void reset_after_external_fd ()
 Function that is call after the finite differencing of the external data. This may be overloaded to reset any slaved variables that may have changed during the finite differencing.
virtual void update_in_external_fd (const unsigned &i)
 Function called within the finite difference loop for external data after a change in any values in the i-th external data object.
virtual void reset_in_external_fd (const unsigned &i)
 Function called within the finite difference loop for external data after the values in the i-th external data object are reset. The default behaviour is to call the update function.
virtual void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 Add the elemental contribution to the jacobian matrix. and the residuals vector. Note that this function will NOT initialise the residuals vector or the jacobian matrix. It must be called after the residuals vector and jacobian matrix have been initialised to zero. The default is to use finite differences to calculate the jacobian.
virtual void fill_in_contribution_to_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix)
 Add the elemental contribution to the mass matrix and the residuals vector. Note that this function will NOT initialise the residuals vector or the mass matrix. It must be called after the residuals vector and jacobian matrix have been initialised to zero. The default is deliberately broken.
virtual void fill_in_contribution_to_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 Add the elemental contribution to the jacobian matrix, mass matrix and the residuals vector. Note that this function should NOT initialise any entries. It must be called after the residuals vector and matrices have been initialised to zero. The default is deliberately broken.

Static Protected Attributes

static DenseMatrix< double > Dummy_matrix
 Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case when only the residuals are being assembled.

Private Attributes

TimeTime_pt
 Pointer to global time.
unsigned long * Eqn_number
 Storage for the global equation numbers represented by the degrees of freedom in the element.
Data ** Data_pt
 Storage for pointers to internal and external data. The data is of the same type and so can be addressed by a single pointer. The idea is that the array is of total size Ninternal_data + Nexternal_data. The internal data are stored at the beginning of the array and the external data are stored at the end of the array.
int ** Data_local_eqn
 Pointer to array storage for local equation numbers associated with internal and external data. Again, we save storage by using a single pointer to access this information. The first index of the array is of dimension Nineternal_data + Nexternal_data and the second index varies with the number of values stored at the data object. In the most general case, however, the scheme is as memory efficient as possible.
unsigned Ndof
 Number of degrees of freedom.
unsigned Ninternal_data
 Number of internal data.
unsigned Nexternal_data
 Number of external data.
std::vector< bool > Data_fd
 Storage for boolean flags of size Ninternal_data + Nexternal_data that correspond to the data used in the element. The flags indicate whether the particular internal or external data should be included in the general finite-difference loop in fill_in_jacobian_from_internal_by_fd() or fill_in_jacobian_from_external_by_fd(). The default is that all data WILL be included in the finite-difference loop, but in many circumstances it is possible to treat certain (external) data analytically. The use of an STL vector is optimal for memory use because the booleans are represented as single-bits.

Constructor & Destructor Documentation

oomph::GeneralisedElement::GeneralisedElement (  )  [inline]

Constructor: Initialise all pointers and all values to zero.

Definition at line 519 of file elements.h.

oomph::GeneralisedElement::~GeneralisedElement (  )  [virtual]

Virtual destructor to clean up any memory allocated by the object.

Destructor for generalised elements: Wipe internal data. Pointers to external data get NULLed but are not deleted because they are (generally) shared by lots of elements.

Definition at line 209 of file elements.cc.

References Data_local_eqn, Data_pt, Eqn_number, and Ninternal_data.

oomph::GeneralisedElement::GeneralisedElement ( const GeneralisedElement  )  [inline]

Broken copy constructor.

Definition at line 530 of file elements.h.

References oomph::BrokenCopy::broken_copy().


Member Function Documentation

unsigned oomph::GeneralisedElement::add_external_data ( Data *const &  data_pt,
const bool &  fd = true 
)

Add a (pointer to an) external data object to the element and return its index (i.e. the index required to obtain it from the access function external_data_pt(...). The optional boolean flag indicates whether the data should be included in the general finite-difference loop when calculating the jacobian. The default value is true, i.e. the data will be included in the finite-differencing

Definition at line 240 of file elements.cc.

References Data_fd, Data_pt, external_data_pt(), Nexternal_data, Ninternal_data, oomph::oomph_info, and Suppress_warning_about_repeated_external_data.

Referenced by oomph::NavierStokesWomersleyPressureControlElement::add_pressure_data(), oomph::TemplateFreeNavierStokesFluxControlElementBase::add_pressure_data(), oomph::DisplacementControlElement::DisplacementControlElement(), oomph::ElasticSurfaceFluidInterfaceElement< ELEMENT >::make_edge_element(), oomph::SpineSurfaceFluidInterfaceElement< ELEMENT >::make_edge_element(), oomph::ElasticLineFluidInterfaceElement< ELEMENT >::make_edge_element(), oomph::SpineLineFluidInterfaceElement< ELEMENT >::make_edge_element(), oomph::ElasticAxisymmetricFluidInterfaceElement< ELEMENT >::make_edge_element(), oomph::SpineAxisymmetricFluidInterfaceElement< ELEMENT >::make_edge_element(), oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::set_boundary_shape_geom_object_pt(), oomph::NavierStokesImpedanceTractionElement< BULK_NAVIER_STOKES_ELEMENT, WOMERSLEY_ELEMENT, DIM >::set_external_data_from_navier_stokes_outflow_mesh(), oomph::FluidInterfaceElement::set_external_pressure_data(), oomph::WomersleyEquations< DIM >::set_pressure_gradient_and_add_as_external_data(), oomph::PseudoBucklingRingElement::set_reference_pressure_pt(), and oomph::DGFaceElement::setup_neighbour_info().

void oomph::GeneralisedElement::add_global_eqn_numbers ( std::deque< unsigned long > const &  global_eqn_numbers  )  [protected]

Add the contents of the queue global_eqn_numbers to the local storage for the local-to-global translation scheme. It is essential that the entries in the queue are added IN ORDER i.e. from the front.

Add the contents of the queue global_eqn_numbers to the local storage for the equation numbers, which represents the local-to-global translation scheme. It is essential that the entries are added in order, i.e. from the front.

Definition at line 143 of file elements.cc.

References Eqn_number, and Ndof.

Referenced by oomph::SpectralElement::assign_all_generic_local_eqn_numbers(), oomph::ElementWithMovingNodes::assign_all_generic_local_eqn_numbers(), oomph::ElementWithExternalElement::assign_external_interaction_data_local_eqn_numbers(), oomph::RefineableElement::assign_hanging_local_eqn_numbers(), assign_internal_and_external_local_eqn_numbers(), oomph::FiniteElement::assign_nodal_local_eqn_numbers(), oomph::RefineableSolidElement::assign_solid_hanging_local_eqn_numbers(), and oomph::SolidFiniteElement::assign_solid_local_eqn_numbers().

unsigned oomph::GeneralisedElement::add_internal_data ( Data *const &  data_pt,
const bool &  fd = true 
) [protected]

Add a (pointer to an) internal data object to the element and return the index required to obtain it from the access function internal_data_pt(). The boolean indicates whether the datum should be included in the general finite-difference loop when calculating the jacobian. The default value is true, i.e. the data will be included in the finite differencing.

Add a (pointer to an) internal data object to the element and return the index required to obtain it from the access function internal_data_pt()

Definition at line 61 of file elements.cc.

References Data_fd, Data_pt, internal_data_pt(), Nexternal_data, Ninternal_data, oomph::oomph_info, and Suppress_warning_about_repeated_internal_data.

Referenced by oomph::AxisymmetricQCrouzeixRaviartElement::AxisymmetricQCrouzeixRaviartElement(), oomph::AxisymQPVDElementWithPressure::AxisymQPVDElementWithPressure(), oomph::DisplacementControlElement::DisplacementControlElement(), oomph::ImposeFluxForWomersleyElement< DIM >::ImposeFluxForWomersleyElement(), oomph::NavierStokesWomersleyPressureControlElement::NavierStokesWomersleyPressureControlElement(), oomph::NetFluxControlElement::NetFluxControlElement(), oomph::PseudoBucklingRingElement::PseudoBucklingRingElement(), oomph::QCrouzeixRaviartElement< DIM >::QCrouzeixRaviartElement(), and oomph::QPVDElementWithPressure< DIM >::QPVDElementWithPressure().

virtual void oomph::GeneralisedElement::assign_additional_local_eqn_numbers (  )  [inline, protected, virtual]

Setup any additional look-up schemes for local equation numbers. Examples of use include using local storage to refer to explicit degrees of freedom. The additional memory cost of such storage may or may not be offset by fast local access.

Reimplemented in oomph::DisplacementControlElement.

Definition at line 244 of file elements.h.

Referenced by assign_local_eqn_numbers().

virtual void oomph::GeneralisedElement::assign_all_generic_local_eqn_numbers (  )  [inline, protected, virtual]

Assign all the local equation numbering schemes that can be applied generically for the element. In most cases, this is the function that will be overloaded by inherited classes. It is required to ensure that assign_additional_local_eqn_numbers() can always be called after ALL other local equation numbering has been performed. The default for the GeneralisedElement is simply to call internal and external local equation numbering.

Reimplemented in oomph::ElementWithMovingNodes, oomph::ElementWithSpecificMovingNodes< ELEMENT, NODE_TYPE >, oomph::FiniteElement, oomph::SolidFiniteElement, oomph::FaceElementAsGeomObject< ELEMENT >, oomph::SpectralElement, oomph::RefineableQSpectralPoissonElement< DIM, NNODE_1D >, oomph::ElementWithSpecificMovingNodes< oomph::FaceGeometry< ELEMENT >, oomph::SpineNode >, oomph::ElementWithSpecificMovingNodes< ELEMENT, oomph::SpineNode >, oomph::ElementWithSpecificMovingNodes< ELEMENT, oomph::AlgebraicNode >, oomph::ElementWithSpecificMovingNodes< oomph::FaceGeometry< oomph::FaceGeometry< ELEMENT > >, oomph::SpineNode >, and oomph::ElementWithSpecificMovingNodes< ELEMENT, oomph::MacroElementNodeUpdateNode >.

Definition at line 235 of file elements.h.

References assign_internal_and_external_local_eqn_numbers().

Referenced by oomph::FiniteElement::assign_all_generic_local_eqn_numbers(), and assign_local_eqn_numbers().

void oomph::GeneralisedElement::assign_internal_and_external_local_eqn_numbers (  )  [protected, virtual]

Assign the local equation numbers for the internal and external Data This must be called after the global equation numbers have all been assigned. It is virtual so that it can be overloaded by ElementWithExternalElements so that any external data from the external elements in included in the numbering scheme.

This function loops over the internal and external data of the element, adds the GLOBAL equation numbers to the local-to-global look-up scheme and fills in the look-up schemes for the local equation numbers

Reimplemented in oomph::ElementWithExternalElement.

Definition at line 528 of file elements.cc.

References add_global_eqn_numbers(), Data_local_eqn, Data_pt, eqn_number(), external_data_pt(), internal_data_pt(), oomph::Data::Is_pinned, oomph::Data::Is_unclassified, local_eqn_number(), ndof(), Nexternal_data, Ninternal_data, and oomph::Data::nvalue().

Referenced by assign_all_generic_local_eqn_numbers(), and oomph::ElementWithExternalElement::assign_internal_and_external_local_eqn_numbers().

void oomph::GeneralisedElement::assign_internal_eqn_numbers ( unsigned long &  global_number,
Vector< double * > &  Dof_pt 
)

Assign the global equation numbers to the internal Data. The arguments are the current highest global equation number (which will be incremented) and a Vector of pointers to the global variables (to which any unpinned values in the internal Data are added).

This function loops over the internal data of the element and assigns GLOBAL equation numbers to the data objects.

Pass:

Definition at line 447 of file elements.cc.

References internal_data_pt(), and Ninternal_data.

void oomph::GeneralisedElement::assign_local_eqn_numbers (  )  [virtual]

Setup the arrays of local equation numbers for the element.

Setup the arrays of local equation numbers for the element. In general, this function should not need to be overloaded. Instead the function assign_all_generic_local_eqn_numbers() will be overloaded by different types of Element. That said, the function is virtual so that it may be overloaded by the user if they *really* know what they are doing.

Reimplemented in oomph::Hijacked< oomph::SpineElement< oomph::FaceGeometry< ELEMENT > > >, and oomph::Hijacked< oomph::SpineElement< oomph::FaceGeometry< oomph::FaceGeometry< ELEMENT > > > >.

Definition at line 464 of file elements.cc.

References assign_additional_local_eqn_numbers(), assign_all_generic_local_eqn_numbers(), clear_global_eqn_numbers(), oomph::FiniteElement::dim(), Eqn_number, Ndof, oomph::FiniteElement::nnode(), oomph::FiniteElement::node_pt(), OOMPH_EXCEPTION_LOCATION, and oomph::Node::x().

void oomph::GeneralisedElement::clear_global_eqn_numbers (  )  [inline, protected]

Clear the storage for the global equation numbers.

Definition at line 205 of file elements.h.

References Eqn_number, and Ndof.

Referenced by assign_local_eqn_numbers().

virtual void oomph::GeneralisedElement::complete_setup_of_dependencies (  )  [inline, virtual]

Complete the setup of any additional dependencies that the element may have. Empty virtual function that may be overloaded for specific derived elements. Used, e.g., for elements with algebraic node update functions to determine the "geometric Data", i.e. the Data that affects the element's shape. This function is called (for all elements) at the very beginning of the equation numbering procedure to ensure that all dependencies are accounted for.

Reimplemented in oomph::ElementWithMovingNodes, oomph::ElementWithSpecificMovingNodes< ELEMENT, NODE_TYPE >, oomph::SpineElement< ELEMENT >, oomph::ElementWithSpecificMovingNodes< oomph::FaceGeometry< ELEMENT >, oomph::SpineNode >, oomph::ElementWithSpecificMovingNodes< ELEMENT, oomph::SpineNode >, oomph::ElementWithSpecificMovingNodes< ELEMENT, oomph::AlgebraicNode >, oomph::ElementWithSpecificMovingNodes< oomph::FaceGeometry< oomph::FaceGeometry< ELEMENT > >, oomph::SpineNode >, oomph::ElementWithSpecificMovingNodes< ELEMENT, oomph::MacroElementNodeUpdateNode >, oomph::SpineElement< oomph::FaceGeometry< ELEMENT > >, and oomph::SpineElement< oomph::FaceGeometry< oomph::FaceGeometry< ELEMENT > > >.

Definition at line 784 of file elements.h.

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

unsigned long oomph::GeneralisedElement::eqn_number ( const unsigned &  ieqn_local  )  [inline]

Return the global equation number corresponding to the ieqn_local-th local equation number.

Definition at line 638 of file elements.h.

References Eqn_number, Ndof, and OOMPH_EXCEPTION_LOCATION.

Referenced by oomph::SpectralElement::assign_all_generic_local_eqn_numbers(), oomph::ElementWithMovingNodes::assign_all_generic_local_eqn_numbers(), oomph::ElementWithExternalElement::assign_external_interaction_data_local_eqn_numbers(), oomph::RefineableElement::assign_hanging_local_eqn_numbers(), 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(), oomph::RefineableAdvectionDiffusionEquations< DIM >::dinterpolated_u_adv_diff_ddata(), oomph::RefineableNavierStokesEquations< DIM >::dinterpolated_u_nst_ddata(), oomph::HopfHandler::eqn_number(), oomph::PitchForkHandler::eqn_number(), oomph::FoldHandler::eqn_number(), oomph::EigenProblemHandler::eqn_number(), oomph::ExplicitTimeStepHandler::eqn_number(), oomph::AssemblyHandler::eqn_number(), oomph::NetFluxControlElementForWomersleyPressureControl::get_dof_numbers_for_unknowns(), oomph::NavierStokesWomersleyPressureControlElement::get_dof_numbers_for_unknowns(), oomph::TTaylorHoodElement< DIM >::get_dof_numbers_for_unknowns(), oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::get_dof_numbers_for_unknowns(), oomph::NetFluxControlElement::get_dof_numbers_for_unknowns(), oomph::QTaylorHoodElement< DIM >::get_dof_numbers_for_unknowns(), oomph::QCrouzeixRaviartElement< DIM >::get_dof_numbers_for_unknowns(), oomph::DisplacementControlElement::get_dof_numbers_for_unknowns(), oomph::HopfHandler::get_jacobian(), oomph::PitchForkHandler::get_jacobian(), oomph::FoldHandler::get_jacobian(), oomph::HopfHandler::get_residuals(), oomph::PitchForkHandler::get_residuals(), and oomph::FoldHandler::get_residuals().

void oomph::GeneralisedElement::exclude_external_data_fd ( const unsigned &  i  )  [inline]

Set the boolean flag to exclude the external datum from the the finite difference loop when computing the jacobian matrix.

Definition at line 710 of file elements.h.

References Data_fd, Nexternal_data, Ninternal_data, and OOMPH_EXCEPTION_LOCATION.

void oomph::GeneralisedElement::exclude_internal_data_fd ( const unsigned &  i  )  [inline, protected]

Set the boolean flag to exclude the internal datum from the finite difference loop when computing the jacobian matrix.

Definition at line 166 of file elements.h.

References Data_fd, Ninternal_data, and OOMPH_EXCEPTION_LOCATION.

bool oomph::GeneralisedElement::external_data_fd ( const unsigned &  i  )  const [inline]

Return the status of the boolean flag indicating whether the external data is included in the finite difference loop.

Definition at line 688 of file elements.h.

References Data_fd, Nexternal_data, Ninternal_data, and OOMPH_EXCEPTION_LOCATION.

Referenced by fill_in_jacobian_from_external_by_fd().

Data* const& oomph::GeneralisedElement::external_data_pt ( const unsigned &  i  )  const [inline]

Return a pointer to i-th external data object (const version).

Definition at line 610 of file elements.h.

References Data_pt, Nexternal_data, Ninternal_data, and OOMPH_EXCEPTION_LOCATION.

Data*& oomph::GeneralisedElement::external_data_pt ( const unsigned &  i  )  [inline]

Return a pointer to i-th external data object.

Definition at line 591 of file elements.h.

References Data_pt, Nexternal_data, Ninternal_data, and OOMPH_EXCEPTION_LOCATION.

Referenced by add_external_data(), oomph::ElementWithExternalElement::assign_external_interaction_data_local_eqn_numbers(), assign_internal_and_external_local_eqn_numbers(), external_local_eqn(), oomph::NavierStokesFluxControlElement< ELEMENT >::fill_in_generic_residual_contribution_fluid_traction(), oomph::NavierStokesWomersleyPressureControlElement::fill_in_generic_residual_contribution_pressure_control(), fill_in_jacobian_from_external_by_fd(), flush_external_data(), oomph::PseudoBucklingRingElement::reference_pressure_pt(), oomph::RefineableNavierStokesFluxControlElement< ELEMENT >::refineable_fill_in_generic_residual_contribution_fluid_traction(), and self_test().

int oomph::GeneralisedElement::external_local_eqn ( const unsigned &  i,
const unsigned &  j 
) [inline, protected]

Return the local equation number corresponding to the j-th value stored at the i-th external data.

Definition at line 293 of file elements.h.

References Data_local_eqn, external_data_pt(), Nexternal_data, Ninternal_data, oomph::Data::nvalue(), and OOMPH_EXCEPTION_LOCATION.

Referenced by oomph::DGFaceElement::add_flux_contributions(), oomph::DisplacementControlElement::assign_additional_local_eqn_numbers(), oomph::NavierStokesWomersleyPressureControlElement::fill_in_generic_residual_contribution_pressure_control(), oomph::WomersleyEquations< DIM >::fill_in_generic_residual_contribution_womersley(), fill_in_jacobian_from_external_by_fd(), oomph::FluidInterfaceElement::pext_local_eqn(), oomph::PseudoBucklingRingElement::reference_pressure_local_eqn(), and oomph::RefineableNavierStokesFluxControlElement< ELEMENT >::refineable_fill_in_generic_residual_contribution_fluid_traction().

virtual void oomph::GeneralisedElement::fill_in_contribution_to_jacobian ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian 
) [inline, protected, virtual]

Add the elemental contribution to the jacobian matrix. and the residuals vector. Note that this function will NOT initialise the residuals vector or the jacobian matrix. It must be called after the residuals vector and jacobian matrix have been initialised to zero. The default is to use finite differences to calculate the jacobian.

Reimplemented in oomph::NavierStokesTractionElement< ELEMENT >, oomph::ImposeParallelOutflowElement< ELEMENT >, oomph::NavierStokesEquations< DIM >, oomph::NetFluxControlElement, oomph::NavierStokesFluxControlElement< ELEMENT >, oomph::RefineableNavierStokesFluxControlElement< ELEMENT >, oomph::ElementWithExternalElement, oomph::FiniteElement, oomph::SolidFiniteElement, oomph::FaceElementAsGeomObject< ELEMENT >, oomph::FSIWallElement, oomph::PoissonEquations< DIM >, oomph::PoissonFluxElement< ELEMENT >, oomph::AxisymmetricNavierStokesEquations, oomph::SpineAxisymmetricFluidInterfaceElement< ELEMENT >, oomph::ElasticAxisymmetricFluidInterfaceElement< ELEMENT >, oomph::SpinePointFluidInterfaceEdgeElement< ELEMENT >, oomph::ElasticPointFluidInterfaceEdgeElement< ELEMENT >, oomph::SpineLineFluidInterfaceEdgeElement< ELEMENT >, oomph::ElasticLineFluidInterfaceEdgeElement< ELEMENT >, oomph::SpineLineFluidInterfaceElement< ELEMENT >, oomph::ElasticLineFluidInterfaceElement< ELEMENT >, oomph::SpineSurfaceFluidInterfaceElement< ELEMENT >, oomph::ElasticSurfaceFluidInterfaceElement< ELEMENT >, oomph::PVDEquations< DIM >, oomph::PVDEquationsWithPressure< DIM >, oomph::SolidTractionElement< ELEMENT >, oomph::FSISolidTractionElement< ELEMENT, DIM >, oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >, oomph::AxisymmetricPVDEquations, oomph::AxisymmetricPVDEquationsWithPressure, oomph::AxisymmetricSolidTractionElement< ELEMENT >, oomph::KirchhoffLoveBeamEquations, oomph::FSIHermiteBeamElement, oomph::KirchhoffLoveShellEquations, oomph::FSIDiagHermiteShellElement, oomph::AdvectionDiffusionEquations< DIM >, oomph::AdvectionDiffusionFluxElement< ELEMENT >, oomph::GeneralisedAdvectionDiffusionEquations< DIM >, oomph::UnsteadyHeatEquations< DIM >, oomph::UnsteadyHeatFluxElement< ELEMENT >, oomph::LinearWaveEquations< DIM >, oomph::LinearWaveFluxElement< ELEMENT >, oomph::WomersleyEquations< DIM >, oomph::NavierStokesImpedanceTractionElement< BULK_NAVIER_STOKES_ELEMENT, WOMERSLEY_ELEMENT, DIM >, oomph::NavierStokesWomersleyPressureControlElement, oomph::LinearElasticityEquations< DIM >, oomph::LinearElasticityTractionElement< ELEMENT >, and oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >.

Definition at line 481 of file elements.h.

References fill_in_contribution_to_residuals(), fill_in_jacobian_from_external_by_fd(), fill_in_jacobian_from_internal_by_fd(), get_residuals(), and Ndof.

Referenced by get_jacobian().

void oomph::GeneralisedElement::fill_in_contribution_to_jacobian_and_mass_matrix ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian,
DenseMatrix< double > &  mass_matrix 
) [protected, virtual]

Add the elemental contribution to the jacobian matrix, mass matrix and the residuals vector. Note that this function should NOT initialise any entries. It must be called after the residuals vector and matrices have been initialised to zero. The default is deliberately broken.

Reimplemented in oomph::NavierStokesEquations< DIM >, oomph::AxisymmetricNavierStokesEquations, oomph::AdvectionDiffusionEquations< DIM >, oomph::GeneralisedAdvectionDiffusionEquations< DIM >, and oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >.

Definition at line 881 of file elements.cc.

References OOMPH_EXCEPTION_LOCATION.

Referenced by get_jacobian_and_mass_matrix().

void oomph::GeneralisedElement::fill_in_contribution_to_mass_matrix ( Vector< double > &  residuals,
DenseMatrix< double > &  mass_matrix 
) [protected, virtual]

Add the elemental contribution to the mass matrix and the residuals vector. Note that this function will NOT initialise the residuals vector or the mass matrix. It must be called after the residuals vector and jacobian matrix have been initialised to zero. The default is deliberately broken.

Definition at line 848 of file elements.cc.

References OOMPH_EXCEPTION_LOCATION.

Referenced by oomph::DGElement::get_inverse_mass_matrix_times_residuals(), get_mass_matrix(), and oomph::DGElement::pre_compute_mass_matrix().

virtual void oomph::GeneralisedElement::fill_in_contribution_to_residuals ( Vector< double > &  residuals  )  [inline, protected, virtual]

Add the elemental contribution to the residuals vector. Note that this function will NOT initialise the residuals vector. It must be called after the residuals vector has been initialised to zero.

Reimplemented in oomph::NavierStokesTractionElement< ELEMENT >, oomph::ImposeParallelOutflowElement< ELEMENT >, oomph::NavierStokesEquations< DIM >, oomph::NetFluxControlElement, oomph::NavierStokesFluxControlElement< ELEMENT >, oomph::RefineableNavierStokesFluxControlElement< ELEMENT >, oomph::DisplacementControlElement, oomph::PoissonEquations< DIM >, oomph::PoissonFluxElement< ELEMENT >, oomph::AxisymmetricNavierStokesEquations, oomph::FluidInterfaceEdgeElement, oomph::FluidInterfaceElement, oomph::PVDEquations< DIM >, oomph::PVDEquationsWithPressure< DIM >, oomph::SolidTractionElement< ELEMENT >, oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >, oomph::AxisymmetricPVDEquations, oomph::AxisymmetricPVDEquationsWithPressure, oomph::AxisymmetricSolidTractionElement< ELEMENT >, oomph::KirchhoffLoveBeamEquations, oomph::ClampedSlidingHermiteBeamBoundaryConditionElement, oomph::KirchhoffLoveShellEquations, oomph::ClampedHermiteShellBoundaryConditionElement, oomph::AdvectionDiffusionEquations< DIM >, oomph::AdvectionDiffusionFluxElement< ELEMENT >, oomph::GeneralisedAdvectionDiffusionEquations< DIM >, oomph::UnsteadyHeatEquations< DIM >, oomph::UnsteadyHeatFluxElement< ELEMENT >, oomph::LinearWaveEquations< DIM >, oomph::LinearWaveFluxElement< ELEMENT >, oomph::WomersleyEquations< DIM >, oomph::NavierStokesImpedanceTractionElement< BULK_NAVIER_STOKES_ELEMENT, WOMERSLEY_ELEMENT, DIM >, oomph::NavierStokesWomersleyPressureControlElement, oomph::LinearElasticityEquations< DIM >, oomph::LinearElasticityTractionElement< ELEMENT >, and oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >.

Definition at line 340 of file elements.h.

References OOMPH_EXCEPTION_LOCATION.

Referenced by oomph::FSIWallElement::fill_in_contribution_to_jacobian(), oomph::SolidFiniteElement::fill_in_contribution_to_jacobian(), oomph::FiniteElement::fill_in_contribution_to_jacobian(), fill_in_contribution_to_jacobian(), oomph::ElementWithExternalElement::fill_in_contribution_to_jacobian(), and get_residuals().

void oomph::GeneralisedElement::fill_in_jacobian_from_external_by_fd ( DenseMatrix< double > &  jacobian,
const bool &  fd_all_data = false 
) [inline, protected]

Calculate the contributions to the jacobian from the external degrees of freedom using finite differences. This version computes the residuals vector before calculating the jacobian terms. If the boolean argument is true, the finite differencing will be performed for all internal data, irrespective of the information in Data_fd. The default value (false) uses the information in Data_fd to selectively difference only certain data.

Definition at line 421 of file elements.h.

References fill_in_jacobian_from_external_by_fd(), get_residuals(), and Ndof.

void oomph::GeneralisedElement::fill_in_jacobian_from_external_by_fd ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian,
const bool &  fd_all_data = false 
) [protected]

Calculate the contributions to the jacobian from the external degrees of freedom using finite differences. This version of the function assumes that the residuals vector has already been calculated. If the boolean argument is true, the finite differencing will be performed for all external data, irrespective of the information in Data_fd. The default value (false) uses the information in Data_fd to selectively difference only certain data.

This function calculates the entries of Jacobian matrix, used in the Newton method, associated with the external degrees of freedom. It does this using finite differences, rather than an analytical formulation, so can be done in total generality. If the boolean argument is true, the finite differencing will be performed for all internal data, irrespective of the information in Data_fd. The default value (false) uses the information in Data_fd to selectively difference only certain data.

Definition at line 758 of file elements.cc.

References Default_fd_jacobian_step, external_data_fd(), external_data_pt(), external_local_eqn(), get_residuals(), ndof(), Nexternal_data, reset_after_external_fd(), reset_in_external_fd(), update_before_external_fd(), and update_in_external_fd().

Referenced by oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::fill_in_contribution_to_jacobian(), oomph::KirchhoffLoveShellEquations::fill_in_contribution_to_jacobian(), oomph::ElasticLineFluidInterfaceEdgeElement< ELEMENT >::fill_in_contribution_to_jacobian(), oomph::SpineLineFluidInterfaceEdgeElement< ELEMENT >::fill_in_contribution_to_jacobian(), oomph::ElasticPointFluidInterfaceEdgeElement< ELEMENT >::fill_in_contribution_to_jacobian(), oomph::SpinePointFluidInterfaceEdgeElement< ELEMENT >::fill_in_contribution_to_jacobian(), oomph::FSIWallElement::fill_in_contribution_to_jacobian(), oomph::SolidFiniteElement::fill_in_contribution_to_jacobian(), oomph::FiniteElement::fill_in_contribution_to_jacobian(), fill_in_contribution_to_jacobian(), oomph::ElementWithExternalElement::fill_in_contribution_to_jacobian(), oomph::KirchhoffLoveBeamEquations::fill_in_contribution_to_jacobian(), and fill_in_jacobian_from_external_by_fd().

void oomph::GeneralisedElement::fill_in_jacobian_from_internal_by_fd ( DenseMatrix< double > &  jacobian,
const bool &  fd_all_data = false 
) [inline, protected]

Calculate the contributions to the jacobian from the internal degrees of freedom using finite differences. This version computes the residuals vector before calculating the jacobian terms. If the boolean argument is true, the finite differencing will be performed for all internal data, irrespective of the information in Data_fd. The default value (false) uses the information in Data_fd to selectively difference only certain data.

Definition at line 387 of file elements.h.

References fill_in_jacobian_from_internal_by_fd(), get_residuals(), and Ndof.

void oomph::GeneralisedElement::fill_in_jacobian_from_internal_by_fd ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian,
const bool &  fd_all_data = false 
) [protected]

Calculate the contributions to the jacobian from the internal degrees of freedom using finite differences. This version of the function assumes that the residuals vector has already been calculated. If the boolean argument is true, the finite differencing will be performed for all internal data, irrespective of the information in Data_fd. The default value (false) uses the information in Data_fd to selectively difference only certain data.

This function calculates the entries of Jacobian matrix, used in the Newton method, associated with the internal degrees of freedom. It does this using finite differences, rather than an analytical formulation, so can be done in total generality. If the boolean argument is true, the finite differencing will be performed for all internal data, irrespective of the information in Data_fd. The default value (false) uses the information in Data_fd to selectively difference only certain data.

Definition at line 664 of file elements.cc.

References Default_fd_jacobian_step, get_residuals(), internal_data_fd(), internal_data_pt(), internal_local_eqn(), ndof(), Ninternal_data, reset_after_internal_fd(), reset_in_internal_fd(), update_before_internal_fd(), and update_in_internal_fd().

Referenced by oomph::FSIWallElement::fill_in_contribution_to_jacobian(), oomph::SolidFiniteElement::fill_in_contribution_to_jacobian(), oomph::FiniteElement::fill_in_contribution_to_jacobian(), fill_in_contribution_to_jacobian(), oomph::ElementWithExternalElement::fill_in_contribution_to_jacobian(), and fill_in_jacobian_from_internal_by_fd().

void oomph::GeneralisedElement::flush_external_data ( Data *const &  data_pt  ) 

Flush the object addressed by data_pt from the external data array.

Remove the object addressed by data_pt from the external data array Note that this could mess up the numbering of other external data

Definition at line 352 of file elements.cc.

References Data_fd, Data_pt, external_data_pt(), Nexternal_data, Ninternal_data, and OOMPH_EXCEPTION_LOCATION.

void oomph::GeneralisedElement::flush_external_data (  ) 

Flush all external data.

Flush all the external data, i.e. clear the pointers to external data from the internal storage.

Definition at line 312 of file elements.cc.

References Data_fd, Data_pt, Nexternal_data, and Ninternal_data.

Referenced by oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::set_boundary_shape_geom_object_pt(), and oomph::PseudoBucklingRingElement::set_reference_pressure_pt().

virtual void oomph::GeneralisedElement::get_dof_numbers_for_unknowns ( std::list< std::pair< unsigned long, unsigned > > &  block_lookup_list  )  [inline, virtual]

Create a list of pairs for the unknowns that this element is "in charge of" -- ignore any unknowns associated with external Data. The first entry in each pair must contain the global equation number of the unknown, while the second one contains the number of the DOF that this unknown is associated with. (The function can obviously only be called if the equation numbering scheme has been set up.).

Reimplemented in oomph::QCrouzeixRaviartElement< DIM >, oomph::QTaylorHoodElement< DIM >, oomph::NetFluxControlElement, oomph::TTaylorHoodElement< DIM >, oomph::DisplacementControlElement, oomph::PVDEquationsBase< DIM >, oomph::FSISolidTractionElement< ELEMENT, DIM >, oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >, oomph::FSIHermiteBeamElement, oomph::FSIDiagHermiteShellElement, oomph::ClampedHermiteShellBoundaryConditionElement, oomph::NavierStokesWomersleyPressureControlElement, oomph::NetFluxControlElementForWomersleyPressureControl, and oomph::LinearElasticityEquationsBase< DIM >.

Definition at line 879 of file elements.h.

References OOMPH_EXCEPTION_LOCATION.

virtual void oomph::GeneralisedElement::get_jacobian ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian 
) [inline, virtual]

Calculate the elemental Jacobian matrix "d equation / d variable".

Reimplemented in oomph::ElementWithSpecificMovingNodes< ELEMENT, NODE_TYPE >, oomph::PseudoBucklingRingElement, oomph::RefineableSolidQElement< 3 >, oomph::RefineableSolidQElement< 2 >, oomph::ImposeFluxForWomersleyElement< DIM >, oomph::ElementWithSpecificMovingNodes< oomph::FaceGeometry< ELEMENT >, oomph::SpineNode >, oomph::ElementWithSpecificMovingNodes< ELEMENT, oomph::SpineNode >, oomph::ElementWithSpecificMovingNodes< ELEMENT, oomph::AlgebraicNode >, oomph::ElementWithSpecificMovingNodes< oomph::FaceGeometry< oomph::FaceGeometry< ELEMENT > >, oomph::SpineNode >, oomph::ElementWithSpecificMovingNodes< ELEMENT, oomph::MacroElementNodeUpdateNode >, oomph::Hijacked< oomph::SpineElement< oomph::FaceGeometry< ELEMENT > > >, and oomph::Hijacked< oomph::SpineElement< oomph::FaceGeometry< oomph::FaceGeometry< ELEMENT > > > >.

Definition at line 799 of file elements.h.

References fill_in_contribution_to_jacobian(), oomph::DenseMatrix< T >::initialise(), and oomph::Vector< _Tp >::initialise().

Referenced by oomph::RefineableSolidQElement< 2 >::get_jacobian(), oomph::RefineableSolidQElement< 3 >::get_jacobian(), oomph::HopfHandler::get_jacobian(), oomph::PitchForkHandler::get_jacobian(), oomph::FoldHandler::get_jacobian(), oomph::AssemblyHandler::get_jacobian(), oomph::PitchForkHandler::get_residuals(), and oomph::FoldHandler::get_residuals().

virtual void oomph::GeneralisedElement::get_jacobian_and_mass_matrix ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian,
DenseMatrix< double > &  mass_matrix 
) [inline, virtual]

Calculate the residuals and jacobian and elemental "mass" matrix, the matrix that multiplies the time derivative terms.

Definition at line 825 of file elements.h.

References fill_in_contribution_to_jacobian_and_mass_matrix(), oomph::DenseMatrix< T >::initialise(), and oomph::Vector< _Tp >::initialise().

Referenced by oomph::EigenProblemHandler::get_all_vectors_and_matrices(), oomph::HopfHandler::get_jacobian(), and oomph::HopfHandler::get_residuals().

virtual void oomph::GeneralisedElement::get_mass_matrix ( Vector< double > &  residuals,
DenseMatrix< double > &  mass_matrix 
) [inline, virtual]

Calculate the residuals and the elemental "mass" matrix, the matrix that multiplies the time derivative terms in a problem.

Definition at line 812 of file elements.h.

References fill_in_contribution_to_mass_matrix(), oomph::DenseMatrix< T >::initialise(), and oomph::Vector< _Tp >::initialise().

Referenced by oomph::ExplicitTimeStepHandler::get_all_vectors_and_matrices(), and oomph::ExplicitTimeStepHandler::get_jacobian().

virtual void oomph::GeneralisedElement::get_residuals ( Vector< double > &  residuals  )  [inline, virtual]

Calculate the vector of residuals of the equations in the element. By default initialise the vector to zero and then call the fill_in_contribution_to_residuals() function. Note that this entire function can be overloaded if desired.

Reimplemented in oomph::PseudoBucklingRingElement, oomph::ImposeFluxForWomersleyElement< DIM >, oomph::Hijacked< oomph::SpineElement< oomph::FaceGeometry< ELEMENT > > >, and oomph::Hijacked< oomph::SpineElement< oomph::FaceGeometry< oomph::FaceGeometry< ELEMENT > > > >.

Definition at line 790 of file elements.h.

References fill_in_contribution_to_residuals(), and oomph::Vector< _Tp >::initialise().

Referenced by oomph::KirchhoffLoveShellEquations::fill_in_contribution_to_jacobian(), oomph::FSIWallElement::fill_in_contribution_to_jacobian(), oomph::SolidFiniteElement::fill_in_contribution_to_jacobian(), oomph::FiniteElement::fill_in_contribution_to_jacobian(), fill_in_contribution_to_jacobian(), oomph::ElementWithExternalElement::fill_in_contribution_to_jacobian(), oomph::KirchhoffLoveBeamEquations::fill_in_contribution_to_jacobian(), fill_in_jacobian_from_external_by_fd(), oomph::ElementWithExternalElement::fill_in_jacobian_from_external_interaction_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::ElementWithMovingNodes::fill_in_jacobian_from_geometric_data(), fill_in_jacobian_from_internal_by_fd(), oomph::RefineableElement::fill_in_jacobian_from_nodal_by_fd(), oomph::FiniteElement::fill_in_jacobian_from_nodal_by_fd(), oomph::RefineableSolidElement::fill_in_jacobian_from_solid_position_by_fd(), oomph::SolidFiniteElement::fill_in_jacobian_from_solid_position_by_fd(), oomph::RefineableElement::get_dresidual_dnodal_coordinates(), oomph::FiniteElement::get_dresidual_dnodal_coordinates(), oomph::PitchForkHandler::get_residuals(), oomph::FoldHandler::get_residuals(), oomph::ExplicitTimeStepHandler::get_residuals(), and oomph::AssemblyHandler::get_residuals().

void oomph::GeneralisedElement::include_external_data_fd ( const unsigned &  i  )  [inline]

Set the boolean flag to include the external datum in the the finite difference loop when computing the jacobian matrix.

Definition at line 731 of file elements.h.

References Data_fd, Nexternal_data, Ninternal_data, and OOMPH_EXCEPTION_LOCATION.

void oomph::GeneralisedElement::include_internal_data_fd ( const unsigned &  i  )  [inline, protected]

Set the boolean flag to include the internal datum in the finite difference loop when computing the jacobian matrix.

Definition at line 186 of file elements.h.

References Data_fd, Ninternal_data, and OOMPH_EXCEPTION_LOCATION.

bool oomph::GeneralisedElement::internal_data_fd ( const unsigned &  i  )  const [inline, protected]

Return the status of the boolean flag indicating whether the internal data is included in the finite difference loop.

Definition at line 145 of file elements.h.

References Data_fd, Ninternal_data, and OOMPH_EXCEPTION_LOCATION.

Referenced by fill_in_jacobian_from_internal_by_fd().

Data* const& oomph::GeneralisedElement::internal_data_pt ( const unsigned &  i  )  const [inline]

Return a pointer to i-th internal data object (const version).

Definition at line 571 of file elements.h.

References Data_pt, Ninternal_data, and OOMPH_EXCEPTION_LOCATION.

Data*& oomph::GeneralisedElement::internal_data_pt ( const unsigned &  i  )  [inline]

Return a pointer to i-th internal data object.

Definition at line 552 of file elements.h.

References Data_pt, Ninternal_data, and OOMPH_EXCEPTION_LOCATION.

Referenced by add_internal_data(), assign_internal_and_external_local_eqn_numbers(), assign_internal_eqn_numbers(), oomph::FSI_functions::doc_fsi(), fill_in_jacobian_from_internal_by_fd(), oomph::QCrouzeixRaviartElement< DIM >::fix_pressure(), oomph::AxisymmetricQCrouzeixRaviartElement::fix_pressure(), oomph::QPVDElementWithPressure< DIM >::fix_solid_pressure(), oomph::AxisymQPVDElementWithPressure::fix_solid_pressure(), oomph::RefineableQCrouzeixRaviartElement< DIM >::further_build(), oomph::RefineableAxisymmetricQCrouzeixRaviartElement::further_build(), internal_local_eqn(), oomph::AxisymmetricQCrouzeixRaviartElement::p_axi_nst(), oomph::QCrouzeixRaviartElement< DIM >::p_nst(), oomph::RefineableQPVDElementWithPressure< DIM >::rebuild_from_sons(), oomph::RefineableQCrouzeixRaviartElement< DIM >::rebuild_from_sons(), oomph::RefineableAxisymmetricQCrouzeixRaviartElement::rebuild_from_sons(), self_test(), oomph::QPVDElementWithPressure< DIM >::set_solid_p(), oomph::QPVDElementWithPressure< DIM >::solid_p(), and oomph::AxisymQPVDElementWithPressure::solid_p().

int oomph::GeneralisedElement::internal_local_eqn ( const unsigned &  i,
const unsigned &  j 
) [inline, protected]

Return the local equation number corresponding to the j-th value stored at the i-th internal data.

Definition at line 248 of file elements.h.

References Data_local_eqn, internal_data_pt(), Ninternal_data, oomph::Data::nvalue(), and OOMPH_EXCEPTION_LOCATION.

Referenced by oomph::DisplacementControlElement::assign_additional_local_eqn_numbers(), oomph::NavierStokesWomersleyPressureControlElement::fill_in_generic_residual_contribution_pressure_control(), fill_in_jacobian_from_internal_by_fd(), oomph::PseudoBucklingRingElement::geometric_local_eqn(), oomph::ImposeFluxForWomersleyElement< DIM >::get_residuals(), oomph::QCrouzeixRaviartElement< DIM >::p_local_eqn(), oomph::AxisymmetricQCrouzeixRaviartElement::p_local_eqn(), oomph::QPVDElementWithPressure< DIM >::solid_p_local_eqn(), and oomph::AxisymQPVDElementWithPressure::solid_p_local_eqn().

int oomph::GeneralisedElement::local_eqn_number ( const unsigned long &  ieqn_global  )  [inline]

Return the local equation number corresponding to the ieqn_global-th global equation number. Returns minus one (-1) if there is no local degree of freedom corresponding to the chosen global equation number.

Definition at line 660 of file elements.h.

References Eqn_number, and Ndof.

Referenced by oomph::SpectralElement::assign_all_generic_local_eqn_numbers(), oomph::ElementWithMovingNodes::assign_all_generic_local_eqn_numbers(), oomph::ElementWithExternalElement::assign_external_interaction_data_local_eqn_numbers(), oomph::RefineableElement::assign_hanging_local_eqn_numbers(), 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(), oomph::TTaylorHoodElement< DIM >::get_dof_numbers_for_unknowns(), oomph::QTaylorHoodElement< DIM >::get_dof_numbers_for_unknowns(), oomph::QCrouzeixRaviartElement< DIM >::get_dof_numbers_for_unknowns(), and oomph::DisplacementControlElement::get_dof_numbers_for_unknowns().

unsigned oomph::GeneralisedElement::ndof (  )  const [inline]

Return the number of equations/dofs in the element.

Definition at line 763 of file elements.h.

References Ndof.

Referenced by oomph::NavierStokesLSCPreconditioner::assemble_velocity_mass_matrix_diagonal(), oomph::SpectralElement::assign_all_generic_local_eqn_numbers(), oomph::ElementWithMovingNodes::assign_all_generic_local_eqn_numbers(), oomph::ElementWithExternalElement::assign_external_interaction_data_local_eqn_numbers(), oomph::RefineableElement::assign_hanging_local_eqn_numbers(), 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(), oomph::HopfHandler::eqn_number(), oomph::PitchForkHandler::eqn_number(), oomph::FoldHandler::eqn_number(), oomph::KirchhoffLoveShellEquations::fill_in_contribution_to_jacobian(), oomph::SolidFiniteElement::fill_in_contribution_to_jacobian(), oomph::FiniteElement::fill_in_contribution_to_jacobian(), oomph::ElementWithExternalElement::fill_in_contribution_to_jacobian(), oomph::KirchhoffLoveBeamEquations::fill_in_contribution_to_jacobian(), oomph::NavierStokesEquations< DIM >::fill_in_generic_residual_contribution_nst(), fill_in_jacobian_from_external_by_fd(), oomph::ElementWithExternalElement::fill_in_jacobian_from_external_interaction_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::ElementWithMovingNodes::fill_in_jacobian_from_geometric_data(), fill_in_jacobian_from_internal_by_fd(), oomph::RefineableElement::fill_in_jacobian_from_nodal_by_fd(), oomph::FiniteElement::fill_in_jacobian_from_nodal_by_fd(), oomph::RefineableSolidElement::fill_in_jacobian_from_solid_position_by_fd(), oomph::SolidFiniteElement::fill_in_jacobian_from_solid_position_by_fd(), oomph::FoldHandler::FoldHandler(), oomph::EigenProblemHandler::get_all_vectors_and_matrices(), oomph::RefineableNavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::RefineableElement::get_dresidual_dnodal_coordinates(), oomph::NavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::FiniteElement::get_dresidual_dnodal_coordinates(), oomph::DGElement::get_inverse_mass_matrix_times_residuals(), oomph::ImposeFluxForWomersleyElement< DIM >::get_jacobian(), oomph::HopfHandler::get_jacobian(), oomph::PitchForkHandler::get_jacobian(), oomph::FoldHandler::get_jacobian(), oomph::HopfHandler::get_residuals(), oomph::PitchForkHandler::get_residuals(), oomph::FoldHandler::get_residuals(), oomph::NavierStokesEquations< DIM >::get_velocity_mass_matrix_diagonal(), oomph::HopfHandler::HopfHandler(), oomph::HopfHandler::ndof(), oomph::PitchForkHandler::ndof(), oomph::FoldHandler::ndof(), oomph::EigenProblemHandler::ndof(), oomph::ExplicitTimeStepHandler::ndof(), oomph::AssemblyHandler::ndof(), oomph::METIS::partition_mesh(), oomph::PitchForkHandler::PitchForkHandler(), oomph::DGElement::pre_compute_mass_matrix(), oomph::BlockHopfLinearSolver::solve(), oomph::AugmentedBlockFoldLinearSolver::solve(), and oomph::BlockHopfLinearSolver::solve_for_two_rhs().

virtual unsigned oomph::GeneralisedElement::ndof_types (  )  [inline, virtual]

The number of types of degrees of freedom in this element are sub-divided into.

Reimplemented in oomph::QCrouzeixRaviartElement< DIM >, oomph::QTaylorHoodElement< DIM >, oomph::NetFluxControlElement, oomph::TTaylorHoodElement< DIM >, oomph::DisplacementControlElement, oomph::PVDEquationsBase< DIM >, oomph::FSISolidTractionElement< ELEMENT, DIM >, oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >, oomph::FSIHermiteBeamElement, oomph::FSIDiagHermiteShellElement, oomph::ClampedHermiteShellBoundaryConditionElement, oomph::NavierStokesWomersleyPressureControlElement, and oomph::NetFluxControlElementForWomersleyPressureControl.

Definition at line 860 of file elements.h.

References OOMPH_EXCEPTION_LOCATION.

Referenced by oomph::BlockPreconditioner< oomph::CRDoubleMatrix >::set_mesh().

unsigned oomph::GeneralisedElement::nexternal_data (  )  const [inline]

Return the number of external data objects.

Definition at line 760 of file elements.h.

References Nexternal_data.

Referenced by oomph::ElementWithExternalElement::assign_external_interaction_data_local_eqn_numbers(), oomph::DGElement::get_inverse_mass_matrix_times_residuals(), oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::set_boundary_shape_geom_object_pt(), and oomph::FluidInterfaceElement::set_external_pressure_data().

unsigned oomph::GeneralisedElement::ninternal_data (  )  const [inline]

Return the number of internal data objects.

Definition at line 757 of file elements.h.

References Ninternal_data.

Referenced by oomph::ElementWithExternalElement::assign_external_interaction_data_local_eqn_numbers(), oomph::Problem::copy(), oomph::FSI_functions::doc_fsi(), oomph::RefineableElement::identify_field_data_for_interactions(), oomph::FiniteElement::identify_field_data_for_interactions(), oomph::QCrouzeixRaviartElement< DIM >::identify_pressure_data(), and oomph::Problem::set_pinned_values_to_zero().

void oomph::GeneralisedElement::operator= ( const GeneralisedElement  )  [inline]

Broken assignment operator.

Definition at line 536 of file elements.h.

References oomph::BrokenCopy::broken_assign().

virtual void oomph::GeneralisedElement::reset_after_external_fd (  )  [inline, protected, virtual]

Function that is call after the finite differencing of the external data. This may be overloaded to reset any slaved variables that may have changed during the finite differencing.

Reimplemented in oomph::FSIWallElement, and oomph::FluidInterfaceEdgeElement.

Definition at line 462 of file elements.h.

Referenced by fill_in_jacobian_from_external_by_fd().

virtual void oomph::GeneralisedElement::reset_after_internal_fd (  )  [inline, protected, virtual]

Function that is call after the finite differencing of the internal data. This may be overloaded to reset any slaved variables that may have changed during the finite differencing.

Reimplemented in oomph::FSIWallElement.

Definition at line 440 of file elements.h.

Referenced by fill_in_jacobian_from_internal_by_fd().

virtual void oomph::GeneralisedElement::reset_in_external_fd ( const unsigned &  i  )  [inline, protected, virtual]

Function called within the finite difference loop for external data after the values in the i-th external data object are reset. The default behaviour is to call the update function.

Reimplemented in oomph::FSIWallElement, and oomph::FluidInterfaceEdgeElement.

Definition at line 472 of file elements.h.

References update_in_external_fd().

Referenced by fill_in_jacobian_from_external_by_fd().

virtual void oomph::GeneralisedElement::reset_in_internal_fd ( const unsigned &  i  )  [inline, protected, virtual]

Function called within the finite difference loop for internal data after the values in the i-th external data object are reset. The default behaviour is to call the update function.

Reimplemented in oomph::FSIWallElement.

Definition at line 450 of file elements.h.

References update_in_internal_fd().

Referenced by fill_in_jacobian_from_internal_by_fd().

unsigned oomph::GeneralisedElement::self_test (  )  [virtual]

Self-test: Have all internal values been classified as pinned/unpinned? Return 0 if OK.

Self-test: Have all internal values been classified as pinned/unpinned? Return 0 if OK.

Reimplemented in oomph::FiniteElement, oomph::PoissonEquations< DIM >, oomph::AdvectionDiffusionEquations< DIM >, oomph::GeneralisedAdvectionDiffusionEquations< DIM >, oomph::UnsteadyHeatEquations< DIM >, oomph::LinearWaveEquations< DIM >, oomph::WomersleyEquations< DIM >, and oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >.

Definition at line 913 of file elements.cc.

References external_data_pt(), internal_data_pt(), Nexternal_data, Ninternal_data, and oomph::oomph_info.

Referenced by oomph::FiniteElement::self_test().

double oomph::GeneralisedElement::time (  )  const

Return the global time, accessed via the time pointer.

Definition at line 197 of file elements.cc.

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

Referenced by oomph::PVDEquationsBase< DIM >::body_force(), oomph::LinearElasticityEquationsBase< DIM >::body_force(), oomph::AxisymmetricSolidTractionElement< ELEMENT >::fill_in_contribution_to_residuals(), oomph::RefineableAxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::AxisymmetricNavierStokesEquations::fill_in_generic_residual_contribution_axi_nst(), oomph::NavierStokesTractionElement< ELEMENT >::fill_in_generic_residual_contribution_fluid_traction(), oomph::RefineableLinearWaveEquations< DIM >::fill_in_generic_residual_contribution_lin_wave(), oomph::LinearWaveEquations< DIM >::fill_in_generic_residual_contribution_lin_wave(), oomph::LinearWaveFluxElement< ELEMENT >::fill_in_generic_residual_contribution_lin_wave_flux(), oomph::RefineableNavierStokesEquations< DIM >::fill_in_generic_residual_contribution_nst(), oomph::NavierStokesEquations< DIM >::fill_in_generic_residual_contribution_nst(), oomph::UnsteadyHeatEquations< DIM >::fill_in_generic_residual_contribution_ust_heat(), oomph::RefineableUnsteadyHeatEquations< DIM >::fill_in_generic_residual_contribution_ust_heat(), oomph::UnsteadyHeatFluxElement< ELEMENT >::fill_in_generic_residual_contribution_ust_heat_flux(), oomph::NavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::get_wind_adv_diff_react(), and oomph::FiniteElement::integrate_fct().

Time* const& oomph::GeneralisedElement::time_pt (  )  const [inline]

Return the pointer to the global time (const version).

Definition at line 546 of file elements.h.

References Time_pt.

Time* & oomph::GeneralisedElement::time_pt (  )  [inline]

Retun the pointer to the global time.

Definition at line 543 of file elements.h.

References Time_pt.

Referenced by oomph::PVDEquationsBase< DIM >::body_force(), and oomph::LinearElasticityEquationsBase< DIM >::body_force().

virtual void oomph::GeneralisedElement::update_before_external_fd (  )  [inline, protected, virtual]

Function that is called before the finite differencing of any external data. This may be overloaded to update any slaved data before finite differencing takes place.

Definition at line 457 of file elements.h.

Referenced by fill_in_jacobian_from_external_by_fd().

virtual void oomph::GeneralisedElement::update_before_internal_fd (  )  [inline, protected, virtual]

Function that is called before the finite differencing of any internal data. This may be overloaded to update any slaved data before finite differencing takes place.

Definition at line 435 of file elements.h.

Referenced by fill_in_jacobian_from_internal_by_fd().

virtual void oomph::GeneralisedElement::update_in_external_fd ( const unsigned &  i  )  [inline, protected, virtual]

Function called within the finite difference loop for external data after a change in any values in the i-th external data object.

Reimplemented in oomph::FSIWallElement, and oomph::FluidInterfaceEdgeElement.

Definition at line 467 of file elements.h.

Referenced by fill_in_jacobian_from_external_by_fd(), and reset_in_external_fd().

virtual void oomph::GeneralisedElement::update_in_internal_fd ( const unsigned &  i  )  [inline, protected, virtual]

Function called within the finite difference loop for internal data after a change in any values in the i-th internal data object.

Reimplemented in oomph::FSIWallElement.

Definition at line 445 of file elements.h.

Referenced by fill_in_jacobian_from_internal_by_fd(), and reset_in_internal_fd().


Member Data Documentation

std::vector<bool> oomph::GeneralisedElement::Data_fd [private]

Storage for boolean flags of size Ninternal_data + Nexternal_data that correspond to the data used in the element. The flags indicate whether the particular internal or external data should be included in the general finite-difference loop in fill_in_jacobian_from_internal_by_fd() or fill_in_jacobian_from_external_by_fd(). The default is that all data WILL be included in the finite-difference loop, but in many circumstances it is possible to treat certain (external) data analytically. The use of an STL vector is optimal for memory use because the booleans are represented as single-bits.

Definition at line 120 of file elements.h.

Referenced by add_external_data(), add_internal_data(), exclude_external_data_fd(), exclude_internal_data_fd(), external_data_fd(), flush_external_data(), include_external_data_fd(), include_internal_data_fd(), and internal_data_fd().

int** oomph::GeneralisedElement::Data_local_eqn [private]

Pointer to array storage for local equation numbers associated with internal and external data. Again, we save storage by using a single pointer to access this information. The first index of the array is of dimension Nineternal_data + Nexternal_data and the second index varies with the number of values stored at the data object. In the most general case, however, the scheme is as memory efficient as possible.

Definition at line 99 of file elements.h.

Referenced by assign_internal_and_external_local_eqn_numbers(), external_local_eqn(), internal_local_eqn(), and ~GeneralisedElement().

Data** oomph::GeneralisedElement::Data_pt [private]

Storage for pointers to internal and external data. The data is of the same type and so can be addressed by a single pointer. The idea is that the array is of total size Ninternal_data + Nexternal_data. The internal data are stored at the beginning of the array and the external data are stored at the end of the array.

Definition at line 90 of file elements.h.

Referenced by add_external_data(), add_internal_data(), assign_internal_and_external_local_eqn_numbers(), external_data_pt(), flush_external_data(), internal_data_pt(), and ~GeneralisedElement().

double oomph::GeneralisedElement::Default_fd_jacobian_step = 1.0e-8 [static]

Double used for the default finite difference step in elemental jacobian calculations.

Default value used as the increment for finite difference calculations of the jacobian matrices

Definition at line 856 of file elements.h.

Referenced by oomph::DGFaceElement::dnumerical_flux_du(), oomph::PVDEquations< DIM >::fill_in_generic_contribution_to_residuals_pvd(), oomph::RefineablePVDEquations< DIM >::fill_in_generic_contribution_to_residuals_pvd(), oomph::PVDEquationsWithPressure< DIM >::fill_in_generic_residual_contribution_pvd_with_pressure(), oomph::RefineablePVDEquationsWithPressure< DIM >::fill_in_generic_residual_contribution_pvd_with_pressure(), 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::ElementWithMovingNodes::fill_in_jacobian_from_geometric_data(), fill_in_jacobian_from_internal_by_fd(), oomph::RefineableElement::fill_in_jacobian_from_nodal_by_fd(), oomph::FiniteElement::fill_in_jacobian_from_nodal_by_fd(), oomph::RefineableSolidElement::fill_in_jacobian_from_solid_position_by_fd(), oomph::SolidFiniteElement::fill_in_jacobian_from_solid_position_by_fd(), oomph::NavierStokesEquations< DIM >::get_body_force_gradient_nst(), oomph::ElementWithMovingNodes::get_dnodal_coordinates_dgeom_dofs(), oomph::RefineablePoissonEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::RefineableNavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::RefineableElement::get_dresidual_dnodal_coordinates(), oomph::PoissonEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::NavierStokesEquations< DIM >::get_dresidual_dnodal_coordinates(), oomph::FiniteElement::get_dresidual_dnodal_coordinates(), oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::get_reaction_deriv_adv_diff_react(), oomph::NavierStokesEquations< DIM >::get_source_gradient_nst(), oomph::PoissonEquations< DIM >::get_source_gradient_poisson(), and oomph::FiniteElement::locate_zeta().

DenseMatrix< double > oomph::GeneralisedElement::Dummy_matrix [static, protected]

Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case when only the residuals are being assembled.

Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case when only the residuals are being assembled

Definition at line 218 of file elements.h.

Referenced by oomph::NavierStokesEquations< DIM >::fill_in_contribution_to_jacobian(), oomph::GeneralisedAdvectionDiffusionEquations< DIM >::fill_in_contribution_to_jacobian(), oomph::AxisymmetricNavierStokesEquations::fill_in_contribution_to_jacobian(), oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::fill_in_contribution_to_jacobian(), oomph::AdvectionDiffusionEquations< DIM >::fill_in_contribution_to_jacobian(), oomph::NavierStokesWomersleyPressureControlElement::fill_in_contribution_to_residuals(), oomph::NavierStokesImpedanceTractionElement< BULK_NAVIER_STOKES_ELEMENT, WOMERSLEY_ELEMENT, DIM >::fill_in_contribution_to_residuals(), oomph::WomersleyEquations< DIM >::fill_in_contribution_to_residuals(), oomph::UnsteadyHeatFluxElement< ELEMENT >::fill_in_contribution_to_residuals(), oomph::UnsteadyHeatEquations< DIM >::fill_in_contribution_to_residuals(), oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::fill_in_contribution_to_residuals(), oomph::PVDEquationsWithPressure< DIM >::fill_in_contribution_to_residuals(), oomph::PVDEquations< DIM >::fill_in_contribution_to_residuals(), oomph::PoissonFluxElement< ELEMENT >::fill_in_contribution_to_residuals(), oomph::PoissonEquations< DIM >::fill_in_contribution_to_residuals(), oomph::RefineableNavierStokesFluxControlElement< ELEMENT >::fill_in_contribution_to_residuals(), oomph::NavierStokesFluxControlElement< ELEMENT >::fill_in_contribution_to_residuals(), oomph::NavierStokesEquations< DIM >::fill_in_contribution_to_residuals(), oomph::LinearWaveFluxElement< ELEMENT >::fill_in_contribution_to_residuals(), oomph::LinearWaveEquations< DIM >::fill_in_contribution_to_residuals(), oomph::LinearElasticityEquations< DIM >::fill_in_contribution_to_residuals(), oomph::FluidInterfaceElement::fill_in_contribution_to_residuals(), oomph::FluidInterfaceEdgeElement::fill_in_contribution_to_residuals(), oomph::ImposeParallelOutflowElement< ELEMENT >::fill_in_contribution_to_residuals(), oomph::GeneralisedAdvectionDiffusionEquations< DIM >::fill_in_contribution_to_residuals(), oomph::NavierStokesTractionElement< ELEMENT >::fill_in_contribution_to_residuals(), oomph::AxisymmetricPVDEquationsWithPressure::fill_in_contribution_to_residuals(), oomph::AxisymmetricNavierStokesEquations::fill_in_contribution_to_residuals(), oomph::AdvectionDiffusionReactionEquations< NREAGENT, DIM >::fill_in_contribution_to_residuals(), oomph::AdvectionDiffusionFluxElement< ELEMENT >::fill_in_contribution_to_residuals(), oomph::AdvectionDiffusionEquations< DIM >::fill_in_contribution_to_residuals(), and oomph::SolidFiniteElement::fill_in_residuals_for_solid_ic().

unsigned long* oomph::GeneralisedElement::Eqn_number [private]

Storage for the global equation numbers represented by the degrees of freedom in the element.

Definition at line 82 of file elements.h.

Referenced by add_global_eqn_numbers(), assign_local_eqn_numbers(), clear_global_eqn_numbers(), eqn_number(), local_eqn_number(), and ~GeneralisedElement().

unsigned oomph::GeneralisedElement::Ndof [private]

Number of degrees of freedom.

Definition at line 102 of file elements.h.

Referenced by add_global_eqn_numbers(), assign_local_eqn_numbers(), clear_global_eqn_numbers(), eqn_number(), fill_in_contribution_to_jacobian(), fill_in_jacobian_from_external_by_fd(), fill_in_jacobian_from_internal_by_fd(), local_eqn_number(), and ndof().

unsigned oomph::GeneralisedElement::Nexternal_data [private]

Number of external data.

Definition at line 108 of file elements.h.

Referenced by add_external_data(), add_internal_data(), assign_internal_and_external_local_eqn_numbers(), exclude_external_data_fd(), external_data_fd(), external_data_pt(), external_local_eqn(), fill_in_jacobian_from_external_by_fd(), flush_external_data(), include_external_data_fd(), nexternal_data(), and self_test().

unsigned oomph::GeneralisedElement::Ninternal_data [private]

Number of internal data.

Definition at line 105 of file elements.h.

Referenced by add_external_data(), add_internal_data(), assign_internal_and_external_local_eqn_numbers(), assign_internal_eqn_numbers(), exclude_external_data_fd(), exclude_internal_data_fd(), external_data_fd(), external_data_pt(), external_local_eqn(), fill_in_jacobian_from_internal_by_fd(), flush_external_data(), include_external_data_fd(), include_internal_data_fd(), internal_data_fd(), internal_data_pt(), internal_local_eqn(), ninternal_data(), self_test(), and ~GeneralisedElement().

bool oomph::GeneralisedElement::Suppress_warning_about_repeated_external_data = true [static]

Static boolean to suppress warnings about repeated external data. Defaults to true.

Static boolean to suppress warnings about repeated external data. Defaults to true

Definition at line 634 of file elements.h.

Referenced by add_external_data().

bool oomph::GeneralisedElement::Suppress_warning_about_repeated_internal_data = false [static]

Static boolean to suppress warnings about repeated internal data. Defaults to false.

Static boolean to suppress warnings about repeated internal data. Defaults to false

Definition at line 630 of file elements.h.

Referenced by add_internal_data().

Time* oomph::GeneralisedElement::Time_pt [private]

Pointer to global time.

Definition at line 78 of file elements.h.

Referenced by time(), and time_pt().


The documentation for this class was generated from the following files:
Generated on Mon Aug 10 11:24:31 2009 by  doxygen 1.4.7