29 #ifndef OOMPH_REFINEABLE_GEN_ADVECTION_DIFFUSION_ELEMENTS_HEADER
30 #define OOMPH_REFINEABLE_GEN_ADVECTION_DIFFUSION_ELEMENTS_HEADER
34 #include <oomph-lib-config.h>
38 #include "../generic/refineable_quad_element.h"
39 #include "../generic/refineable_brick_element.h"
40 #include "../generic/error_estimator.h"
53 template<
unsigned DIM>
103 const unsigned n_node =
nnode();
118 for (
unsigned l = 0; l < n_node; l++)
120 values[0] += this->
nodal_value(l, u_nodal_index) * psi[l];
136 const unsigned n_node =
nnode();
151 for (
unsigned l = 0; l < n_node; l++)
153 values[0] += this->
nodal_value(t, l, u_nodal_index) * psi[l];
162 cast_father_element_pt =
172 this->
Pe_pt = cast_father_element_pt->
pe_pt();
197 template<
unsigned DIM,
unsigned NNODE_1D>
251 return (NNODE_1D - 1);
271 template<
unsigned DIM,
unsigned NNODE_1D>
274 :
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 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.
A class for all elements that solve the Advection Diffusion equations in conservative form using isop...
GeneralisedAdvectionDiffusionWindFctPt & wind_fct_pt()
Access function: Pointer to wind function.
double *& pe_pt()
Pointer to Peclet number.
double * Pe_pt
Pointer to global Peclet number.
GeneralisedAdvectionDiffusionDiffFctPt Diff_fct_pt
Pointer to diffusivity funciton.
double *& pe_st_pt()
Pointer to Peclet number multipled by Strouha number.
GeneralisedAdvectionDiffusionWindFctPt Wind_fct_pt
Pointer to wind function:
GeneralisedAdvectionDiffusionWindFctPt Conserved_wind_fct_pt
Pointer to additional (conservative) wind function:
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: .
GeneralisedAdvectionDiffusionSourceFctPt Source_fct_pt
Pointer to source function:
double * PeSt_pt
Pointer to global Peclet number multiplied by Strouhal number.
virtual unsigned u_index_cons_adv_diff() const
Return the index at which the unknown value is stored. The default value, 0, is appropriate for singl...
GeneralisedAdvectionDiffusionDiffFctPt & diff_fct_pt()
Access function: Pointer to diffusion function.
GeneralisedAdvectionDiffusionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when time-derivatives are computed....
GeneralisedAdvectionDiffusionWindFctPt & conserved_wind_fct_pt()
Access function: Pointer to additional (conservative) wind 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.
A version of the GeneralisedAdvection Diffusion equations that can be used with non-uniform mesh refi...
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: Standard flux.from GeneralisedAdvectionDiffusion equations.
void operator=(const RefineableGeneralisedAdvectionDiffusionEquations< DIM > &)=delete
Broken assignment operator.
RefineableGeneralisedAdvectionDiffusionEquations()
Empty Constructor.
RefineableGeneralisedAdvectionDiffusionEquations(const RefineableGeneralisedAdvectionDiffusionEquations< DIM > &dummy)=delete
Broken copy constructor.
void fill_in_generic_residual_contribution_cons_adv_diff(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Add the element's contribution to the elemental residual vector and/or Jacobian matrix flag=1: comput...
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 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...
A class that is used to template the refineable Q elements by dimension. It's really nothing more tha...
Refineable version of QGeneralisedAdvectionDiffusionElement. Inherit from the standard QGeneralisedAd...
void operator=(const RefineableQGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D > &)=delete
Broken assignment operator.
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.
RefineableQGeneralisedAdvectionDiffusionElement(const RefineableQGeneralisedAdvectionDiffusionElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
unsigned nvertex_node() const
Number of vertex nodes in the element.
RefineableQGeneralisedAdvectionDiffusionElement()
Empty Constructor:
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 1.
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes....
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...