29#ifndef OOMPH_REFINEABLE_TIME_HARMONIC_LINEAR_ELASTICITY_ELEMENTS_HEADER
30#define OOMPH_REFINEABLE_TIME_HARMONIC_LINEAR_ELASTICITY_ELEMENTS_HEADER
43 template<
unsigned DIM>
68 values.resize(2 *
DIM, 0.0);
78 for (
unsigned i = 0;
i <
DIM;
i++)
117 std::ostringstream error_message;
118 error_message <<
"The flux vector has the wrong number of entries, "
135 for (
unsigned i = 0;
i <
DIM;
i++)
144 for (
unsigned i = 0;
i <
DIM;
i++)
146 for (
unsigned j =
i + 1;
j <
DIM;
j++)
175 cast_father_element_pt->elasticity_tensor_pt();
178 this->
Omega_sq_pt = cast_father_element_pt->omega_sq_pt();
196 template<
unsigned DIM,
unsigned NNODE_1D>
245 template<
unsigned NNODE_1D>
248 :
public virtual QElement<1, NNODE_1D>
260 template<
unsigned NNODE_1D>
276 template<
unsigned NNODE_1D>
279 :
public virtual QElement<2, NNODE_1D>
291 template<
unsigned NNODE_1D>
294 :
public virtual QElement<1, NNODE_1D>
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
FaceGeometry class definition: This policy class is used to allow construction of face elements that ...
double nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n. Produces suitably interpolated values for hanging nodes...
double size() const
Calculate the size of the element (length, area, volume,...) in Eulerian computational coordinates....
virtual void shape(const Vector< double > &s, Shape &psi) const =0
Calculate the geometric shape functions at local coordinate s. This function must be overloaded for e...
unsigned nnode() const
Return the number of nodes.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
An OomphLibError object which should be thrown when an run-time error is encountered....
Point element has just a single node and a single shape function which is identically equal to one.
An Element that solves the equations of linear elasticity in Cartesian coordinates,...
RefineableElements are FiniteElements that may be subdivided into children to provide a better local ...
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
A class that is used to template the refineable Q elements by dimension. It's really nothing more tha...
Class for refineable QTimeHarmonicLinearElasticityElement elements.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
void further_setup_hanging_nodes()
No additional hanging node procedures are required.
RefineableQTimeHarmonicLinearElasticityElement()
Constructor:
unsigned nvertex_node() const
Number of vertex nodes in the element.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
void rebuild_from_sons(Mesh *&mesh_pt)
Empty rebuild from sons, no need to reconstruct anything here.
Class for Refineable TimeHarmonicLinearElasticity equations.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 2*DIM.
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Upper triangular entries in strain tensor.
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Get the function value u in Vector. Note: Given the generality of the interface (this function is usu...
void further_build()
Further build function, pass the pointers down to the sons.
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
RefineableTimeHarmonicLinearElasticityEquations()
Constructor.
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Get the current interpolated values (displacements). Note: Given the generality of the interface (thi...
void fill_in_generic_contribution_to_residuals_time_harmonic_linear_elasticity(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Overloaded helper function to take hanging nodes into account.
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...
double * Omega_sq_pt
Square of nondim frequency.
TimeHarmonicElasticityTensor * Elasticity_tensor_pt
Pointer to the elasticity tensor.
BodyForceFctPt Body_force_fct_pt
Pointer to body force function.
virtual std::complex< unsigned > u_index_time_harmonic_linear_elasticity(const unsigned i) const
Return the index at which the i-th real or imag unknown displacement component is stored....
void get_strain(const Vector< double > &s, DenseMatrix< std::complex< double > > &strain) const
Return the strain tensor.
A class for elements that solve the equations of linear elasticity in cartesian coordinates.
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).