27 #ifndef OOMPH_REFINEABLE_SPACE_TIME_UNSTEADY_HEAT_ELEMENTS_HEADER
28 #define OOMPH_REFINEABLE_SPACE_TIME_UNSTEADY_HEAT_ELEMENTS_HEADER
32 #include <oomph-lib-config.h>
50 template<
unsigned SPATIAL_DIM>
75 return SPATIAL_DIM + 1;
85 unsigned n_node =
nnode();
94 DShape dpsidx(n_node, SPATIAL_DIM + 1);
100 for (
unsigned j = 0; j < SPATIAL_DIM + 1; j++)
107 for (
unsigned l = 0; l < n_node; l++)
110 for (
unsigned j = 0; j < SPATIAL_DIM + 1; j++)
113 flux[j] += this->
nodal_value(l, u_nodal_index) * dpsidx(l, j);
130 unsigned n_node =
nnode();
145 for (
unsigned l = 0; l < n_node; l++)
148 values[0] += this->
nodal_value(l, u_nodal_index) * psi[l];
165 unsigned n_node =
nnode();
180 for (
unsigned l = 0; l < n_node; l++)
183 values[0] += this->
nodal_value(t, l, u_nodal_index) * psi[l];
193 cast_father_element_pt =
212 const unsigned& flag);
219 template<
unsigned SPATIAL_DIM,
unsigned NNODE_1D>
281 return (NNODE_1D - 1);
298 template<
unsigned SPATIAL_DIM,
unsigned NNODE_1D>
301 :
public virtual QElement<SPATIAL_DIM, NNODE_1D>
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
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 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.
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Compute the geometric shape functions and also first derivatives w.r.t. global coordinates at local c...
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.
A class that is used to template the refineable Q elements by dimension. It's really nothing more tha...
Refineable version of 2D QUnsteadyHeatSpaceTimeElement elements.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 1.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons (empty)
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes (empt...
unsigned nvertex_node() const
Number of vertex nodes in the element.
RefineableQUnsteadyHeatSpaceTimeElement()
Constructor.
RefineableQUnsteadyHeatSpaceTimeElement(const RefineableQUnsteadyHeatSpaceTimeElement< SPATIAL_DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
Refineable version of Unsteady Heat equations.
RefineableSpaceTimeUnsteadyHeatEquations()
Constructor.
virtual void fill_in_generic_residual_contribution_ust_heat(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Add element's contribution to elemental residual vector and/or Jacobian matrix flag=0: compute residu...
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...
RefineableSpaceTimeUnsteadyHeatEquations(const RefineableSpaceTimeUnsteadyHeatEquations< SPATIAL_DIM > &dummy)=delete
Broken copy constructor.
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 get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Different to the get_flux function in the base class as we also hav...
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...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
A class for all isoparametric elements that solve the SpaceTimeUnsteadyHeat equations.
virtual unsigned u_index_ust_heat() const
Return the index at which the unknown value is stored. The default value, 0, is appropriate for singl...
SpaceTimeUnsteadyHeatSourceFctPt Source_fct_pt
Pointer to source function:
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when time-derivatives are computed....
SpaceTimeUnsteadyHeatSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...