29#ifndef OOMPH_POROELASITICTY_FACE_ELEMENTS_HEADER
30#define OOMPH_POROELASITICTY_FACE_ELEMENTS_HEADER
34#include <oomph-lib-config.h>
47 namespace PoroelasticityFaceElementHelper
85 template<
class ELEMENT>
118 const unsigned&
intpt,
132 const unsigned&
intpt,
159 ELEMENT*
elem_pt =
new ELEMENT;
168 "nodes are hanging\n",
177 Element_pt =
dynamic_cast<ELEMENT*
>(element_pt);
232 const unsigned&
i)
const
286 template<
class ELEMENT>
312 template<
class ELEMENT>
338 template<
class ELEMENT>
351 if (n_position_type != 1)
353 throw OomphLibError(
"Poroelasticity equations are not yet implemented "
354 "for more than one position type",
363 unsigned n_q_basis = Element_pt->nq_basis();
396 for (
unsigned i = 0;
i <
n_dim - 1;
i++)
417 for (
unsigned i = 0;
i <
n_dim;
i++)
424 for (
unsigned j = 0;
j <
n_dim - 1;
j++)
433 for (
unsigned i = 0;
i <
n_dim - 1;
i++)
435 for (
unsigned j = 0;
j <
n_dim - 1;
j++)
441 for (
unsigned k = 0;
k <
n_dim;
k++)
460 Adet = A(0, 0) * A(1, 1) - A(0, 1) * A(1, 0);
463 throw OomphLibError(
"Wrong dimension in PoroelasticityFaceElement",
484 for (
unsigned i = 0;
i <
n_dim;
i++)
506 for (
unsigned i = 0;
i <
n_dim;
i++)
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
FaceElements are elements that coincide with the faces of higher-dimensional "bulk" elements....
int & face_index()
Index of the face (a number that uniquely identifies the face in the element)
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
In a FaceElement, the "global" intrinsic coordinate of the element along the boundary,...
FaceGeometry class definition: This policy class is used to allow construction of face elements that ...
A general Finite Element class.
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
virtual 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...
virtual void output(std::ostream &outfile)
Output the element data — typically the values at the nodes in a format suitable for post-processing.
unsigned nnodal_position_type() const
Return the number of coordinate types that the element requires to interpolate the geometry between t...
double size() const
Calculate the size of the element (length, area, volume,...) in Eulerian computational coordinates....
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
unsigned nnode() const
Return the number of nodes.
double nodal_position(const unsigned &n, const unsigned &i) const
Return the i-th coordinate at local node n. If the node is hanging, the appropriate interpolation is ...
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
virtual void build_face_element(const int &face_index, FaceElement *face_element_pt)
Function for building a lower dimensional FaceElement on the specified face of the FiniteElement....
unsigned nodal_dimension() const
Return the required Eulerian dimension of the nodes in this element.
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
An OomphLibError object which should be thrown when an run-time error is encountered....
An OomphLibWarning object which should be created as a temporary object to issue a warning....
A class for elements that allow the imposition of an applied pressure in the Darcy equations....
void(*&)(const double &time, const Vector< double > &x, const Vector< double > &n, double &pressure) pressure_fct_pt()
Reference to the pressure function pointer.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Fill in contribution from Jacobian.
void traction(const double &time, const Vector< double > &s, Vector< double > &traction)
Compute traction vector at specified local coordinate Should only be used for post-processing; ignore...
PoroelasticityFaceElement(FiniteElement *const &element_pt, const int &face_index)
Constructor, which takes a "bulk" element and the value of the index and its limit.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Return the residuals.
void(* Traction_fct_pt)(const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &result)
Pointer to an imposed traction function. Arguments: Eulerian coordinate; outer unit normal; applied t...
void output(std::ostream &outfile)
Output function.
void(*&)(const double &time, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction) traction_fct_pt()
Reference to the traction function pointer.
ELEMENT * Element_pt
pointer to the bulk element this face element is attached to
void fill_in_contribution_to_residuals_darcy_face(Vector< double > &residuals)
Helper function that actually calculates the residuals.
void output(FILE *file_pt)
C_style output function.
virtual void get_traction(const double &time, const unsigned &intpt, const Vector< double > &x, const Vector< double > &n, Vector< double > &traction)
Get the traction vector: Pass number of integration point (dummy), Eulerlian coordinate and normal ve...
virtual void get_pressure(const double &time, const unsigned &intpt, const Vector< double > &x, const Vector< double > &n, double &pressure)
Get the pressure value: Pass number of integration point (dummy), Eulerlian coordinate and normal vec...
void pressure(const double &time, const Vector< double > &s, double &pressure)
Compute pressure value at specified local coordinate Should only be used for post-processing; ignores...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function.
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Specify the value of nodal zeta from the face geometry The "global" intrinsic coordinate of the eleme...
void(* Pressure_fct_pt)(const double &time, const Vector< double > &x, const Vector< double > &n, double &result)
Pointer to an imposed pressure function. Arguments: Eulerian coordinate; outer unit normal; applied p...
RefineableElements are FiniteElements that may be subdivided into children to provide a better local ...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
TAdvectionDiffusionReactionElement()
Constructor: Call constructors for TElement and AdvectionDiffusionReaction equations.
Time *const & time_pt() const
Access function for the pointer to time (const version)
double & time()
Return the current value of the continuous time.
void Zero_traction_fct(const double &time, const Vector< double > &x, const Vector< double > &N, Vector< double > &load)
Default load function (zero traction)
void Zero_pressure_fct(const double &time, const Vector< double > &x, const Vector< double > &N, double &load)
Default load function (zero pressure)
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).