30#ifndef OOMPH_STORED_SHAPE_FUNCTION_ELEMENTS_HEADER
31#define OOMPH_STORED_SHAPE_FUNCTION_ELEMENTS_HEADER
36#include <oomph-lib-config.h>
45#define OOMPH_STORED_SHAPE_FUNCTIONS_VERBOSE
46#undef OOMPH_STORED_SHAPE_FUNCTIONS_VERBOSE
553 template<
class ELEMENT>
555 public virtual ELEMENT
590 template<
class ELEMENT>
593 public virtual ELEMENT
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
A general Finite Element class.
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Generic class for numerical integration schemes:
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
SolidFiniteElement class.
Base class for elements that allow storage of precomputed shape functions and their derivatives w....
Vector< DShape * > *& dshape_eulerian_stored_pt()
Return a pointer to the vector of pointers to the stored first derivatives of the shape functions w....
Vector< Shape * > *const & shape_stored_pt() const
Return a pointer to the vector of pointers to the stored shape functions (const version)
void shape_at_knot(const unsigned &ipt, Shape &psi) const
Return the geometric shape function at the ipt-th integration point.
virtual void set_integration_scheme(Integral *const &integral_pt)
Set the spatial integration scheme – overloaded from the finite element base class since a change in ...
Shape *& shape_stored_pt(const unsigned &ipt)
Return a pointer to the stored shape function at the ipt-th integration point.
void delete_all_dshape_eulerian_stored()
Delete all storage related to deriatives of shape fcts w.r.t. to global Eulerian coords.
Vector< DShape * > *& d2shape_eulerian_stored_pt()
Return a pointer to the vector of pointers to the stored second derivatives of the shape functions w....
Vector< DShape * > * DShape_local_stored_pt
Pointer to storage for the pointers to the derivatives of the nodal shape functions w....
Vector< DShape * > * D2Shape_local_stored_pt
Pointer to storage for the pointers to the second derivatives of the nodal shape functions w....
virtual ~StorableShapeElementBase()
The destructor cleans up the static memory allocated for shape function storage. Internal and externa...
Vector< Shape * > *& shape_stored_pt()
Return a pointer to the vector of pointers to the stored shape functions.
Vector< DShape * > *const & d2shape_local_stored_pt() const
Return a pointer to the vector of pointers to the stored second derivatives of the shape functions w....
Vector< DShape * > * D2Shape_eulerian_stored_pt
Pointer to storage for the second derivatives of the shape functions w.r.t. global coordinates at int...
void delete_dshape_local_stored()
Delete stored derivatives of shape functions w.r.t. to local coordinates.
void delete_shape_local_stored()
Delete stored shape functions.
Vector< DShape * > *const & dshape_eulerian_stored_pt() const
Return a pointer to the vector of pointers to the stored first derivatives of the shape functions w....
void delete_d2shape_local_stored()
Delete stored 2nd derivatives of shape functions w.r.t. to local coordinates.
StorableShapeElementBase(const StorableShapeElementBase &)=delete
Broken copy constructor.
double dshape_eulerian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidx) const
Return the geometric shape functions and also first derivatives w.r.t. global coordinates at the ipt-...
Vector< Shape * > * Shape_stored_pt
Pointer to storage for the pointers to the nodal shape functions at the integration points (knots)
Vector< double > *& jacobian_eulerian_stored_pt()
Return a pointer to the vector of Jacobians of the mapping between the local and global (eulerian) co...
void pre_compute_d2shape_eulerian_at_knots()
Calculate the first and second derivatives of the shape functions w.r.t global coordinates at the int...
Vector< double > * Jacobian_eulerian_stored_pt
Pointer to storage for the Jacobian of the element w.r.t global coordinates.
void delete_d2shape_eulerian_stored()
Delete stored second derivatives of shape functions w.r.t. to global Eulerian coordinates.
void dshape_local_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsids) const
Return the geometric shape function and its derivative w.r.t. the local coordinates at the ipt-th int...
Vector< DShape * > *const & dshape_local_stored_pt() const
Return a pointer to the vector of pointers to the stored first derivatives of the shape functions w....
void pre_compute_dshape_eulerian_at_knots()
Calculate the first derivatives of the shape functions w.r.t the global coordinates at the integratio...
void set_shape_local_stored_from_element(StorableShapeElementBase *const &element_pt)
Set the shape functions pointed to internally to be those pointed to by the FiniteElement element_pt ...
Vector< DShape * > *const & d2shape_eulerian_stored_pt() const
Return a pointer to the vector of pointers to the stored second derivatives of the shape functions w....
void pre_compute_J_eulerian_at_knots()
Calculate the Jacobian of the mapping from local to global coordinates at the integration points and ...
double d2shape_eulerian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidx, DShape &d2psidx) const
Return the geometric shape functions and also first and second derivatives w.r.t. global coordinates ...
void set_dshape_eulerian_stored_from_element(StorableShapeElementBase *const &element_pt)
Set the derivatives of stored shape functions with respect to the global coordinates to be the same a...
Vector< double > *const & jacobian_eulerian_stored_pt() const
Return a pointer to the vector of Jacobians of the mapping between the local and global (eulerian) co...
Vector< DShape * > *& d2shape_local_stored_pt()
Return a pointer to the vector of pointers to the stored second derivatives of the shape functions w....
void delete_J_eulerian_stored()
Delete stored Jacobian of mapping between local and global (Eulerian) coordinates.
void pre_compute_dshape_local_at_knots()
Calculate the shape functions and first derivatives w.r.t. local coordinatess at the integration poin...
void delete_dshape_eulerian_stored()
Delete stored deriatives of shape fcts w.r.t. to global Eulerian coords.
bool Can_delete_dshape_eulerian_stored
Boolean to determine whether the element can delete the stored derivatives of shape functions w....
StorableShapeElementBase()
Constructor, set most storage pointers to NULL.
bool Can_delete_shape_local_stored
Boolean to determine whether the element can delete the stored local shape functions.
void pre_compute_shape_at_knots()
Calculate the shape functions at the integration points and store the results internally.
double J_eulerian_at_knot(const unsigned &ipt) const
Return the Jacobian of the mapping from local to global coordinates at the ipt-th integration point.
Vector< DShape * > *& dshape_local_stored_pt()
Return a pointer to the vector of pointers to the stored first derivatives of the shape functions w....
Vector< DShape * > * DShape_eulerian_stored_pt
Pointer to storage for the derivatives of the shape functions w.r.t. global coordinates at integratio...
void delete_all_shape_local_stored()
Delete all the objects stored in the [...]_local_stored_pt vectors and delete the vectors themselves.
void pre_compute_d2shape_local_at_knots()
Calculate the second derivatives of the shape functions w.r.t. local coordinates at the integration p...
Shape *const & shape_stored_pt(const unsigned &ipt) const
Return a pointer to the stored shape function at the ipt-th integration point (const version)
void d2shape_local_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsids, DShape &d2psids) const
Return the geometric shape function and its first and second derivatives w.r.t. the local coordinates...
void operator=(const StorableShapeElementBase &)=delete
Broken assignment operator.
Templated wrapper that attaches the ability to store the shape functions and their derivatives w....
StorableShapeElement(const StorableShapeElement &)=delete
Broken copy constructor.
void operator=(const StorableShapeElement &)=delete
Broken assignment operator.
StorableShapeElement()
Constructor, set most storage pointers to zero.
virtual ~StorableShapeElement()
Empty virtual destructor.
Base class for solid elements that allow storage of precomputed shape functions and their derivatives...
Vector< double > * Jacobian_lagrangian_stored_pt
Pointer to storage for the Jacobian of the mapping between the local and the global Lagrangian coordi...
Vector< DShape * > *const & d2shape_lagrangian_stored_pt() const
Return a pointer to the vector of pointers to the stored second derivatives of the shape functions w....
virtual ~StorableShapeSolidElementBase()
Destructor to clean up any allocated memory.
Vector< DShape * > *& dshape_lagrangian_stored_pt()
Return a pointer to the vector of pointers to the stored first derivatives of the shape functions w....
void pre_compute_d2shape_lagrangian_at_knots()
Calculate the first and second derivatives of the shape functions w.r.t Lagrangian coordinates at the...
Vector< DShape * > *const & dshape_lagrangian_stored_pt() const
Return a pointer to the vector of pointers to the stored first derivatives of the shape functions w....
void pre_compute_dshape_lagrangian_at_knots()
Calculate the first derivatives of the shape functions w.r.t Lagrangian coordinates at the integratio...
double d2shape_lagrangian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidxi, DShape &d2psidxi) const
Return the geometric shape functions and also first and second derivatives w.r.t. Lagrangian coordina...
StorableShapeSolidElementBase(const StorableShapeSolidElementBase &)=delete
Broken copy constructor.
Vector< DShape * > * D2Shape_lagrangian_stored_pt
Pointer to storage for the pointers to the second derivatives of the shape functions w....
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-...
Vector< double > *const & jacobian_lagrangian_stored_pt() const
Return a pointer to the vector of Jacobians of the mapping between the local and global (eulerian) co...
void delete_J_lagrangian_stored()
Delete stored Jaocbian of mapping between local and Lagrangian coordinates.
StorableShapeSolidElementBase()
Constructor: Set defaults: Nothing is stored.
Vector< DShape * > *& d2shape_lagrangian_stored_pt()
Return a pointer to the vector of pointers to the stored second derivatives of the shape functions w....
void delete_d2shape_lagrangian_stored()
Delete stored second derivatives of shape functions w.r.t. Lagrangian coordinates.
bool Can_delete_dshape_lagrangian_stored
Boolean to determine whether the element can delete the stored shape function derivatives w....
void set_integration_scheme(Integral *const &integral_pt)
Overload the set_integration_scheme to recompute any stored derivatives w.r.t. Lagrangian coordinates...
void delete_all_dshape_lagrangian_stored()
Delete all the objects stored in the [...]_lagrangian_stored_pt vectors and delete the vectors themse...
void operator=(const StorableShapeSolidElementBase &)=delete
Broken assignment operator.
void delete_dshape_lagrangian_stored()
Delete all the objects stored in the Lagrangian_stored vectors.
Vector< DShape * > * DShape_lagrangian_stored_pt
Pointer to storage for the pointers to the derivatives of the shape functions w.r....
Vector< double > *& jacobian_lagrangian_stored_pt()
Return a pointer to the vector of Jacobians of the mapping between the local and global (eulerian) co...
void set_dshape_lagrangian_stored_from_element(StorableShapeSolidElementBase *const &element_pt)
Set the derivatives of stored shape functions with respect to the lagrangian coordinates to be the sa...
Templated wrapper that attaches the ability to store the shape functions and their derivatives w....
void operator=(const StorableShapeSolidElement &)=delete
Broken assignment operator.
StorableShapeSolidElement(const StorableShapeSolidElement &)=delete
Broken copy constructor.
StorableShapeSolidElement()
Constructor: Set defaults.
virtual ~StorableShapeSolidElement()
Destructor to clean up any allocated memory.
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).