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::SolidFiniteElement Class Reference

#include <elements.h>

Inheritance diagram for oomph::SolidFiniteElement:

oomph::FiniteElement oomph::GeneralisedElement oomph::GeomObject oomph::AxisymmetricPVDEquations oomph::AxisymmetricPVDEquationsWithPressure oomph::ElasticLineFluidInterfaceEdgeElement< ELEMENT > oomph::ElasticPointFluidInterfaceEdgeElement< ELEMENT > oomph::FSIWallElement oomph::KirchhoffLoveBeamEquations oomph::KirchhoffLoveShellEquations oomph::PVDEquationsBase< DIM > oomph::QSolidElementBase oomph::RefineableSolidElement oomph::SolidFaceElement oomph::SolidPointElement oomph::SolidQHermiteElement< DIM > oomph::SolidQHermiteElement< 1 > oomph::SolidQHermiteElement< 2 > oomph::StorableShapeSolidElementBase oomph::TSolidElementBase List of all members.

Detailed Description

SolidFiniteElement class.

Solid elements are elements whose nodal positions are unknowns in the problem -- their nodes are SolidNodes. In such elements, the nodes not only have a variable (Eulerian) but also a fixed (Lagrangian) position. The positional variables have their own local equation numbering scheme which is set up with.

 assign_solid_local_eqn_numbers () 
The derivatives of the `solid equations' (i.e. the equations that determine the nodal positions) with respect to the nodal positions, required in the Jacobian matrix, are determined by finite differencing.

In the present form, the SolidFiniteElement represents a fully functional base class for `proper' solid mechanics elements, but it can also be combined (via inheritance) with elements that solve additional equations. This is particularly useful in cases where the solid equations are merely used to update the nodal positions in a moving mesh problem, although this can prove costly in practice.

Definition at line 2482 of file elements.h.


Public Types

typedef double(*) MultiplierFctPt (const Vector< double > &xi)
 Pointer to function that computes the "multiplier" for the inertia terms in the consistent determination of the initial conditions for Newmark timestepping.

Public Member Functions

void set_lagrangian_dimension (const unsigned &lagrangian_dimension)
 Set the lagrangian dimension of the element --- the number of lagrangian coordinates stored at the nodes in the element.
virtual bool has_internal_solid_data ()
 Return whether there is internal solid data (e.g. discontinuous solid pressure). At present, this is used to report an error in setting initial conditions for ElasticProblems which can't handle such cases. The default is false.
 SolidFiniteElement ()
 Constructor: Set defaults.
virtual ~SolidFiniteElement ()
 Destructor for SolidFiniteElement:.
 SolidFiniteElement (const SolidFiniteElement &)
 Broken copy constructor.
void operator= (const SolidFiniteElement &)
 Broken assignment operator.
unsigned ngeom_data () const
 The number of geometric data affecting a SolidFiniteElemnet is the same as the number of nodes (one variable position data per node).
Datageom_data_pt (const unsigned &j)
 Return pointer to the j-th Data item that the object's shape depends on. (Redirects to the nodes' positional Data).
double zeta_nodal (const unsigned &n, const unsigned &k, const unsigned &i) const
 In a SolidFiniteElement, the "global" intrinsic coordinate of the element when viewed as part of a compound geometric object (a Mesh) is, by default, the Lagrangian coordinate Note the assumption here is that we are always using isoparameteric elements in which the lagrangian coordinate is interpolated by the same shape functions as the eulerian coordinate.
virtual void get_x_and_xi (const Vector< double > &s, Vector< double > &x_fe, Vector< double > &x, Vector< double > &xi_fe, Vector< double > &xi) const
 Eulerian and Lagrangian coordinates as function of the local coordinates: The Eulerian position is returned in FE-interpolated form (x_fe) and then in the form obtained from the "current" MacroElement representation (if it exists -- if not, x is the same as x_fe). This allows the Domain/MacroElement- based representation to be used to apply displacement boundary conditions exactly. Ditto for the Lagrangian coordinates returned in xi_fe and xi. (Broken virtual -- overload in specific geometric element class if you want to use this functionality.).
virtual void set_macro_elem_pt (MacroElement *macro_elem_pt)
 Set pointer to MacroElement -- overloads generic version and uses the MacroElement also as the default for the "undeformed" configuration. This assignment must be overwritten with set_undeformed_macro_elem_pt(...) if the deformation of the solid body is driven by a deformation of the "current" Domain/MacroElement representation of it's boundary. Can be overloaded in derived classes to perform additional tasks.
virtual void set_macro_elem_pt (MacroElement *macro_elem_pt, MacroElement *undeformed_macro_elem_pt)
 Set pointers to "current" and "undeformed" MacroElements. Can be overloaded in derived classes to perform additional tasks.
void set_undeformed_macro_elem_pt (MacroElement *undeformed_macro_elem_pt)
 Set pointer to "undeformed" macro element. Can be overloaded in derived classes to perform additional tasks.
MacroElementundeformed_macro_elem_pt ()
 Access function to pointer to "undeformed" macro element.
double dshape_lagrangian (const Vector< double > &s, Shape &psi, DShape &dpsidxi) const
 Calculate shape functions and derivatives w.r.t. Lagrangian coordinates at local coordinate s. Returns the Jacobian of the mapping from Lagrangian to local coordinates.
virtual double dshape_lagrangian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidxi) const
 Return the geometric shape functions and also first derivatives w.r.t. Lagrangian coordinates at ipt-th integration point.
double d2shape_lagrangian (const Vector< double > &s, Shape &psi, DShape &dpsidxi, DShape &d2psidxi) const
 Compute the geometric shape functions and also first and second derivatives w.r.t. Lagrangian coordinates at local coordinate s; Returns Jacobian of mapping from Lagrangian to local coordinates.

Numbering:
1D:
d2pidxi(i,0) = $ d^2 \psi_j / d \xi^2 $
2D:
d2psidxi(i,0) = $ \partial^2 \psi_j / \partial \xi_0^2 $
d2psidxi(i,1) = $ \partial^2 \psi_j / \partial \xi_1^2 $
d2psidxi(i,2) = $ \partial^2 \psi_j / \partial \xi_0 \partial \xi_1 $

3D:
d2psidxi(i,0) = $ \partial^2 \psi_j / \partial \xi_0^2 $
d2psidxi(i,1) = $ \partial^2 \psi_j / \partial \xi_1^2 $
d2psidxi(i,2) = $ \partial^2 \psi_j / \partial \xi_2^2 $
d2psidxi(i,3) = $ \partial^2 \psi_j/\partial \xi_0 \partial \xi_1 $
d2psidxi(i,4) = $ \partial^2 \psi_j/\partial \xi_0 \partial \xi_2 $
d2psidxi(i,5) = $ \partial^2 \psi_j/\partial \xi_1 \partial \xi_2 $
.
virtual double d2shape_lagrangian_at_knot (const unsigned &ipt, Shape &psi, DShape &dpsidxi, DShape &d2psidxi) const
 Compute the geometric shape functions and also first and second derivatives w.r.t. Lagrangian coordinates at the ipt-th integration point Returns Jacobian of mapping from Lagrangian to local coordinates.

Numbering:
1D:
d2pidxi(i,0) = $ d^2 \psi_j / d \xi^2 $
2D:
d2psidxi(i,0) = $ \partial^2 \psi_j / \partial \xi_0^2 $
d2psidxi(i,1) = $ \partial^2 \psi_j / \partial \xi_1^2 $
d2psidxi(i,2) = $ \partial^2 \psi_j/\partial \xi_0 \partial \xi_1 $

3D:
d2psidxi(i,0) = $ \partial^2 \psi_j / \partial \xi_0^2 $
d2psidxi(i,1) = $ \partial^2 \psi_j / \partial \xi_1^2 $
d2psidxi(i,2) = $ \partial^2 \psi_j / \partial \xi_2^2 $
d2psidxi(i,3) = $ \partial^2 \psi_j/\partial \xi_0 \partial \xi_1 $
d2psidxi(i,4) = $ \partial^2 \psi_j/\partial \xi_0 \partial \xi_2 $
d2psidxi(i,5) = $ \partial^2 \psi_j/\partial \xi_1 \partial \xi_2 $
.
unsigned lagrangian_dimension () const
 Return the number of Lagrangian coordinates that the element requires at all nodes. This is by default the elemental dimension. If we ever need any other case, it can be implemented.
unsigned nnodal_lagrangian_type () const
 Return the number of types of (generalised) nodal Lagrangian coordinates required to interpolate the Lagrangian coordinates in the element. (E.g. 1 for Lagrange-type elements; 2 for Hermite beam elements; 4 for Hermite shell elements). Default value is 1. Needs to be overloaded for any other element.
Nodeconstruct_node (const unsigned &n)
 Construct the local node n and return a pointer to it.
