29 #ifndef OOMPH_REFINEABLE_ADVECTION_DIFFUSION_ELEMENTS_HEADER
30 #define OOMPH_REFINEABLE_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>
101 unsigned n_node =
nnode();
116 for (
unsigned l = 0; l < n_node; l++)
118 values[0] += this->
nodal_value(l, u_nodal_index) * psi[l];
134 const unsigned n_node =
nnode();
149 for (
unsigned l = 0; l < n_node; l++)
151 values[0] += this->
nodal_value(t, l, u_nodal_index) * psi[l];
186 this->
Pe_pt = cast_father_element_pt->
pe_pt();
204 unsigned n_node = this->
nnode();
219 unsigned n_u_dof = 0;
220 for (
unsigned l = 0; l < n_node; l++)
222 unsigned n_master = 1;
231 n_master = hang_info_pt->
nmaster();
240 for (
unsigned m = 0;
m < n_master;
m++)
264 du_ddata.resize(n_u_dof, 0.0);
265 global_eqn_number.resize(n_u_dof, 0);
270 for (
unsigned l = 0; l < n_node; l++)
272 unsigned n_master = 1;
273 double hang_weight = 1.0;
282 n_master = hang_info_pt->
nmaster();
291 for (
unsigned m = 0;
m < n_master;
m++)
321 global_eqn_number[count] = global_eqn;
323 du_ddata[count] = psi[l] * hang_weight;
350 template<
unsigned DIM,
unsigned NNODE_1D>
401 return (NNODE_1D - 1);
420 template<
unsigned DIM,
unsigned NNODE_1D>
422 :
public virtual QElement<DIM - 1, NNODE_1D>
Definition: advection_diffusion_elements.h:52
bool ALE_is_disabled
Definition: advection_diffusion_elements.h:588
double * Pe_pt
Pointer to global Peclet number.
Definition: advection_diffusion_elements.h:574
double *& pe_pt()
Pointer to Peclet number.
Definition: advection_diffusion_elements.h:324
AdvectionDiffusionWindFctPt & wind_fct_pt()
Access function: Pointer to wind function.
Definition: advection_diffusion_elements.h:305
AdvectionDiffusionSourceFctPt Source_fct_pt
Pointer to source function:
Definition: advection_diffusion_elements.h:580
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: .
Definition: advection_diffusion_elements.h:387
AdvectionDiffusionWindFctPt Wind_fct_pt
Pointer to wind function:
Definition: advection_diffusion_elements.h:583
double * PeSt_pt
Pointer to global Peclet number multiplied by Strouhal number.
Definition: advection_diffusion_elements.h:577
virtual unsigned u_index_adv_diff() const
Definition: advection_diffusion_elements.h:91
double *& pe_st_pt()
Pointer to Peclet number multipled by Strouha number.
Definition: advection_diffusion_elements.h:336
AdvectionDiffusionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Definition: advection_diffusion_elements.h:291
long & eqn_number(const unsigned &i)
Return the equation number of the i-th stored variable.
Definition: nodes.h:367
Definition: error_estimator.h:79
FaceGeometry()
Definition: refineable_advection_diffusion_elements.h:427
Definition: elements.h:4998
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
virtual unsigned nvertex_node() const
Definition: elements.h:2491
virtual void shape(const Vector< double > &s, Shape &psi) const =0
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
virtual Node * vertex_node_pt(const unsigned &j) const
Definition: elements.h:2500
double const & master_weight(const unsigned &i) const
Return weight for dofs on i-th master node.
Definition: nodes.h:808
Node *const & master_node_pt(const unsigned &i) const
Return a pointer to the i-th master node.
Definition: nodes.h:791
unsigned nmaster() const
Return the number of master nodes.
Definition: nodes.h:785
HangInfo *const & hanging_pt() const
Definition: nodes.h:1228
bool is_hanging() const
Test whether the node is geometrically hanging.
Definition: nodes.h:1285
Definition: advection_diffusion_elements.h:610
Definition: Qelements.h:459
Definition: refineable_advection_diffusion_elements.h:58
void further_build()
Further build: Copy source function pointer from father element.
Definition: refineable_advection_diffusion_elements.h:177
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Definition: refineable_advection_diffusion_elements.h:126
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
Definition: refineable_advection_diffusion_elements.h:77
RefineableAdvectionDiffusionEquations(const RefineableAdvectionDiffusionEquations< DIM > &dummy)=delete
Broken copy constructor.
void operator=(const RefineableAdvectionDiffusionEquations< DIM > &)=delete
Broken assignment operator.
void fill_in_generic_residual_contribution_adv_diff(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: refineable_advection_diffusion_elements.cc:38
void dinterpolated_u_adv_diff_ddata(const Vector< double > &s, Vector< double > &du_ddata, Vector< unsigned > &global_eqn_number)
Definition: refineable_advection_diffusion_elements.h:199
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Definition: refineable_advection_diffusion_elements.h:84
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: refineable_advection_diffusion_elements.h:94
RefineableAdvectionDiffusionEquations()
Empty Constructor.
Definition: refineable_advection_diffusion_elements.h:61
Definition: refineable_elements.h:97
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
Definition: refineable_elements.h:539
Definition: refineable_advection_diffusion_elements.h:355
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_advection_diffusion_elements.h:383
RefineableQAdvectionDiffusionElement(const RefineableQAdvectionDiffusionElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_advection_diffusion_elements.h:389
void further_setup_hanging_nodes()
Definition: refineable_advection_diffusion_elements.h:406
unsigned nrecovery_order()
Definition: refineable_advection_diffusion_elements.h:399
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 1.
Definition: refineable_advection_diffusion_elements.h:377
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
Definition: refineable_advection_diffusion_elements.h:395
RefineableQAdvectionDiffusionElement()
Empty Constructor:
Definition: refineable_advection_diffusion_elements.h:358
void operator=(const RefineableQAdvectionDiffusionElement< DIM, NNODE_1D > &)=delete
Broken assignment operator.
Definition: Qelements.h:2259
RealScalar s
Definition: level1_cplx_impl.h:130
int * m
Definition: level2_cplx_impl.h:294
#define DIM
Definition: linearised_navier_stokes_elements.h:44
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
t
Definition: plotPSD.py:36
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2