27 #ifndef OOMPH_REFINEABLE_LINEAR_WAVE_ELEMENTS_HEADER
28 #define OOMPH_REFINEABLE_LINEAR_WAVE_ELEMENTS_HEADER
32 #include <oomph-lib-config.h>
36 #include "../generic/refineable_quad_element.h"
37 #include "../generic/refineable_brick_element.h"
38 #include "../generic/error_estimator.h"
46 template<
unsigned DIM>
94 unsigned n_node =
nnode();
109 for (
unsigned l = 0; l < n_node; l++)
111 values[0] += this->
nodal_value(l, u_nodal_index) * psi[l];
131 unsigned n_node =
nnode();
141 for (
unsigned l = 0; l < n_node; l++)
143 values[0] += this->
nodal_value(t, l, u_nodal_index) * psi[l];
172 template<
unsigned DIM,
unsigned NNODE_1D>
221 return (NNODE_1D - 1);
238 template<
unsigned DIM,
unsigned NNODE_1D>
240 :
public virtual QElement<DIM - 1, NNODE_1D>
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
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.
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.
A class for all isoparametric elements that solve the LinearWave equations.
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: flux[i] = du/dx_i.
virtual unsigned u_index_lin_wave() const
Return the index at which the unknown value is stored. The default value, 0, is appropriate for singl...
LinearWaveSourceFctPt Source_fct_pt
Pointer to source function:
LinearWaveSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
//////////////////////////////////////////////////////////////////////// ////////////////////////////...
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.
Refineable version of LinearWave equations.
void fill_in_generic_residual_contribution_lin_wave(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Add element's contribution to elemental residual vector and/or Jacobian matrix flag=1: compute both f...
void get_interpolated_values(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...
RefineableLinearWaveEquations()
Constructor.
RefineableLinearWaveEquations(const RefineableLinearWaveEquations< DIM > &dummy)=delete
Broken copy constructor.
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...
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
void further_build()
Further build: Copy source function pointer from father element.
void operator=(const RefineableLinearWaveEquations< DIM > &)=delete
Broken assignment operator.
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Standard flux.from LinearWave equations.
A class that is used to template the refineable Q elements by dimension. It's really nothing more tha...
Refineable version of 2D QLinearWaveElement elements.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes....
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
unsigned nvertex_node() const
Number of vertex nodes in the element.
void operator=(const RefineableQLinearWaveElement< DIM, NNODE_1D > &)=delete
Broken assignment operator.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 1.
RefineableQLinearWaveElement(const RefineableQLinearWaveElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
RefineableQLinearWaveElement()
Constructor.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...