29 #ifndef OOMPH_REFINEABLE_LINEAR_ELASTICITY_ELEMENTS_HEADER
30 #define OOMPH_REFINEABLE_LINEAR_ELASTICITY_ELEMENTS_HEADER
34 #include "../generic/refineable_quad_element.h"
35 #include "../generic/refineable_brick_element.h"
36 #include "../generic/error_estimator.h"
43 template<
unsigned DIM>
68 values.resize(DIM, 0.0);
71 unsigned n_node = this->
nnode();
78 for (
unsigned i = 0;
i < DIM;
i++)
82 for (
unsigned l = 0; l < n_node; l++)
84 values[
i] += this->
nodal_value(t, l, u_nodal_index) * psi(l);
104 return DIM + DIM * (DIM - 1) / 2;
112 unsigned num_entries = DIM + ((DIM * DIM) - DIM) / 2;
113 if (flux.size() != num_entries)
115 std::ostringstream error_message;
116 error_message <<
"The flux vector has the wrong number of entries, "
117 << flux.size() <<
", whereas it should be " << num_entries
120 OOMPH_CURRENT_FUNCTION,
121 OOMPH_EXCEPTION_LOCATION);
133 for (
unsigned i = 0;
i < DIM;
i++)
135 flux[icount] = strain(
i,
i);
140 for (
unsigned i = 0;
i < DIM;
i++)
142 for (
unsigned j =
i + 1; j < DIM; j++)
144 flux[icount] = strain(
i, j);
192 template<
unsigned DIM,
unsigned NNODE_1D>
238 template<
unsigned DIM>
308 std::ostream& outfile,
318 template<
unsigned NNODE_1D>
320 :
public virtual QElement<1, NNODE_1D>
332 template<
unsigned NNODE_1D>
347 template<
unsigned NNODE_1D>
349 :
public virtual QElement<2, NNODE_1D>
361 template<
unsigned NNODE_1D>
364 :
public virtual QElement<1, NNODE_1D>
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
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...
virtual unsigned nvertex_node() const
Return the number of vertex nodes in this element. Broken virtual function in "pure" finite elements.
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.
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
virtual Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element. Broken virtual function in "pure" finite elements.
virtual void set_integration_scheme(Integral *const &integral_pt)
Set the spatial integration scheme.
Class for multidimensional Gauss Lobatto Legendre integration rules empty - just establishes template...
ElasticityTensor * Elasticity_tensor_pt
Pointer to the elasticity tensor.
double *& lambda_sq_pt()
Access function for pointer to timescale ratio (nondim density)
ElasticityTensor *& elasticity_tensor_pt()
Return the pointer to the elasticity_tensor.
double * Lambda_sq_pt
Timescale ratio (non-dim. density)
void get_strain(const Vector< double > &s, DenseMatrix< double > &strain) const
Return the strain tensor.
BodyForceFctPt Body_force_fct_pt
Pointer to body force function.
bool is_inertia_enabled() const
Access function to flag that switches inertia on/off (const version)
void interpolated_u_linear_elasticity(const Vector< double > &s, Vector< double > &disp) const
Compute vector of FE interpolated displacement u at local coordinate s.
bool Unsteady
Flag that switches inertia on/off.
BodyForceFctPt & body_force_fct_pt()
Access function: Pointer to body force function.
virtual unsigned u_index_linear_elasticity(const unsigned i) const
Return the index at which the i-th unknown displacement component is stored. The default value,...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
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....
A class that is used to template the p-refineable Q elements by dimension. It's really nothing more t...
p-refineable version of 2D QLinearElasticityElement elements
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
PRefineableQLinearElasticityElement()
Constructor, simply call the other constructors.
PRefineableQLinearElasticityElement(const PRefineableQLinearElasticityElement< DIM > &dummy)=delete
Broken copy constructor.
unsigned nvertex_node() const
Number of vertex nodes in the element.
void further_build()
Broken assignment operator.
void compute_energy_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_grad_pt, double &error, double &norm)
Get error against and norm of exact solution.
~PRefineableQLinearElasticityElement()
Destructor (to avoid memory leaks)
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 1.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation:
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
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.
Class for Refineable LinearElasticity equations.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: DIM.
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.
RefineableLinearElasticityEquations()
Constructor.
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 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_linear_elasticity(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Overloaded helper function to take hanging nodes into account.
A class that is used to template the refineable Q elements by dimension. It's really nothing more tha...
/////////////////////////////////////////////////////////////////// /////////////////////////////////...
void further_setup_hanging_nodes()
No additional hanging node procedures are required.
void rebuild_from_sons(Mesh *&mesh_pt)
Empty rebuild from sons, no need to reconstruct anything here.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
RefineableQLinearElasticityElement()
Constructor:
unsigned nvertex_node() const
Number of vertex nodes in the element.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...