Nodeconstruct_node (const unsigned &n, TimeStepper *const &time_stepper_pt)
 Construct the local node n and return a pointer to it. Additionally, create storage for `history' values as required by timestepper.
Nodeconstruct_boundary_node (const unsigned &n)
 Construct the local node n and return a pointer to it. in the case when it is a boundary node; that is it MAY be located on a Mesh boundary.
Nodeconstruct_boundary_node (const unsigned &n, TimeStepper *const &time_stepper_pt)
 Construct the local node n and return a pointer to it, in the case when the node MAY be located on a boundary. Additionally, create storage for `history' values as required by timestepper.
virtual void assign_all_generic_local_eqn_numbers ()
 Overload assign_all_generic_local_equation numbers to include the data associated with solid dofs. It remains virtual so that it can be overloaded by RefineableSolidElements.
double raw_lagrangian_position (const unsigned &n, const unsigned &i) const
 Return i-th Lagrangian coordinate at local node n without using the hanging representation.
double raw_lagrangian_position_gen (const unsigned &n, const unsigned &k, const unsigned &i) const
 Return Generalised Lagrangian coordinate at local node n. `Direction' i, `Type' k. Does not use the hanging node representation.
double lagrangian_position (const unsigned &n, const unsigned &i) const
 Return i-th Lagrangian coordinate at local node n.
double lagrangian_position_gen (const unsigned &n, const unsigned &k, const unsigned &i) const
 Return Generalised Lagrangian coordinate at local node n. `Direction' i, `Type' k.
virtual double interpolated_xi (const Vector< double > &s, const unsigned &i) const
 Return i-th FE-interpolated Lagrangian coordinate xi[i] at local coordinate s.
virtual void interpolated_xi (const Vector< double > &s, Vector< double > &xi) const
 Compute FE interpolated Lagrangian coordinate vector xi[] at local coordinate s as Vector.
virtual void J_lagrangian (const Vector< double > &s) const
 Return the Jacobian of mapping from local to Lagrangian coordinates at local position s. NOT YET IMPLEMENTED.
virtual double J_lagrangian_at_knot (const unsigned &ipt) const
 Return the Jacobian of the mapping from local to Lagrangian coordinates at the ipt-th integration point. NOT YET IMPLEMENTED.
SolidInitialCondition *& solid_ic_pt ()
 Pointer to object that describes the initial condition.
bool & solve_for_consistent_newmark_accel_flag ()
 Access function to flag which indicates which system of equations to solve when assigning initial conditions for time-dependent problems. If true, solve for the history value that corresponds to the acceleration in the Newmark scheme by demanding that the PDE is satisifed at the initial time. In this case the Jacobian is replaced by the mass matrix.
MultiplierFctPtmultiplier_fct_pt ()
 Access function: Pointer to multiplicator function for assignment of consistent assignement of initial conditions for Newmark scheme.
MultiplierFctPt multiplier_fct_pt () const
 Access function: Pointer to multiplicator function for assignment of consistent assignement of initial conditions for Newmark scheme (const version).
virtual void get_residuals_for_solid_ic (Vector< double > &residuals)
 Compute the residuals for the setup of an initial condition. The global equations are:

