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"
80 return 2 * (DIM + ((DIM * DIM) - DIM) / 2);
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()",
98 OOMPH_EXCEPTION_LOCATION);
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>>
void pin(const unsigned &i)
Pin the i-th stored variable.
void unpin(const unsigned &i)
Unpin the i-th stored variable.
void set_value(const unsigned &i, const double &value_)
Set the i-th stored data value to specified value. The only reason that we require an explicit set fu...
double value(const unsigned &i) const
Return i-th stored value. This function is not virtual so that it can be inlined. This means that if ...
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
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 Node * get_node_at_local_coordinate(const Vector< double > &s) const
If there is a node at this local coordinate, return the pointer to the node.
virtual unsigned nvertex_node() const
Return the number of vertex nodes in this element. Broken virtual function in "pure" finite elements.
static const double Node_location_tolerance
Default value for the tolerance to be used when locating nodes via local coordinates.
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.
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.
virtual unsigned nnode_1d() const
Return the number of nodes along one edge of the element Default is to return zero — must be overload...
virtual double local_one_d_fraction_of_node(const unsigned &n1d, const unsigned &i)
Get the local fraction of any node in the n-th position in a one dimensional expansion along the i-th...
Data *& internal_data_pt(const unsigned &i)
Return a pointer to i-th internal data object.
A class for elements that solve the linearised version of the unsteady Navier–Stokes equations in cyl...
virtual unsigned u_index_linearised_nst(const unsigned &i) const
Return the index at which the i-th unknown velocity component is stored. The default value,...
double *& re_pt()
Pointer to Reynolds number.
LinearisedNavierStokesEigenfunctionNormalisationElement * normalisation_element_pt()
Pointer to normalisation element.
void strain_rate(const Vector< double > &s, DenseMatrix< double > &strain_rate, const unsigned &real) const
Strain-rate tensor: where (in that order)
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when the time-derivatives are computed....
void set_eigenfunction_normalisation_element(LinearisedNavierStokesEigenfunctionNormalisationElement *const &normalisation_el_pt)
the boolean flag check_nodal_data is set to false.
double * Viscosity_Ratio_pt
Pointer to the viscosity ratio (relative to the viscosity used in the definition of the Reynolds numb...
double interpolated_p_linearised_nst(const Vector< double > &s, const unsigned &i) const
Return the i-th component of the FE interpolated pressure p[i] at local coordinate s.
double * ReSt_pt
Pointer to global Reynolds number x Strouhal number (=Womersley)
double * Re_pt
Pointer to global Reynolds number.
double *& viscosity_ratio_pt()
Pointer to the viscosity ratio.
double *& re_st_pt()
Pointer to product of Reynolds and Strouhal number (=Womersley number)
double *& density_ratio_pt()
Pointer to the density ratio.
double * Density_Ratio_pt
Pointer to the density ratio (relative to the density used in the definition of the Reynolds number)
double interpolated_u_linearised_nst(const Vector< double > &s, const unsigned &i) const
Compute the element's residual Vector and the jacobian matrix. Virtual function can be overloaded by ...
virtual int p_index_linearised_nst(const unsigned &i) const
Which nodal value represents the pressure?
/////////////////////////////////////////////////////////////////////////// /////////////////////////...
unsigned npres_linearised_nst() const
Return number of pressure values corresponding to a single pressure component.
Vector< unsigned > P_linearised_nst_internal_index
Internal indices that indicate at which internal data the pressure values are stored....
void pshape_linearised_nst(const Vector< double > &s, Shape &psi) const
Compute the pressure shape functions at local coordinate s.
/////////////////////////////////////////////////////////////////////////// /////////////////////////...
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
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....
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 the linearised axisymmetric Navier–Stokes equations.
void further_build()
Further build: pass the pointers down to the sons.
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
static void pin_redundant_nodal_pressures(const Vector< GeneralisedElement * > &element_pt)
Loop over all elements in Vector (which typically contains all the elements in a fluid mesh) and pin ...
virtual void pin_elemental_redundant_nodal_pressure_dofs()
Pin unused nodal pressure dofs (empty by default, because by default pressure dofs are not associated...
virtual Node * pressure_node_pt(const unsigned &n_p)
Pointer to n_p-th pressure node (Default: NULL, indicating that pressure is not based on nodal interp...
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)
Get 'flux' for Z2 error recovery: Upper triangular entries in strain rate tensor.
static void unpin_all_pressure_dofs(const Vector< GeneralisedElement * > &element_pt)
Unpin all pressure dofs in elements listed in Vector.
void fill_in_generic_residual_contribution_linearised_nst(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Add element's contribution to the elemental residual vector and/or Jacobian matrix flag=1: compute bo...
RefineableLinearisedNavierStokesEquations()
Empty Constructor.
/////////////////////////////////////////////////////////////////////////// /////////////////////////...
void further_build()
Further build for Crouzeix_Raviart interpolates the internal pressure dofs from father element: Make ...
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Get all function values [U^C,U^S,...,P^S] at previous timestep t (t=0: present; t>0: previous timeste...
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
RefineableLinearisedQCrouzeixRaviartElement()
Constructor.
unsigned nvertex_node() const
Number of vertex nodes in the element.
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: Reconstruct pressure from the (merged) sons.
void unpin_elemental_pressure_dofs()
Unpin all the internal pressure freedoms.
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...
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 4*DIM (velocities)
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes....
/////////////////////////////////////////////////////////////////////////// /////////////////////////...
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
double local_one_d_fraction_of_interpolating_node(const unsigned &n1d, const unsigned &i, const int &n_value)
The pressure nodes are the corner nodes, so when n_value==6 or 7, the fraction is the same as the 1d ...
unsigned required_nvalue(const unsigned &n) const
Number of values (pinned or dofs) required at node n. Bumped up to 8 so we don't have to worry if a h...
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 8 (6 velocities + 2 pressures)
RefineableLinearisedQTaylorHoodElement()
Constructor:
Node * get_interpolating_node_at_local_coordinate(const Vector< double > &s, const int &n_value)
The velocity nodes are the same as the geometric nodes. The pressure nodes must be calculated by usin...
void unpin_elemental_pressure_dofs()
Unpin all pressure dofs.
Node * pressure_node_pt(const unsigned &n_p)
Pointer to n_p-th pressure node.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes....
Node * interpolating_node_pt(const unsigned &n, const int &n_value)
The velocities are isoparametric and so the "nodes" interpolating the velocities are the geometric no...
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
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 pin_elemental_redundant_nodal_pressure_dofs()
Unpin the proper nodal pressure dofs.
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...
unsigned nvertex_node() const
Number of vertex nodes in the element.
unsigned ninterpolating_node_1d(const int &n_value)
The number of 1d pressure nodes is 2, the number of 1d velocity nodes is the same as the number of 1d...
unsigned ninterpolating_node(const int &n_value)
The number of pressure nodes is 4. The number of velocity nodes is the same as the number of geometri...
void interpolating_basis(const Vector< double > &s, Shape &psi, const int &n_value) const
The basis interpolating the pressure is given by pshape(). / The basis interpolating the velocity is ...
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
QuadTree * quadtree_pt()
Pointer to quadtree representation of this element.
void setup_hang_for_value(const int &value_id)
Internal helper function that is used to construct the hanging node schemes for the value_id-th inter...
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...
RefineableElement * object_pt() const
Return the pointer to the object (RefineableElement) represented by the tree.
Tree * father_pt() const
Return pointer to father: NULL if it's a root node.
int son_type() const
Return son type.
Tree * son_pt(const int &son_index) const
Return pointer to the son for a given index. Note that to aid code readability specific enums have be...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...