28 #ifndef OOMPH_REFINEABLE_LINEARISED_NAVIER_STOKES_ELEMENTS_HEADER
29 #define OOMPH_REFINEABLE_LINEARISED_NAVIER_STOKES_ELEMENTS_HEADER
33 #include <oomph-lib-config.h>
37 #include "../generic/refineable_quad_element.h"
38 #include "../generic/error_estimator.h"
88 unsigned num_entries = 2 * (
DIM + ((
DIM *
DIM) -
DIM) / 2);
89 if (
flux.size() != num_entries)
91 std::ostringstream error_message;
92 error_message <<
"The flux vector has the wrong number of entries, "
93 <<
flux.size() <<
", whereas it should be " << num_entries
97 "RefineableLinearisedNavierStokesEquations::get_Z2_flux()",
113 for (
unsigned i = 0;
i <
DIM;
i++)
115 flux[icount] = real_strainrate(
i,
i);
120 for (
unsigned i = 0;
i <
DIM;
i++)
122 for (
unsigned j =
i + 1;
j <
DIM;
j++)
124 flux[icount] = real_strainrate(
i,
j);
130 for (
unsigned i = 0;
i <
DIM;
i++)
132 flux[icount] = imag_strainrate(
i,
i);
137 for (
unsigned i = 0;
i <
DIM;
i++)
139 for (
unsigned j =
i + 1;
j <
DIM;
j++)
141 flux[icount] = imag_strainrate(
i,
j);
162 this->
Re_pt = cast_father_element_pt->
re_pt();
193 const unsigned n_element = element_pt.size();
194 for (
unsigned e = 0;
e < n_element;
e++)
207 const unsigned n_element = element_pt.size();
208 for (
unsigned e = 0;
e < n_element;
e++)
250 for (
unsigned l = 0; l < n_pres; l++)
253 for (
unsigned i = 0;
i < 2;
i++)
279 using namespace QuadTreeNames;
291 for (
unsigned ison = 0; ison < 4; ison++)
294 for (
unsigned i = 0;
i < 2;
i++)
307 for (
unsigned i = 0;
i < 2;
i++)
316 for (
unsigned i = 0;
i < 2;
i++)
415 const unsigned n_values = 4 *
DIM;
418 values.resize(n_values, 0.0);
421 for (
unsigned i = 0;
i < n_values;
i++)
441 values.resize(4 *
DIM);
444 for (
unsigned i = 0;
i < 4 *
DIM;
i++)
450 const unsigned n_node =
nnode();
457 for (
unsigned i = 0;
i < (4 *
DIM);
i++)
461 for (
unsigned l = 0; l < n_node; l++)
481 using namespace QuadTreeNames;
500 else if (son_type ==
SE)
506 else if (son_type ==
NE)
513 else if (son_type ==
NW)
531 for (
unsigned i = 0;
i < 2;
i++)
572 :
public virtual FaceGeometry<LinearisedQCrouzeixRaviartElement>
586 FaceGeometry<LinearisedQCrouzeixRaviartElement>>
614 return this->
node_pt(this->Pconv[n_p]);
621 const unsigned n_node = this->
nnode();
625 for (
unsigned i = 0;
i < 2;
i++)
631 for (
unsigned n = 0;
n < n_node;
n++)
633 for (
unsigned i = 0;
i < 2;
i++)
644 const unsigned n_node = this->
nnode();
648 for (
unsigned i = 0;
i < 2;
i++)
654 for (
unsigned n = 0;
n < n_node;
n++)
656 for (
unsigned i = 0;
i < 2;
i++)
664 for (
unsigned l = 0; l < n_pres; l++)
667 for (
unsigned i = 0;
i < 2;
i++)
671 nod_pt->
unpin(p_index[
i]);
733 const unsigned n_values = 2 * (
DIM + 1);
736 values.resize(n_values, 0.0);
739 for (
unsigned i = 0;
i < (2 *
DIM);
i++)
745 for (
unsigned i = 0;
i < 2;
i++)
760 values.resize(2 * (
DIM + 1));
763 for (
unsigned i = 0;
i < 2 * (
DIM + 1);
i++)
769 const unsigned n_node =
nnode();
776 for (
unsigned i = 0;
i < (2 *
DIM);
i++)
780 for (
unsigned l = 0; l < n_node; l++)
788 for (
unsigned i = 0;
i < 2;
i++)
799 for (
unsigned i = 0;
i < 2;
i++)
813 if (n_value == (2 *
DIM) || n_value == ((2 *
DIM) + 1))
815 return this->
node_pt(this->Pconv[
n]);
831 if (n_value == (2 *
DIM) || n_value == ((2 *
DIM) + 1))
851 if (n_value ==
static_cast<int>(2 *
DIM) ||
852 n_value ==
static_cast<int>((2 *
DIM) + 1))
855 unsigned total_index = 0;
857 const unsigned NNODE_1D = 2;
862 for (
unsigned i = 0;
i < 2;
i++)
872 else if (
s[
i] == 1.0)
874 index[
i] = NNODE_1D - 1;
881 double float_index = 0.5 * (1.0 +
s[
i]) * (NNODE_1D - 1);
882 index[
i] =
int(float_index);
885 double excess = float_index - index[
i];
896 index[
i] *
static_cast<unsigned>(
pow(
static_cast<float>(NNODE_1D),
897 static_cast<int>(
i)));
900 return this->
node_pt(this->Pconv[total_index]);
914 if (n_value == (2 *
DIM) || n_value == ((2 *
DIM) + 1))
928 if (n_value == (2 *
DIM) || n_value == ((2 *
DIM) + 1))
934 return this->
nnode();
942 const int& n_value)
const
944 if (n_value == (2 *
DIM) || n_value == ((2 *
DIM) + 1))
950 return this->
shape(s, psi);
968 :
public virtual FaceGeometry<LinearisedQTaylorHoodElement>
981 :
public virtual FaceGeometry<FaceGeometry<LinearisedQTaylorHoodElement>>
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Array< double, 1, 3 > e(1./3., 0.5, 2.)
void pin(const unsigned &i)
Pin the i-th stored variable.
Definition: nodes.h:385
void unpin(const unsigned &i)
Unpin the i-th stored variable.
Definition: nodes.h:391
void set_value(const unsigned &i, const double &value_)
Definition: nodes.h:271
double value(const unsigned &i) const
Definition: nodes.h:293
Definition: error_estimator.h:79
FaceGeometry()
Definition: refineable_linearised_navier_stokes_elements.h:589
FaceGeometry()
Definition: refineable_linearised_navier_stokes_elements.h:984
FaceGeometry()
Definition: refineable_linearised_navier_stokes_elements.h:575
FaceGeometry()
Definition: refineable_linearised_navier_stokes_elements.h:971
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
static const double Node_location_tolerance
Definition: elements.h:1374
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 * get_node_at_local_coordinate(const Vector< double > &s) const
Definition: elements.cc:3882
virtual Node * vertex_node_pt(const unsigned &j) const
Definition: elements.h:2500
virtual unsigned nnode_1d() const
Definition: elements.h:2218
virtual double local_one_d_fraction_of_node(const unsigned &n1d, const unsigned &i)
Definition: elements.h:1858
Data *& internal_data_pt(const unsigned &i)
Return a pointer to i-th internal data object.
Definition: elements.h:622
Definition: linearised_navier_stokes_elements.h:54
virtual unsigned u_index_linearised_nst(const unsigned &i) const
Definition: linearised_navier_stokes_elements.h:380
double *& re_pt()
Pointer to Reynolds number.
Definition: linearised_navier_stokes_elements.h:275
LinearisedNavierStokesEigenfunctionNormalisationElement * normalisation_element_pt()
Pointer to normalisation element.
Definition: linearised_navier_stokes_elements.h:299
void strain_rate(const Vector< double > &s, DenseMatrix< double > &strain_rate, const unsigned &real) const
Definition: linearised_navier_stokes_elements.cc:248
bool ALE_is_disabled
Definition: linearised_navier_stokes_elements.h:115
void set_eigenfunction_normalisation_element(LinearisedNavierStokesEigenfunctionNormalisationElement *const &normalisation_el_pt)
the boolean flag check_nodal_data is set to false.
Definition: linearised_navier_stokes_elements.h:305
double * Viscosity_Ratio_pt
Definition: linearised_navier_stokes_elements.h:73
double interpolated_p_linearised_nst(const Vector< double > &s, const unsigned &i) const
Definition: linearised_navier_stokes_elements.h:554
double * ReSt_pt
Pointer to global Reynolds number x Strouhal number (=Womersley)
Definition: linearised_navier_stokes_elements.h:83
double * Re_pt
Pointer to global Reynolds number.
Definition: linearised_navier_stokes_elements.h:80
double *& viscosity_ratio_pt()
Pointer to the viscosity ratio.
Definition: linearised_navier_stokes_elements.h:335
double *& re_st_pt()
Pointer to product of Reynolds and Strouhal number (=Womersley number)
Definition: linearised_navier_stokes_elements.h:281
double *& density_ratio_pt()
Pointer to the density ratio.
Definition: linearised_navier_stokes_elements.h:348
double * Density_Ratio_pt
Definition: linearised_navier_stokes_elements.h:77
double interpolated_u_linearised_nst(const Vector< double > &s, const unsigned &i) const
Definition: linearised_navier_stokes_elements.h:525
virtual int p_index_linearised_nst(const unsigned &i) const
Which nodal value represents the pressure?
Definition: linearised_navier_stokes_elements.h:448
Definition: linearised_navier_stokes_elements.h:597
unsigned npres_linearised_nst() const
Definition: linearised_navier_stokes_elements.h:757
Vector< unsigned > P_linearised_nst_internal_index
Definition: linearised_navier_stokes_elements.h:607
void pshape_linearised_nst(const Vector< double > &s, Shape &psi) const
Compute the pressure shape functions at local coordinate s.
Definition: linearised_navier_stokes_elements.h:866
Definition: linearised_navier_stokes_elements.h:927
bool is_hanging() const
Test whether the node is geometrically hanging.
Definition: nodes.h:1285
Definition: oomph_definitions.h:222
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_linearised_navier_stokes_elements.h:51
void further_build()
Further build: pass the pointers down to the sons.
Definition: refineable_linearised_navier_stokes_elements.h:148
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
Definition: refineable_linearised_navier_stokes_elements.h:77
static void pin_redundant_nodal_pressures(const Vector< GeneralisedElement * > &element_pt)
Definition: refineable_linearised_navier_stokes_elements.h:188
virtual void pin_elemental_redundant_nodal_pressure_dofs()
Definition: refineable_linearised_navier_stokes_elements.h:65
virtual Node * pressure_node_pt(const unsigned &n_p)
Definition: refineable_linearised_navier_stokes_elements.h:55
virtual void unpin_elemental_pressure_dofs()=0
Unpin all pressure dofs in the element.
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Definition: refineable_linearised_navier_stokes_elements.h:85
static void unpin_all_pressure_dofs(const Vector< GeneralisedElement * > &element_pt)
Unpin all pressure dofs in elements listed in Vector.
Definition: refineable_linearised_navier_stokes_elements.h:202
void fill_in_generic_residual_contribution_linearised_nst(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: refineable_linearised_navier_stokes_elements.cc:40
RefineableLinearisedNavierStokesEquations()
Empty Constructor.
Definition: refineable_linearised_navier_stokes_elements.h:69
Definition: refineable_linearised_navier_stokes_elements.h:243
void further_build()
Definition: refineable_linearised_navier_stokes_elements.h:476
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Definition: refineable_linearised_navier_stokes_elements.h:436
unsigned nrecovery_order()
Definition: refineable_linearised_navier_stokes_elements.h:390
RefineableLinearisedQCrouzeixRaviartElement()
Constructor.
Definition: refineable_linearised_navier_stokes_elements.h:262
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_linearised_navier_stokes_elements.h:396
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: Reconstruct pressure from the (merged) sons.
Definition: refineable_linearised_navier_stokes_elements.h:277
void unpin_elemental_pressure_dofs()
Unpin all the internal pressure freedoms.
Definition: refineable_linearised_navier_stokes_elements.h:246
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: refineable_linearised_navier_stokes_elements.h:411
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 4*DIM (velocities)
Definition: refineable_linearised_navier_stokes_elements.h:271
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_linearised_navier_stokes_elements.h:402
void further_setup_hanging_nodes()
Definition: refineable_linearised_navier_stokes_elements.h:470
Definition: refineable_linearised_navier_stokes_elements.h:609
unsigned nrecovery_order()
Definition: refineable_linearised_navier_stokes_elements.h:707
double local_one_d_fraction_of_interpolating_node(const unsigned &n1d, const unsigned &i, const int &n_value)
Definition: refineable_linearised_navier_stokes_elements.h:827
unsigned required_nvalue(const unsigned &n) const
Definition: refineable_linearised_navier_stokes_elements.h:690
unsigned ncont_interpolated_values() const
Definition: refineable_linearised_navier_stokes_elements.h:697
RefineableLinearisedQTaylorHoodElement()
Constructor:
Definition: refineable_linearised_navier_stokes_elements.h:679
Node * get_interpolating_node_at_local_coordinate(const Vector< double > &s, const int &n_value)
Definition: refineable_linearised_navier_stokes_elements.h:847
void unpin_elemental_pressure_dofs()
Unpin all pressure dofs.
Definition: refineable_linearised_navier_stokes_elements.h:618
Node * pressure_node_pt(const unsigned &n_p)
Pointer to n_p-th pressure node.
Definition: refineable_linearised_navier_stokes_elements.h:612
void further_setup_hanging_nodes()
Definition: refineable_linearised_navier_stokes_elements.h:797
Node * interpolating_node_pt(const unsigned &n, const int &n_value)
Definition: refineable_linearised_navier_stokes_elements.h:809
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
Definition: refineable_linearised_navier_stokes_elements.h:703
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: refineable_linearised_navier_stokes_elements.h:728
void pin_elemental_redundant_nodal_pressure_dofs()
Unpin the proper nodal pressure dofs.
Definition: refineable_linearised_navier_stokes_elements.h:641
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Definition: refineable_linearised_navier_stokes_elements.h:755
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_linearised_navier_stokes_elements.h:713
unsigned ninterpolating_node_1d(const int &n_value)
Definition: refineable_linearised_navier_stokes_elements.h:912
unsigned ninterpolating_node(const int &n_value)
Definition: refineable_linearised_navier_stokes_elements.h:926
void interpolating_basis(const Vector< double > &s, Shape &psi, const int &n_value) const
Definition: refineable_linearised_navier_stokes_elements.h:940
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_linearised_navier_stokes_elements.h:719
QuadTree * quadtree_pt()
Pointer to quadtree representation of this element.
Definition: refineable_quad_element.h:152
void setup_hang_for_value(const int &value_id)
Definition: refineable_quad_element.cc:1506
Definition: Qelements.h:2259
RefineableElement * object_pt() const
Definition: tree.h:88
Tree * father_pt() const
Return pointer to father: NULL if it's a root node.
Definition: tree.h:235
int son_type() const
Return son type.
Definition: tree.h:214
Tree * son_pt(const int &son_index) const
Definition: tree.h:103
RealScalar s
Definition: level1_cplx_impl.h:130
#define DIM
Definition: linearised_navier_stokes_elements.h:44
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 pow(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:625
squared absolute value
Definition: GlobalFunctions.h:87
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
@ SE
Definition: quadtree.h:57
@ NW
Definition: quadtree.h:58
@ NE
Definition: quadtree.h:59
@ SW
Definition: quadtree.h:56
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
t
Definition: plotPSD.py:36
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2