29 #ifndef OOMPH_REFINEABLE_AXISYM_ADVECTION_DIFFUSION_ELEMENTS_HEADER
30 #define OOMPH_REFINEABLE_AXISYM_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"
107 unsigned n_node =
nnode();
122 for (
unsigned l = 0; l < n_node; l++)
124 values[0] += this->
nodal_value(l, u_nodal_index) * psi[l];
140 const unsigned n_node =
nnode();
155 for (
unsigned l = 0; l < n_node; l++)
157 values[0] += this->
nodal_value(t, l, u_nodal_index) * psi[l];
164 "Time-dependent version of get_interpolated_values() ";
165 error_message +=
"not implemented for this element \n";
167 "RefineableAxisymAdvectionDiffusionEquations::get_"
168 "interpolated_values()",
169 OOMPH_EXCEPTION_LOCATION);
196 this->
Pe_pt = cast_father_element_pt->
pe_pt();
214 unsigned n_node = this->
nnode();
229 unsigned n_u_dof = 0;
230 for (
unsigned l = 0; l < n_node; l++)
232 unsigned n_master = 1;
241 n_master = hang_info_pt->
nmaster();
250 for (
unsigned m = 0; m < n_master; m++)
274 du_ddata.resize(n_u_dof, 0.0);
275 global_eqn_number.resize(n_u_dof, 0);
280 for (
unsigned l = 0; l < n_node; l++)
282 unsigned n_master = 1;
283 double hang_weight = 1.0;
292 n_master = hang_info_pt->
nmaster();
301 for (
unsigned m = 0; m < n_master; m++)
331 global_eqn_number[count] = global_eqn;
333 du_ddata[count] = psi[l] * hang_weight;
360 template<
unsigned NNODE_1D>
412 return (NNODE_1D - 1);
432 template<
unsigned NNODE_1D>
434 :
public virtual QElement<1, NNODE_1D>
A class for all elements that solve the Advection Diffusion equations in a cylindrical polar coordina...
double * PeSt_pt
Pointer to global Peclet number multiplied by Strouhal number.
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: [du/dr,du/dz].
AxisymAdvectionDiffusionSourceFctPt Source_fct_pt
Pointer to source function:
virtual unsigned u_index_axi_adv_diff() const
Broken assignment operator.
AxisymAdvectionDiffusionWindFctPt Wind_fct_pt
Pointer to wind function:
double *& pe_st_pt()
Pointer to Peclet number multipled by Strouha number.
AxisymAdvectionDiffusionWindFctPt & wind_fct_pt()
Access function: Pointer to wind function.
double *& pe_pt()
Pointer to Peclet number.
AxisymAdvectionDiffusionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
double * Pe_pt
Pointer to global Peclet number.
long & eqn_number(const unsigned &i)
Return the equation number of the i-th stored variable.
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.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
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.
Class that contains data for hanging nodes.
double const & master_weight(const unsigned &i) const
Return weight for dofs on i-th master node.
Node *const & master_node_pt(const unsigned &i) const
Return a pointer to the i-th master node.
unsigned nmaster() const
Return the number of master nodes.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
HangInfo *const & hanging_pt() const
Return pointer to hanging node data (this refers to the geometric hanging node status) (const version...
bool is_hanging() const
Test whether the node is geometrically hanging.
An OomphLibError object which should be thrown when an run-time error is encountered....
//////////////////////////////////////////////////////////////////////// ////////////////////////////...
unsigned nvertex_node() const
Number of vertex nodes in the element.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
A version of the Advection Diffusion in axisym coordinates equations that can be used with non-unifor...
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Standard flux.from AdvectionDiffusion equations.
RefineableAxisymAdvectionDiffusionEquations(const RefineableAxisymAdvectionDiffusionEquations &dummy)=delete
Broken copy constructor.
double geometric_jacobian(const Vector< double > &x)
Fill in the geometric Jacobian, which in this case is r.
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...
void fill_in_generic_residual_contribution_axi_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...
unsigned num_Z2_flux_terms()
Broken assignment operator.
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...
void dinterpolated_u_adv_diff_ddata(const Vector< double > &s, Vector< double > &du_ddata, Vector< unsigned > &global_eqn_number)
Compute the derivatives of the i-th component of velocity at point s with respect to all data that ca...
void further_build()
Further build: Copy source function pointer from father element.
RefineableAxisymAdvectionDiffusionEquations()
Empty Constructor.
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 QAxisymAdvectionDiffusionElement. Inherit from the standard QAxisymAdvectionDif...
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.
unsigned nvertex_node() const
Number of vertex nodes in the element.
unsigned ncont_interpolated_values() const
Broken assignment operator.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes....
RefineableQAxisymAdvectionDiffusionElement(const RefineableQAxisymAdvectionDiffusionElement< NNODE_1D > &dummy)=delete
Broken copy constructor.
RefineableQAxisymAdvectionDiffusionElement()
Empty Constructor:
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
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...
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...