28#ifndef OOMPH_GENERALISED_NEWTONIAN_REFINEABLE_NAVIER_STOKES_ELEMENTS_HEADER
29#define OOMPH_GENERALISED_NEWTONIAN_REFINEABLE_NAVIER_STOKES_ELEMENTS_HEADER
33#include <oomph-lib-config.h>
55 template<
unsigned DIM>
149 std::ostringstream error_message;
150 error_message <<
"The flux vector has the wrong number of entries, "
151 << flux.
size() <<
", whereas it should be at least "
167 for (
unsigned i = 0;
i <
DIM;
i++)
174 for (
unsigned i = 0;
i <
DIM;
i++)
176 for (
unsigned j =
i + 1;
j <
DIM;
j++)
198 this->
Re_pt = cast_father_element_pt->re_pt();
200 this->
ReSt_pt = cast_father_element_pt->re_st_pt();
202 this->
ReInvFr_pt = cast_father_element_pt->re_invfr_pt();
204 this->
G_pt = cast_father_element_pt->g_pt();
290 global_eqn_number.resize(
n_u_dof, 0);
381 template<
unsigned DIM>
481 values.resize(
DIM + 1, 0.0);
484 for (
unsigned i = 0;
i <
DIM;
i++)
502 values.resize(
DIM + 1);
505 for (
unsigned i = 0;
i <
DIM + 1;
i++)
518 for (
unsigned i = 0;
i <
DIM;
i++)
602 for (
unsigned i = 0;
i <
DIM;
i++)
611 else if (
s[
i] == 1.0)
634 index[
i] *
static_cast<unsigned>(
pow(
static_cast<float>(
NNODE_1D),
635 static_cast<int>(
i)));
668 return static_cast<unsigned>(
pow(2.0,
static_cast<int>(
DIM)));
672 return this->
nnode();
706 unsigned u_index[
DIM];
707 for (
unsigned i = 0;
i <
DIM;
i++)
723 unsigned nmaster =
nod_pt->hanging_pt()->nmaster();
726 for (
unsigned j = 0;
j < nmaster;
j++)
732 for (
unsigned i = 0;
i <
DIM;
i++)
744 for (
unsigned i = 0;
i <
DIM;
i++)
747 std::make_pair(this->
node_pt(
n), u_index[
i]));
772 for (
unsigned m = 0;
m < nmaster;
m++)
796 template<
unsigned DIM>
798 :
public virtual FaceGeometry<GeneralisedNewtonianQTaylorHoodElement<DIM>>
812 template<
unsigned DIM>
816 FaceGeometry<GeneralisedNewtonianQTaylorHoodElement<DIM>>>
835 template<
unsigned DIM>
917 values.resize(
DIM, 0.0);
920 for (
unsigned i = 0;
i <
DIM;
i++)
943 for (
unsigned i = 0;
i <
DIM;
i++)
956 for (
unsigned i = 0;
i <
DIM;
i++)
991 unsigned u_index[
DIM];
992 for (
unsigned i = 0;
i <
DIM;
i++)
1005 if (
nod_pt->is_hanging())
1008 unsigned nmaster =
nod_pt->hanging_pt()->nmaster();
1011 for (
unsigned j = 0;
j < nmaster;
j++)
1017 for (
unsigned i = 0;
i <
DIM;
i++)
1029 for (
unsigned i = 0;
i <
DIM;
i++)
1032 std::make_pair(this->
node_pt(
n), u_index[
i]));
1055 template<
unsigned DIM>
1083 this->p_order() = 3;
1131 double p_nst(
const unsigned&
t,
const unsigned&
i)
const
1139 return (this->p_order() - 2) * (this->p_order() - 2);
1237 values.resize(
DIM, 0.0);
1240 for (
unsigned i = 0;
i <
DIM;
i++)
1263 for (
unsigned i = 0;
i <
DIM;
i++)
1276 for (
unsigned i = 0;
i <
DIM;
i++)
1302 template<
unsigned DIM>
1305 GeneralisedNewtonianQCrouzeixRaviartElement<DIM>>
1320 template<
unsigned DIM>
1324 FaceGeometry<GeneralisedNewtonianQCrouzeixRaviartElement<DIM>>>
1341 inline void RefineableGeneralisedNewtonianQCrouzeixRaviartElement<
1342 2>::rebuild_from_sons(
Mesh*& mesh_pt)
1344 using namespace QuadTreeNames;
1369 internal_data_pt(this->P_nst_internal_index)->set_value(0, 0.25 *
av_press);
1380 double slope1 = quadtree_pt()
1391 double slope2 = quadtree_pt()
1404 internal_data_pt(this->P_nst_internal_index)
1440 internal_data_pt(this->P_nst_internal_index)
1450 3>::rebuild_from_sons(
Mesh*& mesh_pt)
1452 using namespace OcTreeNames;
1475 internal_data_pt(this->P_nst_internal_index)
1487 double slope1 = octree_pt()
1498 double slope2 = octree_pt()
1509 double slope3 = octree_pt()
1520 double slope4 = octree_pt()
1533 internal_data_pt(this->P_nst_internal_index)
1591 internal_data_pt(this->P_nst_internal_index)
1648 internal_data_pt(this->P_nst_internal_index)
1666 using namespace QuadTreeNames;
1669 int son_type = quadtree_pt()->son_type();
1685 else if (son_type == SE)
1691 else if (son_type == NE)
1698 else if (son_type == NW)
1713 internal_data_pt(this->P_nst_internal_index)->set_value(0,
press);
1716 for (
unsigned i = 1;
i < 3;
i++)
1723 internal_data_pt(this->P_nst_internal_index)
1741 using namespace OcTreeNames;
1744 int son_type = octree_pt()->son_type();
1748 octree_pt()->father_pt()->object_pt());
1753 for (
unsigned i = 0;
i < 3;
i++)
1767 internal_data_pt(this->P_nst_internal_index)->set_value(0,
press);
1770 for (
unsigned i = 1;
i < 4;
i++)
1777 internal_data_pt(this->P_nst_internal_index)
1801 for (
unsigned i = 0;
i < nnode_1d() * nnode_1d();
i++)
1820 2>::dshape_and_dtest_eulerian_at_knot_nst(
const unsigned&
ipt,
1831 for (
unsigned i = 0;
i < nnode_1d() * nnode_1d();
i++)
1861 for (
unsigned i = 0;
i < nnode_1d() * nnode_1d() * nnode_1d();
i++)
1881 3>::dshape_and_dtest_eulerian_at_knot_nst(
const unsigned&
ipt,
1892 for (
unsigned i = 0;
i < nnode_1d() * nnode_1d() * nnode_1d();
i++)
1912 unsigned npres = this->npres_nst();
1949 if (this->npres_nst() == 1)
1955 for (
unsigned i = 0;
i < this->npres_nst();
i++)
test[
i] =
psi[
i];
1967 unsigned npres = this->npres_nst();
2009 if (this->npres_nst() == 1)
2015 for (
unsigned i = 0;
i < this->npres_nst();
i++)
test[
i] =
psi[
i];
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when time-derivatives are computed....
AdvectionDiffusionReactionSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
A class that represents a collection of data; each Data object may contain many different individual ...
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...
unsigned nvalue() const
Return number of values stored in data object (incl pinned ones).
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 ...
long & eqn_number(const unsigned &i)
Return the equation number of the i-th stored variable.
virtual void resize(const unsigned &n_value)
Change (increase) the number of values that may be stored.
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
FaceGeometry class definition: This policy class is used to allow construction of face elements that ...
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
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.
double size() const
Calculate the size of the element (length, area, volume,...) in Eulerian computational coordinates....
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.
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
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 void set_integration_scheme(Integral *const &integral_pt)
Set the spatial integration scheme.
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 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...
unsigned long eqn_number(const unsigned &ieqn_local) const
Return the global equation number corresponding to the ieqn_local-th local equation number.
Data *& internal_data_pt(const unsigned &i)
Return a pointer to i-th internal data object.
A class for elements that solve the cartesian Navier–Stokes equations, templated by the dimension DIM...
double * ReInvFr_pt
Pointer to global Reynolds number x inverse Froude number (= Bond number / Capillary number)
double * Viscosity_Ratio_pt
Pointer to the viscosity ratio (relative to the viscosity used in the definition of the Reynolds numb...
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when time-derivatives are computed....
void strain_rate(const Vector< double > &s, DenseMatrix< double > &strain_rate) const
Strain-rate tensor: 1/2 (du_i/dx_j + du_j/dx_i)
double * Density_Ratio_pt
Pointer to the density ratio (relative to the density used in the definition of the Reynolds number)
virtual unsigned u_index_nst(const unsigned &i) const
Return the index at which the i-th unknown velocity component is stored. The default value,...
void interpolated_u_nst(const Vector< double > &s, Vector< double > &veloc) const
Compute vector of FE interpolated velocity u at local coordinate s.
NavierStokesBodyForceFctPt Body_force_fct_pt
Pointer to body force function.
NavierStokesSourceFctPt Source_fct_pt
Pointer to volumetric source function.
double * Re_pt
Pointer to global Reynolds number.
double * ReSt_pt
Pointer to global Reynolds number x Strouhal number (=Womersley)
virtual double interpolated_p_nst(const Vector< double > &s) const
Return FE interpolated pressure at local coordinate s.
Vector< double > * G_pt
Pointer to global gravity Vector.
Crouzeix_Raviart elements are Navier–Stokes elements with quadratic interpolation for velocities and ...
unsigned P_nst_internal_index
Internal index that indicates at which internal data the pressure is stored.
unsigned npres_nst() const
Return number of pressure values.
Taylor–Hood elements are Navier–Stokes elements with quadratic interpolation for velocities and posit...
virtual int p_nodal_index_nst() const
Set the value at which the pressure is stored in the nodes.
void pshape_nst(const Vector< double > &s, Shape &psi) const
Pressure shape functions at local coordinate s.
static const unsigned Pconv[]
Static array of ints to hold conversion from pressure node numbers to actual node numbers.
unsigned npres_nst() const
Return number of pressure values.
Class that contains data for hanging nodes.
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.
HangInfo *const & hanging_pt() const
Return pointer to hanging node data (this refers to the geometric hanging node status) (const version...
static Vector< Vector< int > > Direction_to_vector
For each direction, i.e. a son_type (vertex), a face or an edge, this defines a vector that indicates...
Non-templated class that returns modal hierachical shape functions based on Legendre polynomials.
An OomphLibError object which should be thrown when an run-time error is encountered....
p-refineable version of Crouzeix Raviart elements. Generic class definitions
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation:
void pshape_nst(const Vector< double > &s, Shape &psi) const
Pressure shape functions at local coordinate s.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes....
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
PRefineableGeneralisedNewtonianQCrouzeixRaviartElement()
Constructor.
double dshape_and_dtest_eulerian_at_knot_nst(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Velocity shape and test functions and their derivs w.r.t. to global coords at ipt-th integation point...
double p_nst(const unsigned &i) const
Broken assignment operator.
double p_nst(const unsigned &t, const unsigned &i) const
Return the i-th pressure value (Discontinous pressure interpolation – no need to cater for hanging no...
void unpin_elemental_pressure_dofs()
Unpin all internal pressure dofs.
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: Reconstruct pressure from the (merged) sons This must be specialised for each dime...
void further_build()
Further build for Crouzeix_Raviart interpolates the internal pressure dofs from father element: Make ...
double dshape_and_dtest_eulerian_nst(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Velocity shape and test functions and their derivs w.r.t. to global coords at local coordinate s (tak...
void pshape_nst(const Vector< double > &s, Shape &psi, Shape &test) const
Pressure shape and test functions at local coordinte s.
unsigned npres_nst() const
// Return number of pressure values
~PRefineableGeneralisedNewtonianQCrouzeixRaviartElement()
Destructor.
unsigned required_nvalue(const unsigned &n) const
Number of values (pinned or dofs) required at local node n.
void fix_pressure(const unsigned &p_dof, const double &p_value)
Pin p_dof-th pressure dof and set it to value specified by p_value.
PRefineableGeneralisedNewtonianQCrouzeixRaviartElement(const PRefineableGeneralisedNewtonianQCrouzeixRaviartElement< DIM > &dummy)=delete
Broken copy constructor.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: DIM (velocities)
unsigned nvertex_node() const
Number of vertex nodes in the element.
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Get all function values [u,v..,p] at previous timestep t (t=0: present; t>0: previous timestep).
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...
A class that is used to template the p-refineable Q elements by dimension. It's really nothing more t...
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 Navier–Stokes equations.
void dinterpolated_u_nst_ddata(const Vector< double > &s, const unsigned &i, 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...
virtual void get_dresidual_dnodal_coordinates(RankThreeTensor< double > &dresidual_dnodal_coordinates)
Compute derivatives of elemental residual vector with respect to nodal coordinates....
static void unpin_all_pressure_dofs(const Vector< GeneralisedElement * > &element_pt)
Unpin all pressure dofs in elements listed in vector.
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 ...
void further_build()
Further build, pass the pointers down to the sons.
virtual void unpin_elemental_pressure_dofs()=0
Unpin all pressure dofs in the element.
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...
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get 'flux' for Z2 error recovery: Upper triangular entries in strain rate tensor.
RefineableGeneralisedNewtonianNavierStokesEquations()
Constructor.
virtual void pin_elemental_redundant_nodal_pressure_dofs()
Pin unused nodal pressure dofs (empty by default, because by default pressure dofs are not associated...
void fill_in_generic_residual_contribution_nst(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Add element's contribution to elemental residual vector and/or Jacobian matrix flag=1: compute both f...
void get_pressure_and_velocity_mass_matrix_diagonal(Vector< double > &press_mass_diag, Vector< double > &veloc_mass_diag, const unsigned &which_one=0)
Compute the diagonal of the velocity/pressure mass matrices. If which one=0, both are computed,...
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
Refineable version of Crouzeix Raviart elements. Generic class definitions.
RefineableGeneralisedNewtonianQCrouzeixRaviartElement(const RefineableGeneralisedNewtonianQCrouzeixRaviartElement< DIM > &dummy)=delete
Broken copy constructor.
void identify_load_data(std::set< std::pair< Data *, unsigned > > &paired_load_data)
Add to the set paired_load_data pairs containing.
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: Reconstruct pressure from the (merged) sons This must be specialised for each dime...
void unpin_elemental_pressure_dofs()
Unpin all internal pressure dofs.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Get all function values [u,v..,p] at previous timestep t (t=0: present; t>0: previous timestep).
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes....
RefineableGeneralisedNewtonianQCrouzeixRaviartElement()
Constructor.
unsigned ncont_interpolated_values() const
Broken assignment operator.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
void further_build()
Further build for Crouzeix_Raviart interpolates the internal pressure dofs from father element: Make ...
unsigned nvertex_node() const
Number of vertex nodes in the element.
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...
Refineable version of Taylor Hood elements. These classes can be written in total generality.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: (DIM velocities + 1 pressure)
unsigned ninterpolating_node(const int &value_id)
The number of pressure nodes is 2^DIM. The number of velocity nodes is the same as the number of geom...
void identify_load_data(std::set< std::pair< Data *, unsigned > > &paired_load_data)
Add to the set paired_load_data pairs containing.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes....
unsigned nvertex_node() const
Number of vertex nodes in the element.
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...
Node * interpolating_node_pt(const unsigned &n, const int &value_id)
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 pin_elemental_redundant_nodal_pressure_dofs()
Pin all nodal pressure dofs that are not required.
Node * get_interpolating_node_at_local_coordinate(const Vector< double > &s, const int &value_id)
The velocity nodes are the same as the geometric nodes. The pressure nodes must be calculated by usin...
unsigned required_nvalue(const unsigned &n) const
Number of values required at local node n. In order to simplify matters, we allocate storage for pres...
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
RefineableGeneralisedNewtonianQTaylorHoodElement()
Constructor.
void unpin_elemental_pressure_dofs()
Unpin all pressure dofs.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
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 ninterpolating_node_1d(const int &value_id)
The number of 1d pressure nodes is 2, the number of 1d velocity nodes is the same as the number of 1d...
Node * pressure_node_pt(const unsigned &n_p)
Pointer to n_p-th pressure node.
void interpolating_basis(const Vector< double > &s, Shape &psi, const int &value_id) const
The basis interpolating the pressure is given by pshape(). / The basis interpolating the velocity is ...
double local_one_d_fraction_of_interpolating_node(const unsigned &n1d, const unsigned &i, const int &value_id)
The pressure nodes are the corner nodes, so when n_value==DIM, the fraction is the same as the 1d nod...
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...
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).