29 #ifndef OOMPH_REFINEABLE_ADVECTION_DIFFUSION_REACTION_ELEMENTS_HEADER
30 #define OOMPH_REFINEABLE_ADVECTION_DIFFUSION_REACTION_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 NREAGENT,
unsigned DIM>
82 return NREAGENT * DIM;
101 values.resize(NREAGENT);
104 unsigned n_node =
nnode();
113 for (
unsigned r = 0; r < NREAGENT; r++)
121 for (
unsigned l = 0; l < n_node; l++)
123 values[r] += this->
nodal_value(l, c_nodal_index) * psi[l];
137 values.resize(NREAGENT);
140 const unsigned n_node =
nnode();
147 for (
unsigned r = 0; r < NREAGENT; r++)
156 for (
unsigned l = 0; l < n_node; l++)
158 values[r] += this->
nodal_value(t, l, c_nodal_index) * psi[l];
169 cast_father_element_pt =
dynamic_cast<
205 template<
unsigned NREAGENT,
unsigned DIM,
unsigned NNODE_1D>
234 NNODE_1D>&) =
delete;
263 return (NNODE_1D - 1);
282 template<
unsigned NREAGENT,
unsigned DIM,
unsigned NNODE_1D>
285 :
public virtual QElement<DIM - 1, NNODE_1D>
A class for all elements that solve the Advection Diffusion Reaction equations using isoparametric el...
AdvectionDiffusionReactionSourceFctPt 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....
AdvectionDiffusionReactionReactionDerivFctPt Reaction_deriv_fct_pt
Pointer to reaction derivatives.
AdvectionDiffusionReactionReactionDerivFctPt & reaction_deriv_fct_pt()
Access function: Pointer to reaction derivatives function.
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: .
AdvectionDiffusionReactionReactionFctPt & reaction_fct_pt()
Access function: Pointer to reaction function.
Vector< double > * Tau_pt
Pointer to global timescales.
AdvectionDiffusionReactionWindFctPt & wind_fct_pt()
Access function: Pointer to wind function.
Vector< double > *& diff_pt()
Pointer to vector of diffusion coefficients.
AdvectionDiffusionReactionWindFctPt Wind_fct_pt
Pointer to wind function:
AdvectionDiffusionReactionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
virtual unsigned c_index_adv_diff_react(const unsigned &i) const
Return the index at which the unknown i-th reagent is stored. The default value, i,...
Vector< double > *& tau_pt()
Pointer to vector of dimensionless timescales.
AdvectionDiffusionReactionReactionFctPt Reaction_fct_pt
Pointer to reaction function.
Vector< double > * Diff_pt
Pointer to global diffusion coefficients.
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.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
//////////////////////////////////////////////////////////////////////// ////////////////////////////...
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
A version of the Advection Diffusion Reaction equations that can be used with non-uniform mesh refine...
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
RefineableAdvectionDiffusionReactionEquations(const RefineableAdvectionDiffusionReactionEquations< NREAGENT, DIM > &dummy)=delete
Broken copy constructor.
void operator=(const RefineableAdvectionDiffusionReactionEquations< NREAGENT, DIM > &)=delete
Broken assignment operator.
RefineableAdvectionDiffusionReactionEquations()
Empty Constructor.
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Get the function values c in Vector. Note: Given the generality of the interface (this function is us...
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Standard flux.from AdvectionDiffusionReaction equations.
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Get the function values c in Vector. Note: Given the generality of the interface (this function is us...
void fill_in_generic_residual_contribution_adv_diff_react(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 further_build()
Further build: Copy all pointers from the father element.
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 QAdvectionDiffusionReactionElement. Inherit from the standard QAdvectionDiffusi...
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: NREAGENT.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes....
unsigned nvertex_node() const
Number of vertex nodes in the element.
void operator=(const RefineableQAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D > &)=delete
Broken assignment operator.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
RefineableQAdvectionDiffusionReactionElement()
Empty Constructor:
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
RefineableQAdvectionDiffusionReactionElement(const RefineableQAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
A class that is used to template the refineable Q elements by dimension. It's really nothing more tha...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...