\[ /// 0 = \int \left( \sum_{j=1}^N \sum_{k=1}^K X_{ijk} \psi_{jk}(\xi_n) /// - \frac{\partial^D R^{(IC)}_i(\xi_n)}{\partial t^D} /// \right) \psi_{lm}(\xi_n) \ dv /// \mbox{ \ \ \ \ for \ \ \ $l=1,...,N, \ \ m=1,...,K$} /// \]

where $ N $ is the number of nodes in the mesh and $ K $ the number of generalised nodal coordinates. The initial shape of the solid body, $ {\bf R}^{(IC)},$ and its time-derivatives are specified via the GeomObject that is stored in the SolidFiniteElement::SolidInitialCondition object. The latter also stores the order of the time-derivative $ D $ to be assigned.

void fill_in_residuals_for_solid_ic (Vector< double > &residuals)
 Fill in the residuals for the setup of an initial condition. The global equations are:

\[ /// 0 = \int \left( \sum_{j=1}^N \sum_{k=1}^K X_{ijk} \psi_{jk}(\xi_n) /// - \frac{\partial^D R^{(IC)}_i(\xi_n)}{\partial t^D} /// \right) \psi_{lm}(\xi_n) \ dv /// \mbox{ \ \ \ \ for \ \ \ $l=1,...,N, \ \ m=1,...,K$} /// \]

where $ N $ is the number of nodes in the mesh and $ K $ the number of generalised nodal coordinates. The initial shape of the solid body, $ {\bf R}^{(IC)},$ and its time-derivatives are specified via the GeomObject that is stored in the SolidFiniteElement::SolidInitialCondition object. The latter also stores the order of the time-derivative $ D $ to be assigned.

void fill_in_jacobian_for_solid_ic (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 Fill in the residuals and Jacobian for the setup of an initial condition. The global equations are:

\[ /// 0 = \int \left( \sum_{j=1}^N \sum_{k=1}^K X_{ijk} \psi_{jk}(\xi_n) /// - \frac{\partial^D R^{(IC)}_i(\xi_n)}{\partial t^D} /// \right) \psi_{lm}(\xi_n) \ dv /// \mbox{ \ \ \ \ for \ \ \ $l=1,...,N, \ \ m=1,...,K$} /// \]

where $ N $ is the number of nodes in the mesh and $ K $ the number of generalised nodal coordinates. The initial shape of the solid body, $ {\bf R}^{(IC)},$ and its time-derivatives are specified via the GeomObject that is stored in the SolidFiniteElement::SolidInitialCondition object. The latter also stores the order of the time-derivative $ D $ to be assigned.

void fill_in_jacobian_for_newmark_accel (DenseMatrix< double > &jacobian)
 Fill in the contributions of the Jacobian matrix for the consistent assignment of the initial "accelerations" in Newmark scheme. In this case the Jacobian is the mass matrix.

Protected Member Functions

void fill_in_generic_jacobian_for_solid_ic (Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
 Helper function to fill in the residuals and (if flag==1) the Jacobian for the setup of an initial condition. The global equations are:

\[ /// 0 = \int \left( \sum_{j=1}^N \sum_{k=1}^K X_{ijk} \psi_{jk}(\xi_n) /// - \frac{\partial^D R^{(IC)}_i(\xi_n)}{\partial t^D} /// \right) \psi_{lm}(\xi_n) \ dv /// \mbox{ \ \ \ \ for \ \ \ $l=1,...,N, \ \ m=1,...,K$} /// \]

where $ N $ is the number of nodes in the mesh and $ K $ the number of generalised nodal coordinates. The initial shape of the solid body, $ {\bf R}^{(IC)},$ and its time-derivatives are specified via the GeomObject that is stored in the SolidFiniteElement::SolidInitialCondition object. The latter also stores the order of the time-derivative $ D $ to be assigned.

void set_nnodal_lagrangian_type (const unsigned &nlagrangian_type)
 Set the number of types required to interpolate the Lagrangian coordinates.
virtual double local_to_lagrangian_mapping (const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 Calculate the mapping from local to lagrangian coordinates, given the derivatives of the shape functions w.r.t. local coorindates. Return the determinant of the jacobian, the jacobian and inverse jacobian.
double local_to_lagrangian_mapping (const DShape &dpsids, DenseMatrix< double > &inverse_jacobian) const
 Calculate the mapping from local to lagrangian coordinates, given the derivatives of the shape functions w.r.t. local coordinates, Return only the determinant of the jacobian and the inverse of the mapping (ds/dx).
virtual double local_to_lagrangian_mapping_diagonal (const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
 Calculate the mapping from local to Lagrangian coordinates given the derivatives of the shape functions w.r.t the local coorindates. assuming that the coordinates are aligned in the direction of the local coordinates, i.e. there are no cross terms and the jacobian is diagonal. This function returns the determinant of the jacobian, the jacobian and the inverse jacobian.
virtual void assign_solid_local_eqn_numbers ()
 Assign local equation numbers for the solid equations in the element.
int position_local_eqn (const unsigned &n, const unsigned &k, const unsigned &j)
 Access function that returns the local equation number that corresponds to the j-th coordinate of the k-th position-type at the n-th local node.
void fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 Overload the fill_in_contribution_to_jacobian() function to use finite differences to calculate the solid residuals by default.
virtual void fill_in_jacobian_from_solid_position_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 Use finite differences to calculate the Jacobian entries corresponding to the solid positions. This version assumes that the residuals vector has already been computed.
void fill_in_jacobian_from_solid_position_by_fd (DenseMatrix< double > &jacobian)
 Use finite differences to calculate the Jacobian entries corresponding to the solid positions.
virtual void update_before_solid_position_fd ()
 Function that is called before the finite differencing of any solid position data. This may be overloaded to update any slaved data before finite differencing takes place.
virtual void reset_after_solid_position_fd ()
 Function that is call after the finite differencing of the solid position data. This may be overloaded to reset any slaved variables that may have changed during the finite differencing.
virtual void update_in_solid_position_fd (const unsigned &i)
 Function called within the finite difference loop for the solid position dat after a change in any values in the n-th node.
virtual void reset_in_solid_position_fd (const unsigned &i)
 Function called within the finite difference loop for solid position data after the values in the i-th node are reset. The default behaviour is to call the update function.

Protected Attributes

MacroElementUndeformed_macro_elem_pt
 Pointer to the element's "undeformed" macro element (NULL by default).
SolidInitialConditionSolid_ic_pt
 Pointer to object that specifies the initial condition.
bool Solve_for_consistent_newmark_accel_flag
 Flag to indicate which system of equations to solve when assigning initial conditions for time-dependent problems. If true, solve for the history value that corresponds to the acceleration in the Newmark scheme by demanding that the PDE is satisifed at the initial time. In this case the Jacobian is replaced by the mass matrix.

Private Member Functions

virtual void assemble_local_to_lagrangian_jacobian (const DShape &dpsids, DenseMatrix< double > &jacobian) const
 Assemble the jacobian matrix for the mapping from local to lagrangian coordinates, given the derivatives of the shape function.
virtual void assemble_local_to_lagrangian_jacobian2 (const DShape &d2psids, DenseMatrix< double > &jacobian2) const
 Assemble the the "jacobian" matrix of second derivatives, given the second derivatives of the shape functions w.r.t. local coordinates.
double multiplier (const Vector< double > &xi)
 Access to the "multiplier" for the inertia terms in the consistent determination of the initial conditions for Newmark timestepping.

Private Attributes

MultiplierFctPt Multiplier_fct_pt
 Pointer to function that computes the "multiplier" for the inertia terms in the consistent determination of the initial conditions for Newmark timestepping.
int * Position_local_eqn
 Array to hold the local equation number information for the solid equations, whatever they may be.
unsigned Lagrangian_dimension
 The Lagrangian dimension of the nodes stored in the element,.
unsigned Nnodal_lagrangian_type
 The number of coordinate types requried to intepolate the Lagrangian coordinates in the element. For Lagrange elements it is 1 (the default). It must be over-ridden by using the set_nlagrangian_type() function in the constructors of elements that use generalised coordinate, e.g. for 1D Hermite elements Nnodal_position_types =2.

Member Typedef Documentation

typedef double(*) oomph::SolidFiniteElement::MultiplierFctPt(const Vector< double > &xi)

Pointer to function that computes the "multiplier" for the inertia terms in the consistent determination of the initial conditions for Newmark timestepping.

Definition at line 2501 of file elements.h.


Constructor & Destructor Documentation

oomph::SolidFiniteElement::SolidFiniteElement (  )  [inline]

Constructor: Set defaults.

Definition at line 2504 of file elements.h.

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

Destructor for SolidFiniteElement:.

Definition at line 4254 of file elements.cc.

References Position_local_eqn.

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

Broken copy constructor.

Definition at line 2515 of file elements.h.

References oomph::BrokenCopy::broken_copy().


Member Function Documentation

void oomph::SolidFiniteElement::assemble_local_to_lagrangian_jacobian ( const DShape dpsids,
DenseMatrix< double > &  jacobian 
) const [private, virtual]

Assemble the jacobian matrix for the mapping from local to lagrangian coordinates, given the derivatives of the shape function.

Internal function that is used to assemble the jacobian of the mapping from local coordinates (s) to the lagrangian coordinates (xi), given the derivatives of the shape functions.

Reimplemented in oomph::RefineableSolidElement.

Definition at line 4153 of file elements.cc.

References oomph::FiniteElement::dim(), Lagrangian_dimension, nnodal_lagrangian_type(), oomph::FiniteElement::nnode(), OOMPH_EXCEPTION_LOCATION, and raw_lagrangian_position_gen().

Referenced by local_to_lagrangian_mapping().

void oomph::SolidFiniteElement::assemble_local_to_lagrangian_jacobian2 ( const DShape d2psids,
DenseMatrix< double > &  jacobian2 
) const [private, virtual]

Assemble the the "jacobian" matrix of second derivatives, given the second derivatives of the shape functions w.r.t. local coordinates.

Internal function that is used to assemble the jacobian of second derivatives of the the mapping from local coordinates (s) to the lagrangian coordinates (xi), given the second derivatives of the shape functions.

Reimplemented in oomph::RefineableSolidElement.

Definition at line 4212 of file elements.cc.

References oomph::FiniteElement::dim(), oomph::FiniteElement::N2deriv, nnodal_lagrangian_type(), oomph::FiniteElement::nnode(), and raw_lagrangian_position_gen().

Referenced by d2shape_lagrangian(), and d2shape_lagrangian_at_knot().

virtual void oomph::SolidFiniteElement::assign_all_generic_local_eqn_numbers (  )  [inline, virtual]

Overload assign_all_generic_local_equation numbers to include the data associated with solid dofs. It remains virtual so that it can be overloaded by RefineableSolidElements.

Reimplemented from oomph::FiniteElement.

Definition at line 2753 of file elements.h.

References oomph::FiniteElement::assign_all_generic_local_eqn_numbers(), and assign_solid_local_eqn_numbers().

void oomph::SolidFiniteElement::assign_solid_local_eqn_numbers (  )  [protected, virtual]

Assign local equation numbers for the solid equations in the element.

Reimplemented in oomph::RefineableSolidElement.

Definition at line 4484 of file elements.cc.

References oomph::GeneralisedElement::add_global_eqn_numbers(), oomph::GeneralisedElement::eqn_number(), oomph::Data::Is_pinned, oomph::GeneralisedElement::local_eqn_number(), oomph::GeneralisedElement::ndof(), oomph::FiniteElement::nnodal_position_type(), oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_dimension(), oomph::FiniteElement::node_pt(), and Position_local_eqn.

Referenced by assign_all_generic_local_eqn_numbers(), and oomph::RefineableSolidElement::assign_solid_local_eqn_numbers().

Node* oomph::SolidFiniteElement::construct_boundary_node ( const unsigned &  n,
TimeStepper *const &  time_stepper_pt 
) [inline, virtual]

Construct the local node n and return a pointer to it, in the case when the node MAY be located on a boundary. Additionally, create storage for `history' values as required by timestepper.

Reimplemented from oomph::FiniteElement.

Definition at line 2732 of file elements.h.

References lagrangian_dimension(), nnodal_lagrangian_type(), oomph::FiniteElement::nnodal_position_type(), oomph::FiniteElement::nodal_dimension(), oomph::FiniteElement::node_pt(), oomph::FiniteElement::required_nvalue(), and oomph::GeomObject::time_stepper_pt().

Node* oomph::SolidFiniteElement::construct_boundary_node ( const unsigned &  n  )  [inline, virtual]

Construct the local node n and return a pointer to it. in the case when it is a boundary node; that is it MAY be located on a Mesh boundary.

Reimplemented from oomph::FiniteElement.

Definition at line 2712 of file elements.h.

References lagrangian_dimension(), nnodal_lagrangian_type(), oomph::FiniteElement::nnodal_position_type(), oomph::FiniteElement::nodal_dimension(), oomph::FiniteElement::node_pt(), and oomph::FiniteElement::required_nvalue().

Node* oomph::SolidFiniteElement::construct_node ( const unsigned &  n,
TimeStepper *const &  time_stepper_pt 
) [inline, virtual]

Construct the local node n and return a pointer to it. Additionally, create storage for `history' values as required by timestepper.

Reimplemented from oomph::FiniteElement.

Definition at line 2693 of file elements.h.

References lagrangian_dimension(), nnodal_lagrangian_type(), oomph::FiniteElement::nnodal_position_type(), oomph::FiniteElement::nodal_dimension(), oomph::FiniteElement::node_pt(), oomph::FiniteElement::required_nvalue(), and oomph::GeomObject::time_stepper_pt().

Node* oomph::SolidFiniteElement::construct_node ( const unsigned &  n  )  [inline, virtual]

Construct the local node n and return a pointer to it.

Reimplemented from oomph::FiniteElement.

Definition at line 2675 of file elements.h.

References lagrangian_dimension(), nnodal_lagrangian_type(), oomph::FiniteElement::nnodal_position_type(), oomph::FiniteElement::nodal_dimension(), oomph::FiniteElement::node_pt(), and oomph::FiniteElement::required_nvalue().

double oomph::SolidFiniteElement::d2shape_lagrangian ( const Vector< double > &  s,
Shape psi,
DShape dpsidxi,
DShape d2psidxi 
) const

Compute the geometric shape functions and also first and second derivatives w.r.t. Lagrangian coordinates at local coordinate s; Returns Jacobian of mapping from Lagrangian to local coordinates.

Numbering:
1D:
d2pidxi(i,0) = $ d^2 \psi_j / d \xi^2 $
2D:
d2psidxi(i,0) = $ \partial^2 \psi_j / \partial \xi_0^2 $
d2psidxi(i,1) = $ \partial^2 \psi_j / \partial \xi_1^2 $
d2psidxi(i,2) = $ \partial^2 \psi_j / \partial \xi_0 \partial \xi_1 $

3D:
d2psidxi(i,0) = $ \partial^2 \psi_j / \partial \xi_0^2 $
d2psidxi(i,1) = $ \partial^2 \psi_j / \partial \xi_1^2 $
d2psidxi(i,2) = $ \partial^2 \psi_j / \partial \xi_2^2 $
d2psidxi(i,3) = $ \partial^2 \psi_j/\partial \xi_0 \partial \xi_1 $
d2psidxi(i,4) = $ \partial^2 \psi_j/\partial \xi_0 \partial \xi_2 $
d2psidxi(i,5) = $ \partial^2 \psi_j/\partial \xi_1 \partial \xi_2 $
.

Definition at line 4398 of file elements.cc.

References assemble_local_to_lagrangian_jacobian2(), oomph::FiniteElement::d2shape_local(), oomph::FiniteElement::dim(), local_to_lagrangian_mapping(), oomph::FiniteElement::N2deriv, and oomph::FiniteElement::transform_second_derivatives().

double oomph::SolidFiniteElement::d2shape_lagrangian_at_knot ( const unsigned &  ipt,
Shape psi,
DShape dpsidxi,
DShape d2psidxi 
) const [virtual]

Compute the geometric shape functions and also first and second derivatives w.r.t. Lagrangian coordinates at the ipt-th integration point Returns Jacobian of mapping from Lagrangian to local coordinates.

Numbering:
1D:
d2pidxi(i,0) = $ d^2 \psi_j / d \xi^2 $
2D:
d2psidxi(i,0) = $ \partial^2 \psi_j / \partial \xi_0^2 $
d2psidxi(i,1) = $ \partial^2 \psi_j / \partial \xi_1^2 $
d2psidxi(i,2) = $ \partial^2 \psi_j/\partial \xi_0 \partial \xi_1 $

3D:
d2psidxi(i,0) = $ \partial^2 \psi_j / \partial \xi_0^2 $
d2psidxi(i,1) = $ \partial^2 \psi_j / \partial \xi_1^2 $
d2psidxi(i,2) = $ \partial^2 \psi_j / \partial \xi_2^2 $
d2psidxi(i,3) = $ \partial^2 \psi_j/\partial \xi_0 \partial \xi_1 $
d2psidxi(i,4) = $ \partial^2 \psi_j/\partial \xi_0 \partial \xi_2 $
d2psidxi(i,5) = $ \partial^2 \psi_j/\partial \xi_1 \partial \xi_2 $
.

Reimplemented in oomph::StorableShapeSolidElementBase.

Definition at line 4447 of file elements.cc.

References assemble_local_to_lagrangian_jacobian2(), oomph::FiniteElement::d2shape_local_at_knot(), oomph::FiniteElement::dim(), local_to_lagrangian_mapping(), oomph::FiniteElement::N2deriv, and oomph::FiniteElement::transform_second_derivatives().

Referenced by oomph::StorableShapeSolidElementBase::d2shape_lagrangian_at_knot(), oomph::KirchhoffLoveBeamEquations::fill_in_contribution_to_residuals_beam(), oomph::KirchhoffLoveShellEquations::fill_in_contribution_to_residuals_shell(), oomph::KirchhoffLoveShellEquations::get_energy(), oomph::KirchhoffLoveBeamEquations::get_energy(), and oomph::StorableShapeSolidElementBase::pre_compute_d2shape_lagrangian_at_knots().

double oomph::SolidFiniteElement::dshape_lagrangian ( const Vector< double > &  s,
Shape psi,
DShape dpsi 
) const

Calculate shape functions and derivatives w.r.t. Lagrangian coordinates at local coordinate s. Returns the Jacobian of the mapping from Lagrangian to local coordinates.

Calculate shape functions and derivatives w.r.t. Lagrangian coordinates at local coordinate s. Returns the Jacobian of the mapping from Lagrangian to local coordinates. General case, may be overloaded

Definition at line 4328 of file elements.cc.

References oomph::FiniteElement::dim(), oomph::FiniteElement::dshape_local(), local_to_lagrangian_mapping(), and oomph::FiniteElement::transform_derivatives().

Referenced by oomph::FSIDiagHermiteShellElement::dposition_dlagrangian_at_local_coordinate(), oomph::FSIHermiteBeamElement::dposition_dlagrangian_at_local_coordinate(), fill_in_jacobian_for_newmark_accel(), oomph::PVDEquationsBase< DIM >::get_deformed_covariant_basis_vectors(), oomph::KirchhoffLoveShellEquations::get_normal(), oomph::KirchhoffLoveBeamEquations::get_normal(), oomph::PVDEquationsBase< DIM >::get_strain(), and oomph::KirchhoffLoveShellEquations::load_rate_of_work().

double oomph::SolidFiniteElement::dshape_lagrangian_at_knot ( const unsigned &  ipt,
Shape psi,
DShape dpsi 
) const [virtual]

Return the geometric shape functions and also first derivatives w.r.t. Lagrangian coordinates at ipt-th integration point.

Compute the geometric shape functions and also first derivatives w.r.t. Lagrangian coordinates at integration point ipt. Most general form of function, but may be over-loaded if desired

Reimplemented in oomph::StorableShapeSolidElementBase.

Definition at line 4353 of file elements.cc.

References oomph::FiniteElement::dim(), oomph::FiniteElement::dshape_local_at_knot(), local_to_lagrangian_mapping(), and oomph::FiniteElement::transform_derivatives().

Referenced by oomph::StorableShapeSolidElementBase::dshape_lagrangian_at_knot(), oomph::AxisymmetricPVDEquations::fill_in_contribution_to_residuals_axisym_pvd(), oomph::PVDEquations< DIM >::fill_in_generic_contribution_to_residuals_pvd(), oomph::RefineablePVDEquations< DIM >::fill_in_generic_contribution_to_residuals_pvd(), oomph::AxisymmetricPVDEquationsWithPressure::fill_in_generic_residual_contribution_axisym_pvd_with_pressure(), oomph::PVDEquationsWithPressure< DIM >::fill_in_generic_residual_contribution_pvd_with_pressure(), oomph::RefineablePVDEquationsWithPressure< DIM >::fill_in_generic_residual_contribution_pvd_with_pressure(), oomph::StorableShapeSolidElementBase::pre_compute_dshape_lagrangian_at_knots(), oomph::AxisymmetricPVDEquationsWithPressure::size(), and oomph::AxisymmetricPVDEquations::size().

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

Overload the fill_in_contribution_to_jacobian() function to use finite differences to calculate the solid residuals by default.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::FSIWallElement, oomph::ElasticPointFluidInterfaceEdgeElement< ELEMENT >, oomph::ElasticLineFluidInterfaceEdgeElement< ELEMENT >, oomph::PVDEquations< DIM >, oomph::PVDEquationsWithPressure< DIM >, oomph::SolidTractionElement< ELEMENT >, oomph::FSISolidTractionElement< ELEMENT, DIM >, oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >, oomph::AxisymmetricPVDEquations, oomph::AxisymmetricPVDEquationsWithPressure, oomph::KirchhoffLoveBeamEquations, oomph::FSIHermiteBeamElement, oomph::KirchhoffLoveShellEquations, and oomph::FSIDiagHermiteShellElement.

Definition at line 3046 of file elements.h.

References oomph::GeneralisedElement::fill_in_contribution_to_residuals(), fill_in_jacobian_for_newmark_accel(), 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(), fill_in_jacobian_from_solid_position_by_fd(), oomph::GeneralisedElement::get_residuals(), oomph::GeneralisedElement::ndof(), and Solve_for_consistent_newmark_accel_flag.

void oomph::SolidFiniteElement::fill_in_generic_jacobian_for_solid_ic ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian,
const unsigned &  flag 
) [protected]

Helper function to fill in the residuals and (if flag==1) the Jacobian for the setup of an initial condition. The global equations are:

\[ /// 0 = \int \left( \sum_{j=1}^N \sum_{k=1}^K X_{ijk} \psi_{jk}(\xi_n) /// - \frac{\partial^D R^{(IC)}_i(\xi_n)}{\partial t^D} /// \right) \psi_{lm}(\xi_n) \ dv /// \mbox{ \ \ \ \ for \ \ \ $l=1,...,N, \ \ m=1,...,K$} /// \]

where $ N $ is the number of nodes in the mesh and $ K $ the number of generalised nodal coordinates. The initial shape of the solid body, $ {\bf R}^{(IC)},$ and its time-derivatives are specified via the GeomObject that is stored in the SolidFiniteElement::SolidInitialCondition object. The latter also stores the order of the time-derivative $ D $ to be assigned.

Helper function to fill in the residuals and (if flag==1) the Jacobian for the setup of an initial condition. The global equations are:

\[ /// 0 = \int \left( \sum_{j=1}^N \sum_{k=1}^K X_{ijk} \psi_{jk}(\xi_n) /// - \frac{\partial^D R^{(IC)}_i(\xi_n)}{\partial t^D} /// \right) \psi_{lm}(\xi_n) \ dv /// \mbox{ \ \ \ \ for \ \ \ $l=1,...,N, \ \ m=1,...,K$} /// \]

where $ N $ is the number of nodes in the mesh and $ K $ the number of generalised nodal coordinates. The initial shape of the solid body, $ {\bf R}^{(IC)},$ and its time-derivatives are specified via the GeomObject that is stored in the SolidFiniteElement::SolidInitialCondition object. The latter also stores the order of the time-derivative $ D $ to be assigned.

Definition at line 4908 of file elements.cc.

References oomph::FiniteElement::dim(), oomph::FiniteElement::dposition_dt(), oomph::SolidInitialCondition::geom_object_pt(), oomph::SolidInitialCondition::ic_time_deriv(), oomph::FiniteElement::integral_pt(), oomph::FiniteElement::interpolated_x(), lagrangian_position_gen(), nnodal_lagrangian_type(), oomph::FiniteElement::nnodal_position_type(), oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_dimension(), oomph::Integral::nweight(), oomph::oomph_info, position_local_eqn(), oomph::FiniteElement::shape(), Solid_ic_pt, and oomph::Integral::weight().

Referenced by fill_in_jacobian_for_solid_ic(), and fill_in_residuals_for_solid_ic().

void oomph::SolidFiniteElement::fill_in_jacobian_for_newmark_accel ( DenseMatrix< double > &  jacobian  ) 

Fill in the contributions of the Jacobian matrix for the consistent assignment of the initial "accelerations" in Newmark scheme. In this case the Jacobian is the mass matrix.

Add jacobian and residuals for consistent assignment of initial "accelerations" in Newmark scheme. Jacobian is the mass matrix.

Definition at line 4750 of file elements.cc.

References oomph::FiniteElement::dim(), dshape_lagrangian(), oomph::FiniteElement::integral_pt(), interpolated_xi(), multiplier(), oomph::FiniteElement::nnodal_position_type(), oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_dimension(), oomph::FiniteElement::node_pt(), oomph::TimeStepper::ntstorage(), oomph::Integral::nweight(), OOMPH_EXCEPTION_LOCATION, position_local_eqn(), oomph::Node::position_time_stepper_pt(), Solid_ic_pt, Solve_for_consistent_newmark_accel_flag, oomph::TimeStepper::type(), oomph::QuadTreeNames::W, oomph::Integral::weight(), and oomph::TimeStepper::weight().

Referenced by oomph::PVDEquations< DIM >::fill_in_contribution_to_jacobian(), oomph::KirchhoffLoveShellEquations::fill_in_contribution_to_jacobian(), oomph::FSIWallElement::fill_in_contribution_to_jacobian(), fill_in_contribution_to_jacobian(), oomph::KirchhoffLoveBeamEquations::fill_in_contribution_to_jacobian(), and oomph::AxisymmetricPVDEquations::fill_in_contribution_to_jacobian().

void oomph::SolidFiniteElement::fill_in_jacobian_for_solid_ic ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian 
) [inline]

Fill in the residuals and Jacobian for the setup of an initial condition. The global equations are:

\[ /// 0 = \int \left( \sum_{j=1}^N \sum_{k=1}^K X_{ijk} \psi_{jk}(\xi_n) /// - \frac{\partial^D R^{(IC)}_i(\xi_n)}{\partial t^D} /// \right) \psi_{lm}(\xi_n) \ dv /// \mbox{ \ \ \ \ for \ \ \ $l=1,...,N, \ \ m=1,...,K$} /// \]

where $ N $ is the number of nodes in the mesh and $ K $ the number of generalised nodal coordinates. The initial shape of the solid body, $ {\bf R}^{(IC)},$ and its time-derivatives are specified via the GeomObject that is stored in the SolidFiniteElement::SolidInitialCondition object. The latter also stores the order of the time-derivative $ D $ to be assigned.

Definition at line 2905 of file elements.h.

References fill_in_generic_jacobian_for_solid_ic().

Referenced by oomph::PVDEquations< DIM >::fill_in_contribution_to_jacobian().

void oomph::SolidFiniteElement::fill_in_jacobian_from_solid_position_by_fd ( DenseMatrix< double > &  jacobian  )  [inline, protected]

Use finite differences to calculate the Jacobian entries corresponding to the solid positions.

Definition at line 3088 of file elements.h.

References fill_in_jacobian_from_solid_position_by_fd(), oomph::GeneralisedElement::get_residuals(), and oomph::GeneralisedElement::ndof().

void oomph::SolidFiniteElement::fill_in_jacobian_from_solid_position_by_fd ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian 
) [protected, virtual]

Use finite differences to calculate the Jacobian entries corresponding to the solid positions. This version assumes that the residuals vector has already been computed.

This function calculates the entries of Jacobian matrix, used in the Newton method, associated with the elastic problem in which the nodal position is a variable. It does this using finite differences, rather than an analytical formulation, so can be done in total generality.

Reimplemented in oomph::RefineableSolidElement.

Definition at line 4558 of file elements.cc.

References oomph::GeneralisedElement::Default_fd_jacobian_step, oomph::GeneralisedElement::get_residuals(), oomph::GeneralisedElement::ndof(), oomph::FiniteElement::nnodal_position_type(), oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_dimension(), oomph::FiniteElement::node_pt(), oomph::Node::perform_auxiliary_node_update_fct(), position_local_eqn(), reset_after_solid_position_fd(), reset_in_solid_position_fd(), update_before_solid_position_fd(), and update_in_solid_position_fd().

Referenced by oomph::PVDEquationsWithPressure< DIM >::fill_in_contribution_to_jacobian(), oomph::PVDEquations< DIM >::fill_in_contribution_to_jacobian(), oomph::KirchhoffLoveShellEquations::fill_in_contribution_to_jacobian(), oomph::ElasticLineFluidInterfaceEdgeElement< ELEMENT >::fill_in_contribution_to_jacobian(), oomph::ElasticPointFluidInterfaceEdgeElement< ELEMENT >::fill_in_contribution_to_jacobian(), oomph::FSIWallElement::fill_in_contribution_to_jacobian(), fill_in_contribution_to_jacobian(), oomph::KirchhoffLoveBeamEquations::fill_in_contribution_to_jacobian(), oomph::AxisymmetricPVDEquationsWithPressure::fill_in_contribution_to_jacobian(), oomph::AxisymmetricPVDEquations::fill_in_contribution_to_jacobian(), and fill_in_jacobian_from_solid_position_by_fd().

void oomph::SolidFiniteElement::fill_in_residuals_for_solid_ic ( Vector< double > &  residuals  )  [inline]

Fill in the residuals for the setup of an initial condition. The global equations are:

\[ /// 0 = \int \left( \sum_{j=1}^N \sum_{k=1}^K X_{ijk} \psi_{jk}(\xi_n) /// - \frac{\partial^D R^{(IC)}_i(\xi_n)}{\partial t^D} /// \right) \psi_{lm}(\xi_n) \ dv /// \mbox{ \ \ \ \ for \ \ \ $l=1,...,N, \ \ m=1,...,K$} /// \]

where $ N $ is the number of nodes in the mesh and $ K $ the number of generalised nodal coordinates. The initial shape of the solid body, $ {\bf R}^{(IC)},$ and its time-derivatives are specified via the GeomObject that is stored in the SolidFiniteElement::SolidInitialCondition object. The latter also stores the order of the time-derivative $ D $ to be assigned.

Definition at line 2883 of file elements.h.

References oomph::GeneralisedElement::Dummy_matrix, and fill_in_generic_jacobian_for_solid_ic().

Referenced by oomph::KirchhoffLoveBeamEquations::fill_in_contribution_to_residuals_beam(), oomph::PVDEquations< DIM >::fill_in_generic_contribution_to_residuals_pvd(), and get_residuals_for_solid_ic().

Data* oomph::SolidFiniteElement::geom_data_pt ( const unsigned &  j  )  [inline, virtual]

Return pointer to the j-th Data item that the object's shape depends on. (Redirects to the nodes' positional Data).

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::RefineableSolidElement.

Definition at line 2533 of file elements.h.

References oomph::FiniteElement::node_pt().

Referenced by oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::set_boundary_shape_geom_object_pt().

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

Compute the residuals for the setup of an initial condition. The global equations are:

\[ /// 0 = \int \left( \sum_{j=1}^N \sum_{k=1}^K X_{ijk} \psi_{jk}(\xi_n) /// - \frac{\partial^D R^{(IC)}_i(\xi_n)}{\partial t^D} /// \right) \psi_{lm}(\xi_n) \ dv /// \mbox{ \ \ \ \ for \ \ \ $l=1,...,N, \ \ m=1,...,K$} /// \]

where $ N $ is the number of nodes in the mesh and $ K $ the number of generalised nodal coordinates. The initial shape of the solid body, $ {\bf R}^{(IC)},$ and its time-derivatives are specified via the GeomObject that is stored in the SolidFiniteElement::SolidInitialCondition object. The latter also stores the order of the time-derivative $ D $ to be assigned.

Definition at line 2863 of file elements.h.

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

Referenced by oomph::RefineablePVDEquations< DIM >::fill_in_generic_contribution_to_residuals_pvd(), oomph::PVDEquationsWithPressure< DIM >::fill_in_generic_residual_contribution_pvd_with_pressure(), and oomph::RefineablePVDEquationsWithPressure< DIM >::fill_in_generic_residual_contribution_pvd_with_pressure().

virtual void oomph::SolidFiniteElement::get_x_and_xi ( const Vector< double > &  s,
Vector< double > &  x_fe,
Vector< double > &  x,
Vector< double > &  xi_fe,
Vector< double > &  xi 
) const [inline, virtual]

Eulerian and Lagrangian coordinates as function of the local coordinates: The Eulerian position is returned in FE-interpolated form (x_fe) and then in the form obtained from the "current" MacroElement representation (if it exists -- if not, x is the same as x_fe). This allows the Domain/MacroElement- based representation to be used to apply displacement boundary conditions exactly. Ditto for the Lagrangian coordinates returned in xi_fe and xi. (Broken virtual -- overload in specific geometric element class if you want to use this functionality.).

Reimplemented in oomph::QSolidElementBase.

Definition at line 2556 of file elements.h.

References OOMPH_EXCEPTION_LOCATION.

virtual bool oomph::SolidFiniteElement::has_internal_solid_data (  )  [inline, virtual]

Return whether there is internal solid data (e.g. discontinuous solid pressure). At present, this is used to report an error in setting initial conditions for ElasticProblems which can't handle such cases. The default is false.

Reimplemented in oomph::QPVDElementWithPressure< DIM >, and oomph::AxisymQPVDElementWithPressure.

Definition at line 2496 of file elements.h.

void oomph::SolidFiniteElement::interpolated_xi ( const Vector< double > &  s,
Vector< double > &  xi 
) const [virtual]

Compute FE interpolated Lagrangian coordinate vector xi[] at local coordinate s as Vector.

Compute FE-interpolated Lagrangian coordinate vector xi[] at local coordinate s.

Reimplemented in oomph::SolidFaceElement.

Definition at line 4708 of file elements.cc.

References lagrangian_dimension(), lagrangian_position_gen(), nnodal_lagrangian_type(), oomph::FiniteElement::nnode(), and oomph::FiniteElement::shape().

double oomph::SolidFiniteElement::interpolated_xi ( const Vector< double > &  s,
const unsigned &  i 
) const [virtual]

Return i-th FE-interpolated Lagrangian coordinate xi[i] at local coordinate s.

Return i-th FE-interpolated Lagrangian coordinate at local coordinate s.

Reimplemented in oomph::SolidFaceElement.

Definition at line 4673 of file elements.cc.

References lagrangian_position_gen(), nnodal_lagrangian_type(), oomph::FiniteElement::nnode(), and oomph::FiniteElement::shape().

Referenced by oomph::AxisymmetricPVDEquations::fill_in_contribution_to_residuals_axisym_pvd(), oomph::KirchhoffLoveBeamEquations::fill_in_contribution_to_residuals_beam(), oomph::KirchhoffLoveShellEquations::fill_in_contribution_to_residuals_shell(), oomph::PVDEquations< DIM >::fill_in_generic_contribution_to_residuals_pvd(), oomph::RefineablePVDEquations< DIM >::fill_in_generic_contribution_to_residuals_pvd(), oomph::AxisymmetricPVDEquationsWithPressure::fill_in_generic_residual_contribution_axisym_pvd_with_pressure(), 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_for_newmark_accel(), oomph::KirchhoffLoveShellEquations::get_energy(), oomph::KirchhoffLoveBeamEquations::get_energy(), oomph::PVDEquationsBase< DIM >::get_strain(), oomph::PVDEquationsWithPressure< DIM >::get_stress(), oomph::PVDEquations< DIM >::get_stress(), oomph::QSolidElementBase::get_x_and_xi(), oomph::PVDEquationsWithPressure< DIM >::output(), oomph::PVDEquations< DIM >::output(), oomph::ClampedHermiteShellBoundaryConditionElement::output(), oomph::HermiteShellElement::output(), oomph::SolidQHermiteElement< 2 >::output(), oomph::SolidQHermiteElement< 1 >::output(), oomph::HermiteBeamElement::output(), oomph::AxisymmetricPVDEquationsWithPressure::output(), oomph::AxisymDiagHermitePVDElement::output(), oomph::AxisymmetricPVDEquations::output(), oomph::AxisymmetricPVDEquationsWithPressure::size(), and oomph::AxisymmetricPVDEquations::size().

virtual void oomph::SolidFiniteElement::J_lagrangian ( const Vector< double > &  s  )  const [inline, virtual]

Return the Jacobian of mapping from local to Lagrangian coordinates at local position s. NOT YET IMPLEMENTED.

Definition at line 2795 of file elements.h.

References OOMPH_EXCEPTION_LOCATION.

virtual double oomph::SolidFiniteElement::J_lagrangian_at_knot ( const unsigned &  ipt  )  const [inline, virtual]

Return the Jacobian of the mapping from local to Lagrangian coordinates at the ipt-th integration point. NOT YET IMPLEMENTED.

Definition at line 2805 of file elements.h.

References OOMPH_EXCEPTION_LOCATION.

unsigned oomph::SolidFiniteElement::lagrangian_dimension (  )  const [inline]

Return the number of Lagrangian coordinates that the element requires at all nodes. This is by default the elemental dimension. If we ever need any other case, it can be implemented.

Definition at line 2664 of file elements.h.

References Lagrangian_dimension.

Referenced by construct_boundary_node(), construct_node(), and interpolated_xi().

double oomph::SolidFiniteElement::lagrangian_position ( const unsigned &  n,
const unsigned &  i 
) const [inline]

Return i-th Lagrangian coordinate at local node n.

Definition at line 2774 of file elements.h.

References oomph::FiniteElement::node_pt().

Referenced by oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::fill_in_generic_contribution_to_residuals_displ_lagr_multiplier().

double oomph::SolidFiniteElement::lagrangian_position_gen ( const unsigned &  n,
const unsigned &  k,
const unsigned &  i 
) const [inline]

Return Generalised Lagrangian coordinate at local node n. `Direction' i, `Type' k.

Definition at line 2779 of file elements.h.

References oomph::FiniteElement::node_pt().

Referenced by oomph::RefineableSolidElement::assemble_local_to_lagrangian_jacobian(), oomph::RefineableSolidElement::assemble_local_to_lagrangian_jacobian2(), oomph::AxisymmetricPVDEquations::fill_in_contribution_to_residuals_axisym_pvd(), oomph::SolidTractionElement< ELEMENT >::fill_in_contribution_to_residuals_solid_traction(), oomph::RefineablePVDEquations< DIM >::fill_in_generic_contribution_to_residuals_pvd(), fill_in_generic_jacobian_for_solid_ic(), oomph::AxisymmetricPVDEquationsWithPressure::fill_in_generic_residual_contribution_axisym_pvd_with_pressure(), oomph::RefineablePVDEquationsWithPressure< DIM >::fill_in_generic_residual_contribution_pvd_with_pressure(), oomph::KirchhoffLoveShellEquations::get_energy(), interpolated_xi(), oomph::KirchhoffLoveShellEquations::load_rate_of_work(), oomph::RefineableSolidElement::local_to_lagrangian_mapping_diagonal(), oomph::AxisymmetricPVDEquationsWithPressure::size(), oomph::AxisymmetricPVDEquations::size(), and zeta_nodal().

double oomph::SolidFiniteElement::local_to_lagrangian_mapping ( const DShape dpsids,
DenseMatrix< double > &  inverse_jacobian 
) const [inline, protected]

Calculate the mapping from local to lagrangian coordinates, given the derivatives of the shape functions w.r.t. local coordinates, Return only the determinant of the jacobian and the inverse of the mapping (ds/dx).

Definition at line 2964 of file elements.h.

References oomph::FiniteElement::dim(), and local_to_lagrangian_mapping().

virtual double oomph::SolidFiniteElement::local_to_lagrangian_mapping ( const DShape dpsids,
DenseMatrix< double > &  jacobian,
DenseMatrix< double > &  inverse_jacobian 
) const [inline, protected, virtual]

Calculate the mapping from local to lagrangian coordinates, given the derivatives of the shape functions w.r.t. local coorindates. Return the determinant of the jacobian, the jacobian and inverse jacobian.

Reimplemented in oomph::SolidDiagQHermiteElement< DIM >, oomph::SolidDiagQHermiteElement< 1 >, and oomph::SolidDiagQHermiteElement< 2 >.

Definition at line 2950 of file elements.h.

References assemble_local_to_lagrangian_jacobian(), and oomph::FiniteElement::invert_jacobian_mapping().

Referenced by d2shape_lagrangian(), d2shape_lagrangian_at_knot(), dshape_lagrangian(), dshape_lagrangian_at_knot(), and local_to_lagrangian_mapping().

double oomph::SolidFiniteElement::local_to_lagrangian_mapping_diagonal ( const DShape dpsids,
DenseMatrix< double > &  jacobian,
DenseMatrix< double > &  inverse_jacobian 
) const [protected, virtual]

Calculate the mapping from local to Lagrangian coordinates given the derivatives of the shape functions w.r.t the local coorindates. assuming that the coordinates are aligned in the direction of the local coordinates, i.e. there are no cross terms and the jacobian is diagonal. This function returns the determinant of the jacobian, the jacobian and the inverse jacobian.

Calculate the mapping from local to lagrangian coordinates assuming that the coordinates are aligned in the direction of the local coordinates, i.e. there are no cross terms and the jacobian is diagonal. The local derivatives are passed as dpsids and the jacobian and inverse jacobian are returned.

Reimplemented in oomph::RefineableSolidElement.

Definition at line 4269 of file elements.cc.

References oomph::FiniteElement::check_jacobian(), oomph::FiniteElement::dim(), nnodal_lagrangian_type(), oomph::FiniteElement::nnode(), and raw_lagrangian_position_gen().

Referenced by oomph::SolidDiagQHermiteElement< 2 >::local_to_lagrangian_mapping().

double oomph::SolidFiniteElement::multiplier ( const Vector< double > &  xi  )  [inline, private]

Access to the "multiplier" for the inertia terms in the consistent determination of the initial conditions for Newmark timestepping.

Definition at line 3171 of file elements.h.

References Multiplier_fct_pt.

Referenced by fill_in_jacobian_for_newmark_accel().

MultiplierFctPt oomph::SolidFiniteElement::multiplier_fct_pt (  )  const [inline]

Access function: Pointer to multiplicator function for assignment of consistent assignement of initial conditions for Newmark scheme (const version).

Definition at line 2843 of file elements.h.

References Multiplier_fct_pt.

MultiplierFctPt& oomph::SolidFiniteElement::multiplier_fct_pt (  )  [inline]

Access function: Pointer to multiplicator function for assignment of consistent assignement of initial conditions for Newmark scheme.

Definition at line 2834 of file elements.h.

References Multiplier_fct_pt.

unsigned oomph::SolidFiniteElement::ngeom_data (  )  const [inline, virtual]

The number of geometric data affecting a SolidFiniteElemnet is the same as the number of nodes (one variable position data per node).

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::RefineableSolidElement.

Definition at line 2529 of file elements.h.

References oomph::FiniteElement::nnode().

unsigned oomph::SolidFiniteElement::nnodal_lagrangian_type (  )  const [inline]

Return the number of types of (generalised) nodal Lagrangian coordinates required to interpolate the Lagrangian coordinates in the element. (E.g. 1 for Lagrange-type elements; 2 for Hermite beam elements; 4 for Hermite shell elements). Default value is 1. Needs to be overloaded for any other element.

Definition at line 2672 of file elements.h.

References Nnodal_lagrangian_type.

Referenced by oomph::RefineableSolidElement::assemble_local_to_lagrangian_jacobian(), assemble_local_to_lagrangian_jacobian(), oomph::RefineableSolidElement::assemble_local_to_lagrangian_jacobian2(), assemble_local_to_lagrangian_jacobian2(), construct_boundary_node(), construct_node(), fill_in_generic_jacobian_for_solid_ic(), interpolated_xi(), oomph::RefineableSolidElement::local_to_lagrangian_mapping_diagonal(), and local_to_lagrangian_mapping_diagonal().

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

Broken assignment operator.

Definition at line 2522 of file elements.h.

References oomph::BrokenCopy::broken_assign().

int oomph::SolidFiniteElement::position_local_eqn ( const unsigned &  n,
const unsigned &  k,
const unsigned &  j 
) [inline, protected]

Access function that returns the local equation number that corresponds to the j-th coordinate of the k-th position-type at the n-th local node.

Definition at line 2996 of file elements.h.

References oomph::FiniteElement::nnodal_position_type(), oomph::FiniteElement::nnode(), oomph::FiniteElement::nodal_dimension(), OOMPH_EXCEPTION_LOCATION, and Position_local_eqn.

Referenced by oomph::RefineableSolidElement::assign_solid_hanging_local_eqn_numbers(), oomph::ClampedHermiteShellBoundaryConditionElement::fill_in_contribution_to_residuals(), oomph::ClampedSlidingHermiteBeamBoundaryConditionElement::fill_in_contribution_to_residuals(), oomph::AxisymmetricPVDEquations::fill_in_contribution_to_residuals_axisym_pvd(), oomph::KirchhoffLoveBeamEquations::fill_in_contribution_to_residuals_beam(), oomph::KirchhoffLoveShellEquations::fill_in_contribution_to_residuals_shell(), oomph::ImposeDisplacementByLagrangeMultiplierElement< ELEMENT >::fill_in_generic_contribution_to_residuals_displ_lagr_multiplier(), oomph::RefineablePVDEquations< DIM >::fill_in_generic_contribution_to_residuals_pvd(), fill_in_generic_jacobian_for_solid_ic(), oomph::AxisymmetricPVDEquationsWithPressure::fill_in_generic_residual_contribution_axisym_pvd_with_pressure(), oomph::RefineablePVDEquationsWithPressure< DIM >::fill_in_generic_residual_contribution_pvd_with_pressure(), fill_in_jacobian_for_newmark_accel(), oomph::RefineableSolidElement::fill_in_jacobian_from_solid_position_by_fd(), fill_in_jacobian_from_solid_position_by_fd(), oomph::FSISolidTractionElement< ELEMENT, DIM >::get_dof_numbers_for_unknowns(), oomph::PVDEquationsBase< DIM >::get_dof_numbers_for_unknowns(), oomph::ClampedHermiteShellBoundaryConditionElement::get_dof_numbers_for_unknowns(), oomph::FSIDiagHermiteShellElement::get_dof_numbers_for_unknowns(), and oomph::FSIHermiteBeamElement::get_dof_numbers_for_unknowns().

double oomph::SolidFiniteElement::raw_lagrangian_position ( const unsigned &  n,
const unsigned &  i 
) const [inline]

Return i-th Lagrangian coordinate at local node n without using the hanging representation.

Definition at line 2764 of file elements.h.

References oomph::FiniteElement::node_pt().

Referenced by oomph::FSIDiagHermiteShellElement::locate_zeta(), and oomph::FSIHermiteBeamElement::locate_zeta().

double oomph::SolidFiniteElement::raw_lagrangian_position_gen ( const unsigned &  n,
const unsigned &  k,
const unsigned &  i 
) const [inline]

Return Generalised Lagrangian coordinate at local node n. `Direction' i, `Type' k. Does not use the hanging node representation.

Definition at line 2769 of file elements.h.

References oomph::FiniteElement::node_pt().

Referenced by assemble_local_to_lagrangian_jacobian(), assemble_local_to_lagrangian_jacobian2(), oomph::KirchhoffLoveBeamEquations::fill_in_contribution_to_residuals_beam(), oomph::KirchhoffLoveShellEquations::fill_in_contribution_to_residuals_shell(), oomph::KirchhoffLoveBeamEquations::get_energy(), local_to_lagrangian_mapping_diagonal(), and oomph::HermiteBeamElement::output().

virtual void oomph::SolidFiniteElement::reset_after_solid_position_fd (  )  [inline, protected, virtual]

Function that is call after the finite differencing of the solid position 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 3107 of file elements.h.

Referenced by oomph::RefineableSolidElement::fill_in_jacobian_from_solid_position_by_fd(), and fill_in_jacobian_from_solid_position_by_fd().

virtual void oomph::SolidFiniteElement::reset_in_solid_position_fd ( const unsigned &  i  )  [inline, protected, virtual]

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

Reimplemented in oomph::FSIWallElement.

Definition at line 3117 of file elements.h.

References update_in_solid_position_fd().

Referenced by oomph::RefineableSolidElement::fill_in_jacobian_from_solid_position_by_fd(), and fill_in_jacobian_from_solid_position_by_fd().

void oomph::SolidFiniteElement::set_lagrangian_dimension ( const unsigned &  lagrangian_dimension  )  [inline]

Set the lagrangian dimension of the element --- the number of lagrangian coordinates stored at the nodes in the element.

Definition at line 2489 of file elements.h.

References Lagrangian_dimension.

Referenced by oomph::SolidQHermiteElement< DIM >::build_face_element(), and oomph::SolidQHermiteElement< 2 >::SolidQHermiteElement().

virtual void oomph::SolidFiniteElement::set_macro_elem_pt ( MacroElement macro_elem_pt,
MacroElement undeformed_macro_elem_pt 
) [inline, virtual]

Set pointers to "current" and "undeformed" MacroElements. Can be overloaded in derived classes to perform additional tasks.

Reimplemented in oomph::QSolidElementBase, oomph::RefineableSolidQElement< 3 >, and oomph::RefineableSolidQElement< 2 >.

Definition at line 2586 of file elements.h.

References oomph::FiniteElement::macro_elem_pt(), oomph::FiniteElement::Macro_elem_pt, undeformed_macro_elem_pt(), and Undeformed_macro_elem_pt.

virtual void oomph::SolidFiniteElement::set_macro_elem_pt ( MacroElement macro_elem_pt  )  [inline, virtual]

Set pointer to MacroElement -- overloads generic version and uses the MacroElement also as the default for the "undeformed" configuration. This assignment must be overwritten with set_undeformed_macro_elem_pt(...) if the deformation of the solid body is driven by a deformation of the "current" Domain/MacroElement representation of it's boundary. Can be overloaded in derived classes to perform additional tasks.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::QSolidElementBase, oomph::RefineableSolidQElement< 3 >, and oomph::RefineableSolidQElement< 2 >.

Definition at line 2577 of file elements.h.

References oomph::FiniteElement::macro_elem_pt(), oomph::FiniteElement::Macro_elem_pt, and Undeformed_macro_elem_pt.

void oomph::SolidFiniteElement::set_nnodal_lagrangian_type ( const unsigned &  nlagrangian_type  )  [inline, protected]

Set the number of types required to interpolate the Lagrangian coordinates.

Definition at line 2941 of file elements.h.

References Nnodal_lagrangian_type.

Referenced by oomph::SolidQHermiteElement< 2 >::SolidQHermiteElement().

void oomph::SolidFiniteElement::set_undeformed_macro_elem_pt ( MacroElement undeformed_macro_elem_pt  )  [inline]

Set pointer to "undeformed" macro element. Can be overloaded in derived classes to perform additional tasks.

Definition at line 2596 of file elements.h.

References undeformed_macro_elem_pt(), and Undeformed_macro_elem_pt.

Referenced by oomph::QSolidElementBase::set_macro_elem_pt().

SolidInitialCondition*& oomph::SolidFiniteElement::solid_ic_pt (  )  [inline]

Pointer to object that describes the initial condition.

Definition at line 2814 of file elements.h.

References Solid_ic_pt.

bool& oomph::SolidFiniteElement::solve_for_consistent_newmark_accel_flag (  )  [inline]

Access function to flag which indicates which system of equations to solve when assigning initial conditions for time-dependent problems. If true, solve for the history value that corresponds to the acceleration in the Newmark scheme by demanding that the PDE is satisifed at the initial time. In this case the Jacobian is replaced by the mass matrix.

Definition at line 2826 of file elements.h.

References Solve_for_consistent_newmark_accel_flag.

MacroElement* oomph::SolidFiniteElement::undeformed_macro_elem_pt (  )  [inline]

Access function to pointer to "undeformed" macro element.

Definition at line 2602 of file elements.h.

References Undeformed_macro_elem_pt.

Referenced by oomph::QSolidElementBase::set_macro_elem_pt(), set_macro_elem_pt(), and set_undeformed_macro_elem_pt().

virtual void oomph::SolidFiniteElement::update_before_solid_position_fd (  )  [inline, protected, virtual]

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

Definition at line 3102 of file elements.h.

Referenced by oomph::RefineableSolidElement::fill_in_jacobian_from_solid_position_by_fd(), and fill_in_jacobian_from_solid_position_by_fd().

virtual void oomph::SolidFiniteElement::update_in_solid_position_fd ( const unsigned &  i  )  [inline, protected, virtual]

Function called within the finite difference loop for the solid position dat after a change in any values in the n-th node.

Reimplemented in oomph::FSIWallElement.

Definition at line 3112 of file elements.h.

Referenced by oomph::RefineableSolidElement::fill_in_jacobian_from_solid_position_by_fd(), fill_in_jacobian_from_solid_position_by_fd(), and reset_in_solid_position_fd().

double oomph::SolidFiniteElement::zeta_nodal ( const unsigned &  n,
const unsigned &  k,
const unsigned &  i 
) const [inline, virtual]

In a SolidFiniteElement, the "global" intrinsic coordinate of the element when viewed as part of a compound geometric object (a Mesh) is, by default, the Lagrangian coordinate Note the assumption here is that we are always using isoparameteric elements in which the lagrangian coordinate is interpolated by the same shape functions as the eulerian coordinate.

Reimplemented from oomph::FiniteElement.

Reimplemented in oomph::SolidFaceElement, oomph::ElasticPointFluidInterfaceEdgeElement< ELEMENT >, and oomph::ElasticLineFluidInterfaceEdgeElement< ELEMENT >.

Definition at line 2542 of file elements.h.

References lagrangian_position_gen().


Member Data Documentation

unsigned oomph::SolidFiniteElement::Lagrangian_dimension [private]

The Lagrangian dimension of the nodes stored in the element,.

Definition at line 3146 of file elements.h.

Referenced by assemble_local_to_lagrangian_jacobian(), lagrangian_dimension(), and set_lagrangian_dimension().

MultiplierFctPt oomph::SolidFiniteElement::Multiplier_fct_pt [private]

Pointer to function that computes the "multiplier" for the inertia terms in the consistent determination of the initial conditions for Newmark timestepping.

Definition at line 3137 of file elements.h.

Referenced by multiplier(), and multiplier_fct_pt().

unsigned oomph::SolidFiniteElement::Nnodal_lagrangian_type [private]

The number of coordinate types requried to intepolate the Lagrangian coordinates in the element. For Lagrange elements it is 1 (the default). It must be over-ridden by using the set_nlagrangian_type() function in the constructors of elements that use generalised coordinate, e.g. for 1D Hermite elements Nnodal_position_types =2.

Definition at line 3154 of file elements.h.

Referenced by nnodal_lagrangian_type(), and set_nnodal_lagrangian_type().

int* oomph::SolidFiniteElement::Position_local_eqn [private]

Array to hold the local equation number information for the solid equations, whatever they may be.

Definition at line 3142 of file elements.h.

Referenced by assign_solid_local_eqn_numbers(), position_local_eqn(), and ~SolidFiniteElement().

SolidInitialCondition* oomph::SolidFiniteElement::Solid_ic_pt [protected]

Pointer to object that specifies the initial condition.

Definition at line 2991 of file elements.h.

Referenced by oomph::KirchhoffLoveBeamEquations::fill_in_contribution_to_residuals_beam(), oomph::RefineablePVDEquations< DIM >::fill_in_generic_contribution_to_residuals_pvd(), fill_in_generic_jacobian_for_solid_ic(), oomph::RefineablePVDEquationsWithPressure< DIM >::fill_in_generic_residual_contribution_pvd_with_pressure(), fill_in_jacobian_for_newmark_accel(), and solid_ic_pt().

bool oomph::SolidFiniteElement::Solve_for_consistent_newmark_accel_flag [protected]

Flag to indicate which system of equations to solve when assigning initial conditions for time-dependent problems. If true, solve for the history value that corresponds to the acceleration in the Newmark scheme by demanding that the PDE is satisifed at the initial time. In this case the Jacobian is replaced by the mass matrix.

Definition at line 3164 of file elements.h.

Referenced by oomph::KirchhoffLoveShellEquations::fill_in_contribution_to_jacobian(), oomph::FSIWallElement::fill_in_contribution_to_jacobian(), fill_in_contribution_to_jacobian(), oomph::KirchhoffLoveBeamEquations::fill_in_contribution_to_jacobian(), oomph::AxisymmetricPVDEquations::fill_in_contribution_to_jacobian(), fill_in_jacobian_for_newmark_accel(), and solve_for_consistent_newmark_accel_flag().

MacroElement* oomph::SolidFiniteElement::Undeformed_macro_elem_pt [protected]

Pointer to the element's "undeformed" macro element (NULL by default).

Definition at line 2945 of file elements.h.

Referenced by oomph::QSolidElementBase::get_x_and_xi(), set_macro_elem_pt(), set_undeformed_macro_elem_pt(), and undeformed_macro_elem_pt().


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