33 #ifndef OOMPH_REFINEABLE_ADVECTION_DIFFUSION_REACTION_ELEMENTS_HEADER 34 #define OOMPH_REFINEABLE_ADVECTION_DIFFUSION_REACTION_ELEMENTS_HEADER 38 #include <oomph-lib-config.h> 42 #include "../generic/refineable_quad_element.h" 43 #include "../generic/refineable_brick_element.h" 44 #include "../generic/error_estimator.h" 58 template <
unsigned NREAGENT,
unsigned DIM>
104 values.resize(NREAGENT);
107 unsigned n_node =
nnode();
116 for(
unsigned r=0;r<NREAGENT;r++)
124 for(
unsigned l=0;l<n_node;l++)
126 values[r] += this->
nodal_value(l,c_nodal_index)*psi[l];
139 values.resize(NREAGENT);
142 const unsigned n_node =
nnode();
149 for(
unsigned r=0;r<NREAGENT;r++)
158 for(
unsigned l=0;l<n_node;l++)
160 values[r] += this->
nodal_value(t,l,c_nodal_index)*psi[l];
171 cast_father_element_pt
208 template <
unsigned NREAGENT,
unsigned DIM,
unsigned NNODE_1D>
278 template<
unsigned NREAGENT,
unsigned DIM,
unsigned NNODE_1D>
281 public virtual QElement<DIM-1,NNODE_1D>
A version of the Advection Diffusion Reaction equations that can be used with non-uniform mesh refine...
void broken_copy(const std::string &class_name)
Issue error message and terminate execution.
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
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.
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
void operator=(const RefineableQAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D > &)
Broken assignment operator.
Refineable version of QAdvectionDiffusionReactionElement. Inherit from the standard QAdvectionDiffusi...
void operator=(const RefineableAdvectionDiffusionReactionEquations< NREAGENT, DIM > &)
Broken assignment operator.
QAdvectionDiffusionReactionElement elements are linear/quadrilateral/brick-shaped Advection Diffusion...
AdvectionDiffusionReactionReactionDerivFctPt Reaction_deriv_fct_pt
Pointer to reaction derivatives.
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when time-derivatives are computed...
AdvectionDiffusionReactionReactionFctPt & reaction_fct_pt()
Access function: Pointer to reaction function.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
RefineableAdvectionDiffusionReactionEquations()
Empty Constructor.
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, is appropriate for single-physics problems, when there are only i variables, the values that satisfies the advection-diffusion-reaction equation. In derived multi-physics elements, this function should be overloaded to reflect the chosen storage scheme. Note that these equations require that the unknown is always stored at the same index at each node.
RefineableQAdvectionDiffusionReactionElement(const RefineableQAdvectionDiffusionReactionElement< NREAGENT, DIM, NNODE_1D > &dummy)
Broken copy constructor.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes...
AdvectionDiffusionReactionReactionDerivFctPt & reaction_deriv_fct_pt()
Access function: Pointer to reaction derivatives function.
AdvectionDiffusionReactionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
A class for all elements that solve the Advection Diffusion Reaction equations using isoparametric el...
AdvectionDiffusionReactionWindFctPt & wind_fct_pt()
Access function: Pointer to wind function.
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
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...
AdvectionDiffusionReactionWindFctPt Wind_fct_pt
Pointer to wind function:
virtual unsigned nvertex_node() const
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement. ...
Vector< double > *& tau_pt()
Pointer to vector of dimensionless timescales.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: .
unsigned nvertex_node() const
Number of vertex nodes in the element.
Vector< double > * Diff_pt
Pointer to global diffusion coefficients.
AdvectionDiffusionReactionSourceFctPt Source_fct_pt
Pointer to source function:
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...
void broken_assign(const std::string &class_name)
Issue error message and terminate execution.
Vector< double > * Tau_pt
Pointer to global timescales.
RefineableQAdvectionDiffusionReactionElement()
Empty Constructor:
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Standard flux.from AdvectionDiffusionReaction equations.
Vector< double > *& diff_pt()
Pointer to vector of diffusion coefficients.
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...
RefineableAdvectionDiffusionReactionEquations(const RefineableAdvectionDiffusionReactionEquations< NREAGENT, DIM > &dummy)
Broken copy constructor.
AdvectionDiffusionReactionReactionFctPt Reaction_fct_pt
Pointer to reaction function.
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...
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: NREAGENT.
unsigned nnode() const
Return the number of 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...
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.