33 #ifndef OOMPH_REFINEABLE_POLAR_NAVIER_STOKES_HEADER
34 #define OOMPH_REFINEABLE_POLAR_NAVIER_STOKES_HEADER
38 #include <oomph-lib-config.h>
43 #include "../generic/refineable_quad_element.h"
44 #include "../generic/error_estimator.h"
99 unsigned n_element = element_pt.size();
100 for (
unsigned e = 0;
e < n_element;
e++)
113 unsigned n_element = element_pt.size();
114 for (
unsigned e = 0;
e < n_element;
e++)
134 unsigned num_entries = 3;
135 if (
flux.size() != num_entries)
137 std::ostringstream error_message;
138 error_message <<
"The flux vector has the wrong number of entries, "
139 <<
flux.size() <<
", whereas it should be " << num_entries
156 for (
unsigned i = 0;
i < 2;
i++)
158 flux[icount] = strainrate(
i,
i);
163 for (
unsigned i = 0;
i < 2;
i++)
165 for (
unsigned j =
i + 1;
j < 2;
j++)
167 flux[icount] = strainrate(
i,
j);
186 this->
Re_pt = cast_father_element_pt->
re_pt();
192 this->
G_pt = cast_father_element_pt->
g_pt();
237 unsigned n_node = this->
nnode();
239 for (
unsigned n = 0;
n < n_node;
n++)
251 unsigned n_node = this->
nnode();
253 for (
unsigned n = 0;
n < n_node;
n++)
260 for (
unsigned l = 0; l < n_pres; l++)
265 nod_pt->
unpin(p_index);
325 values.resize(3, 0.0);
328 for (
unsigned i = 0;
i < 2;
i++)
351 unsigned N_prev_time =
356 std::ostringstream error_message;
358 <<
"The value of t in get_interpolated_values(...), " <<
t
360 <<
"is greater than the number of previous stored timesteps";
369 values.resize(2 + 1);
372 for (
unsigned i = 0;
i < 2 + 1;
i++)
378 unsigned n_node = this->
nnode();
382 this->
shape(s, psif);
385 for (
unsigned i = 0;
i < 2;
i++)
389 for (
unsigned l = 0; l < n_node; l++)
391 values[
i] += this->
nodal_value(t, l, u_nodal_index) * psif[l];
456 unsigned total_index = 0;
458 unsigned NNODE_1D = 2;
462 for (
unsigned i = 0;
i < 2;
i++)
471 else if (
s[
i] == 1.0)
473 index[
i] = NNODE_1D - 1;
479 double float_index = 0.5 * (1.0 +
s[
i]) * (NNODE_1D - 1);
480 index[
i] =
int(float_index);
483 double excess = float_index - index[
i];
494 index[
i] *
static_cast<unsigned>(
pow(
static_cast<float>(NNODE_1D),
495 static_cast<int>(
i)));
528 return static_cast<unsigned>(
pow(2.0,
static_cast<int>(2)));
532 return this->
nnode();
540 const int& value_id)
const
548 return this->
shape(s, psi);
563 std::set<std::pair<Data*, unsigned>>& paired_load_data)
567 for (
unsigned i = 0;
i < 2;
i++)
573 unsigned n_node = this->
nnode();
574 for (
unsigned n = 0;
n < n_node;
n++)
586 for (
unsigned j = 0;
j < nmaster;
j++)
592 for (
unsigned i = 0;
i < 2;
i++)
594 paired_load_data.insert(
595 std::make_pair(master_nod_pt, u_index[
i]));
604 for (
unsigned i = 0;
i < 2;
i++)
606 paired_load_data.insert(
607 std::make_pair(this->
node_pt(
n), u_index[
i]));
617 for (
unsigned l = 0; l < n_pres; l++)
629 unsigned nmaster = hang_info_pt->
nmaster();
632 for (
unsigned m = 0;
m < nmaster;
m++)
636 paired_load_data.insert(
645 paired_load_data.insert(std::make_pair(pres_node_pt, p_index));
684 for (
unsigned l = 0; l < n_pres; l++)
737 values.resize(2, 0.0);
740 for (
unsigned i = 0;
i < 2;
i++)
764 unsigned N_prev_time =
769 std::ostringstream error_message;
771 <<
"The value of t in get_interpolated_values(...), " <<
t
773 <<
"is greater than the number of previous stored timesteps";
785 for (
unsigned i = 0;
i < 2;
i++)
791 unsigned n_node = this->
nnode();
795 this->
shape(s, psif);
798 for (
unsigned i = 0;
i < 2;
i++)
802 for (
unsigned l = 0; l < n_node; l++)
804 values[
i] += this->
nodal_value(t, l, u_nodal_index) * psif[l];
830 std::set<std::pair<Data*, unsigned>>& paired_load_data)
834 for (
unsigned i = 0;
i < 2;
i++)
840 unsigned n_node = this->
nnode();
841 for (
unsigned n = 0;
n < n_node;
n++)
853 for (
unsigned j = 0;
j < nmaster;
j++)
859 for (
unsigned i = 0;
i < 2;
i++)
861 paired_load_data.insert(
862 std::make_pair(master_nod_pt, u_index[
i]));
871 for (
unsigned i = 0;
i < 2;
i++)
873 paired_load_data.insert(
874 std::make_pair(this->
node_pt(
n), u_index[
i]));
882 for (
unsigned l = 0; l < n_pres; l++)
886 paired_load_data.insert(std::make_pair(
898 :
public virtual FaceGeometry<PolarCrouzeixRaviartElement>
911 using namespace QuadTreeNames;
920 double av_press = 0.0;
923 for (
unsigned ison = 0; ison < 4; ison++)
1009 ->
set_value(2, 0.5 * (slope1 + slope2));
1023 using namespace QuadTreeNames;
1042 else if (son_type ==
SE)
1048 else if (son_type ==
NE)
1055 else if (son_type ==
NW)
1071 for (
unsigned i = 1;
i < 3;
i++)
1073 double half_father_slope =
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
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
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_polar_navier_stokes_elements.h:901
FaceGeometry()
Definition: refineable_polar_navier_stokes_elements.h:661
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
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: oomph_definitions.h:222
Definition: polar_navier_stokes_elements.h:772
unsigned npres_pnst() const
Return number of pressure values.
Definition: polar_navier_stokes_elements.h:852
unsigned P_pnst_internal_index
Definition: polar_navier_stokes_elements.h:780
Definition: polar_navier_stokes_elements.h:108
double *& density_ratio_pt()
Pointer to Density ratio.
Definition: polar_navier_stokes_elements.h:322
double *& re_invfr_pt()
Pointer to global inverse Froude number.
Definition: polar_navier_stokes_elements.h:334
NavierStokesBodyForceFctPt & body_force_fct_pt()
Access function for the body-force pointer.
Definition: polar_navier_stokes_elements.h:352
double * ReSt_pt
Pointer to global Reynolds number x Strouhal number (=Womersley)
Definition: polar_navier_stokes_elements.h:157
double * Viscosity_Ratio_pt
Definition: polar_navier_stokes_elements.h:142
Vector< double > * G_pt
Pointer to global gravity Vector.
Definition: polar_navier_stokes_elements.h:164
NavierStokesBodyForceFctPt Body_force_fct_pt
Pointer to body force function.
Definition: polar_navier_stokes_elements.h:167
double *& alpha_pt()
Pointer to Alpha.
Definition: polar_navier_stokes_elements.h:290
NavierStokesSourceFctPt & source_fct_pt()
Access function for the source-function pointer.
Definition: polar_navier_stokes_elements.h:364
double * Re_pt
Pointer to global Reynolds number.
Definition: polar_navier_stokes_elements.h:154
double *& re_st_pt()
Pointer to product of Reynolds and Strouhal number (=Womersley number)
Definition: polar_navier_stokes_elements.h:296
void strain_rate_by_r(const Vector< double > &s, DenseMatrix< double > &strain_rate) const
Function to return polar strain multiplied by r.
Definition: polar_navier_stokes_elements.cc:909
double *& re_pt()
Pointer to Reynolds number.
Definition: polar_navier_stokes_elements.h:284
virtual unsigned u_index_pnst(const unsigned &i) const
Definition: polar_navier_stokes_elements.h:394
double * Density_Ratio_pt
Definition: polar_navier_stokes_elements.h:146
NavierStokesSourceFctPt Source_fct_pt
Pointer to volumetric source function.
Definition: polar_navier_stokes_elements.h:170
double * ReInvFr_pt
Definition: polar_navier_stokes_elements.h:161
double *& viscosity_ratio_pt()
Pointer to Viscosity Ratio.
Definition: polar_navier_stokes_elements.h:309
void interpolated_u_pnst(const Vector< double > &s, Vector< double > &veloc) const
Compute vector of FE interpolated velocity u at local coordinate s.
Definition: polar_navier_stokes_elements.h:662
double * Alpha_pt
Pointer to the angle alpha.
Definition: polar_navier_stokes_elements.h:151
double interpolated_p_pnst(const Vector< double > &s) const
Return FE interpolated pressure at local coordinate s.
Definition: polar_navier_stokes_elements.h:706
Vector< double > *& g_pt()
Pointer to Vector of gravitational components.
Definition: polar_navier_stokes_elements.h:346
Definition: polar_navier_stokes_elements.h:1014
void pshape_pnst(const Vector< double > &s, Shape &psi) const
Pressure shape functions at local coordinate s.
Definition: polar_navier_stokes_elements.h:1210
static const unsigned Pconv[]
Definition: polar_navier_stokes_elements.h:1022
unsigned npres_pnst() const
Return number of pressure values.
Definition: polar_navier_stokes_elements.h:1103
virtual int p_nodal_index_pnst()
Which nodal value represents the pressure?
Definition: polar_navier_stokes_elements.h:1068
Definition: refineable_elements.h:97
Tree * tree_pt()
Access function: Pointer to quadtree representation of this element.
Definition: refineable_elements.h:211
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
Definition: refineable_elements.h:539
Refineable version of Crouzeix Raviart elements. Generic class definitions.
Definition: refineable_polar_navier_stokes_elements.h:677
RefineablePolarCrouzeixRaviartElement()
Constructor.
Definition: refineable_polar_navier_stokes_elements.h:692
void further_build()
Definition: refineable_polar_navier_stokes_elements.h:1018
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_polar_navier_stokes_elements.h:724
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_polar_navier_stokes_elements.h:718
void further_setup_hanging_nodes()
Definition: refineable_polar_navier_stokes_elements.h:811
void unpin_elemental_pressure_dofs()
Unpin all internal pressure dofs.
Definition: refineable_polar_navier_stokes_elements.h:680
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Definition: refineable_polar_navier_stokes_elements.h:753
unsigned nrecovery_order()
Definition: refineable_polar_navier_stokes_elements.h:712
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: DIM (velocities)
Definition: refineable_polar_navier_stokes_elements.h:701
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: refineable_polar_navier_stokes_elements.h:733
void rebuild_from_sons(Mesh *&mesh_pt)
2D Rebuild from sons: Reconstruct pressure from the (merged) sons
Definition: refineable_polar_navier_stokes_elements.h:908
void insert_load_data(std::set< std::pair< Data *, unsigned >> &paired_load_data)
Definition: refineable_polar_navier_stokes_elements.h:829
Definition: refineable_polar_navier_stokes_elements.h:58
virtual void pin_elemental_redundant_nodal_pressure_dofs()
Definition: refineable_polar_navier_stokes_elements.h:72
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
Definition: refineable_polar_navier_stokes_elements.h:123
virtual void fill_in_generic_residual_contribution(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Start of what would've been refineable_navier_stokes_elements.cc //.
Definition: refineable_polar_navier_stokes_elements.cc:48
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Definition: refineable_polar_navier_stokes_elements.h:131
virtual Node * pressure_node_pt(const unsigned &n_p)
Definition: refineable_polar_navier_stokes_elements.h:62
static void pin_redundant_nodal_pressures(const Vector< GeneralisedElement * > &element_pt)
Definition: refineable_polar_navier_stokes_elements.h:94
void further_build()
Further build, pass the pointers down to the sons.
Definition: refineable_polar_navier_stokes_elements.h:174
static void unpin_all_pressure_dofs(const Vector< GeneralisedElement * > &element_pt)
Unpin all pressure dofs in elements listed in vector.
Definition: refineable_polar_navier_stokes_elements.h:108
RefineablePolarNavierStokesEquations()
Constructor.
Definition: refineable_polar_navier_stokes_elements.h:76
virtual void unpin_elemental_pressure_dofs()=0
Unpin all pressure dofs in the element.
Definition: refineable_polar_navier_stokes_elements.h:224
void unpin_elemental_pressure_dofs()
Unpin all pressure dofs.
Definition: refineable_polar_navier_stokes_elements.h:233
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
Definition: refineable_polar_navier_stokes_elements.h:296
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_polar_navier_stokes_elements.h:306
unsigned ninterpolating_node_1d(const int &value_id)
Definition: refineable_polar_navier_stokes_elements.h:510
unsigned ncont_interpolated_values() const
Definition: refineable_polar_navier_stokes_elements.h:290
Node * get_interpolating_node_at_local_coordinate(const Vector< double > &s, const int &value_id)
Definition: refineable_polar_navier_stokes_elements.h:449
Node * pressure_node_pt(const unsigned &n_p)
Pointer to n_p-th pressure node.
Definition: refineable_polar_navier_stokes_elements.h:227
Node * interpolating_node_pt(const unsigned &n, const int &value_id)
Definition: refineable_polar_navier_stokes_elements.h:412
void interpolating_basis(const Vector< double > &s, Shape &psi, const int &value_id) const
Definition: refineable_polar_navier_stokes_elements.h:538
void pin_elemental_redundant_nodal_pressure_dofs()
Pin all nodal pressure dofs that are not required.
Definition: refineable_polar_navier_stokes_elements.h:246
unsigned ninterpolating_node(const int &value_id)
Definition: refineable_polar_navier_stokes_elements.h:524
unsigned required_nvalue(const unsigned &n) const
Definition: refineable_polar_navier_stokes_elements.h:283
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_polar_navier_stokes_elements.h:312
RefineablePolarTaylorHoodElement()
Constructor.
Definition: refineable_polar_navier_stokes_elements.h:272
void insert_load_data(std::set< std::pair< Data *, unsigned >> &paired_load_data)
Definition: refineable_polar_navier_stokes_elements.h:562
unsigned nrecovery_order()
Definition: refineable_polar_navier_stokes_elements.h:300
double local_one_d_fraction_of_interpolating_node(const unsigned &n1d, const unsigned &i, const int &value_id)
Definition: refineable_polar_navier_stokes_elements.h:429
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: refineable_polar_navier_stokes_elements.h:321
void further_setup_hanging_nodes()
Definition: refineable_polar_navier_stokes_elements.h:403
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Definition: refineable_polar_navier_stokes_elements.h:341
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
virtual unsigned nprev_values() const =0
Number of previous values available: 0 for static, 1 for BDF<1>,...
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
TreeRoot *& root_pt()
Return pointer to root of the tree.
Definition: tree.h:141
RealScalar s
Definition: level1_cplx_impl.h:130
int * m
Definition: level2_cplx_impl.h:294
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 pow(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:625
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
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
void set(Container &c, Position position, const Value &value)
Definition: stdlist_overload.cpp:36
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2