28 #ifndef OOMPH_SPHERICAL_NAVIER_STOKES_ELEMENTS_HEADER
29 #define OOMPH_SPHERICAL_NAVIER_STOKES_ELEMENTS_HEADER
33 #include <oomph-lib-config.h>
38 #include "../generic/Qelements.h"
39 #include "../generic/fsi.h"
139 DShape& dtestdx)
const = 0;
149 DShape& dtestdx)
const = 0;
153 Shape& psi)
const = 0;
159 Shape& test)
const = 0;
176 for (
unsigned i = 0;
i < 3;
i++)
184 (*Body_force_fct_pt)(time, x, result);
218 inline double cot(
const double& th)
const
220 return (1 / tan(th));
249 const double&
re()
const
376 const unsigned&
i)
const
412 for (
unsigned t = 0;
t < n_time;
t++)
443 virtual void fix_pressure(
const unsigned& p_dof,
const double& p_value) = 0;
512 const unsigned& which_one = 0);
525 void output(std::ostream& outfile,
const unsigned& nplot);
537 void output(FILE* file_pt,
const unsigned& nplot);
551 void full_output(std::ostream& outfile,
const unsigned& nplot);
558 const unsigned& nplot,
570 const unsigned& nplot,
577 const unsigned& nplot,
605 std::ostream& outfile,
623 const double Re = this->
re();
629 ans[2] = r * sin(theta);
630 ans[3] = 0.5 * Re * r * r * sin(theta) * sin(theta);
639 const double Re = this->
re();
646 ans[3] = Re * r * sin(theta) * sin(theta);
655 const double Re = this->
re();
661 ans[2] = r * cos(theta);
662 ans[3] = Re * r * r * sin(theta) * cos(theta);
698 residuals, jacobian, mass_matrix, 2);
706 unsigned n_node =
nnode();
712 for (
unsigned i = 0;
i < 3;
i++)
719 for (
unsigned l = 0; l < n_node; l++)
728 const unsigned&
i)
const
731 unsigned n_node =
nnode();
741 double interpolated_u = 0.0;
743 for (
unsigned l = 0; l < n_node; l++)
745 interpolated_u +=
nodal_value(l, u_nodal_index) * psi[l];
748 return (interpolated_u);
755 const unsigned& j)
const
758 unsigned n_node =
nnode();
769 double interpolated_dudx = 0.0;
771 for (
unsigned l = 0; l < n_node; l++)
773 interpolated_dudx +=
nodal_value(l, u_nodal_index) * dpsidx(l, j);
776 return (interpolated_dudx);
790 double interpolated_p = 0.0;
792 for (
unsigned l = 0; l < n_pres; l++)
797 return (interpolated_p);
907 std::set<std::pair<Data*, unsigned>>& paired_load_data);
917 std::set<std::pair<Data*, unsigned>>& paired_pressure_data);
926 void output(std::ostream& outfile,
const unsigned& nplot)
939 void output(FILE* file_pt,
const unsigned& nplot)
976 std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list)
const;
997 for (
unsigned i = 0;
i < 9;
i++)
1000 dtestdx(
i, 0) = dpsidx(
i, 0);
1001 dtestdx(
i, 1) = dpsidx(
i, 1);
1023 for (
unsigned i = 0;
i < 9;
i++)
1026 dtestdx(
i, 0) = dpsidx(
i, 0);
1027 dtestdx(
i, 1) = dpsidx(
i, 1);
1052 for (
unsigned i = 0;
i < 3;
i++) test[
i] = psi[
i];
1116 const unsigned& ipt,
1187 std::set<std::pair<Data*, unsigned>>& paired_load_data);
1198 std::set<std::pair<Data*, unsigned>>& paired_pressure_data);
1207 void output(std::ostream& outfile,
const unsigned& nplot)
1219 void output(FILE* file_pt,
const unsigned& nplot)
1239 std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list)
const;
1261 for (
unsigned i = 0;
i < 9;
i++)
1264 dtestdx(
i, 0) = dpsidx(
i, 0);
1265 dtestdx(
i, 1) = dpsidx(
i, 1);
1288 for (
unsigned i = 0;
i < 9;
i++)
1291 dtestdx(
i, 0) = dpsidx(
i, 0);
1292 dtestdx(
i, 1) = dpsidx(
i, 1);
1305 double psi1[2], psi2[2];
1312 for (
unsigned i = 0;
i < 2;
i++)
1314 for (
unsigned j = 0; j < 2; j++)
1317 psi[2 *
i + j] = psi2[
i] * psi1[j];
1332 for (
unsigned i = 0;
i < 4;
i++) test[
i] = psi[
i];
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.
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
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 ...
/////////////////////////////////////////////////////////////////////////
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
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 double dshape_eulerian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidx) const
Return the geometric shape functions and also first derivatives w.r.t. global coordinates at the ipt-...
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...
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
unsigned nnode() const
Return the number of nodes.
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Compute the geometric shape functions and also first derivatives w.r.t. global coordinates at local c...
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as .
unsigned add_internal_data(Data *const &data_pt, const bool &fd=true)
Add a (pointer to an) internal data object to the element and return the index required to obtain it ...
Data *& internal_data_pt(const unsigned &i)
Return a pointer to i-th internal data object.
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
int internal_local_eqn(const unsigned &i, const unsigned &j) const
Return the local equation number corresponding to the j-th value stored at the i-th internal data.
TimeStepper *& time_stepper_pt()
Access function for pointer to time stepper: Null if object is not time-dependent.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
/////////////////////////////////////////////////////////////////////////// /////////////////////////...
void identify_load_data(std::set< std::pair< Data *, unsigned >> &paired_load_data)
Add to the set paired_load_data pairs containing.
void full_output(std::ostream &outfile, const unsigned &nplot)
Full output function: x,y,[z],u,v,[w],p,du/dt,dv/dt,[dw/dt],dissipation in tecplot format....
void output(std::ostream &outfile, const unsigned &nplot)
Redirect output to SphericalNavierStokesEquations output.
void output(std::ostream &outfile)
Redirect output to SphericalNavierStokesEquations output.
void output(FILE *file_pt, const unsigned &nplot)
Redirect output to SphericalNavierStokesEquations output.
int p_local_eqn(const unsigned &n) const
Return the local equation numbers for the pressure values.
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.
unsigned required_nvalue(const unsigned &n) const
Number of values (pinned or dofs) required at local node n.
double dshape_and_dtest_eulerian_at_knot_spherical_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...
void full_output(std::ostream &outfile)
Full output function: x,y,[z],u,v,[w],p,du/dt,dv/dt,[dw/dt],dissipation in tecplot format....
double p_spherical_nst(const unsigned &i) const
Return the i-th pressure value (Discontinous pressure interpolation – no need to cater for hanging no...
void output(FILE *file_pt)
Redirect output to SphericalNavierStokesEquations output.
QSphericalCrouzeixRaviartElement()
Constructor, there are 3 internal values (for the pressure)
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contain...
void identify_pressure_data(std::set< std::pair< Data *, unsigned >> &paired_pressure_data)
Add to the set paired_pressure_data pairs containing.
static const unsigned Initial_Nvalue[]
Static array of ints to hold required number of variables at nodes.
unsigned npres_spherical_nst() const
Return number of pressure values.
void pshape_spherical_nst(const Vector< double > &s, Shape &psi) const
Pressure shape functions at local coordinate s.
unsigned ndof_types() const
The number of "DOF types" that degrees of freedom in this element are sub-divided into: Velocity (thr...
unsigned P_spherical_nst_internal_index
Internal index that indicates at which internal data the pressure is stored.
double dshape_and_dtest_eulerian_spherical_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...
/////////////////////////////////////////////////////////////////////////
QSphericalTaylorHoodElement()
Constructor, no internal data points.
double p_spherical_nst(const unsigned &n_p) const
Access function for the pressure values at local pressure node n_p (const version)
unsigned ndof_types() const
The number of "DOF types" that degrees of freedom in this element are sub-divided into: Velocity (3 c...
int p_local_eqn(const unsigned &n) const
Return the local equation numbers for the pressure values.
int p_nodal_index_spherical_nst() const
Set the value at which the pressure is stored in the nodes In this case the third index because there...
virtual unsigned required_nvalue(const unsigned &n) const
Number of values (pinned or dofs) required at node n. Can be overwritten for hanging node version.
void output(FILE *file_pt)
Redirect output to SphericalNavierStokesEquations output.
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contain...
static const unsigned Initial_Nvalue[]
Static array of ints to hold number of variables at node.
void output(std::ostream &outfile, const unsigned &nplot)
Redirect output to SphericalNavierStokesEquations output.
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.
void pshape_spherical_nst(const Vector< double > &s, Shape &psi) const
Pressure shape functions at local coordinate s.
double dshape_and_dtest_eulerian_at_knot_spherical_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 local coordinate s (tak...
double dshape_and_dtest_eulerian_spherical_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 output(FILE *file_pt, const unsigned &nplot)
Redirect output to SphericalNavierStokesEquations output.
void identify_load_data(std::set< std::pair< Data *, unsigned >> &paired_load_data)
Add to the set paired_load_data pairs containing.
static const unsigned Pconv[]
Static array of ints to hold conversion from pressure node numbers to actual node numbers.
void output(std::ostream &outfile)
Redirect output to SphericalNavierStokesEquations output.
void identify_pressure_data(std::set< std::pair< Data *, unsigned >> &paired_pressure_data)
Add to the set paired_pressure_data pairs containing.
unsigned npres_spherical_nst() const
Return number of pressure values.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
A class for elements that solve the Navier–Stokes equations, in axisymmetric spherical polar coordina...
SphericalNavierStokesSourceFctPt & source_fct_pt()
Access function for the source-function pointer.
double d_kin_energy_dt() const
Get integral of time derivative of kinetic energy over element.
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Validate against exact solution at given time Solution is provided via function pointer....
void get_vorticity(const Vector< double > &s, Vector< double > &vorticity) const
Compute the vorticity vector at local coordinate s.
SphericalNavierStokesSourceFctPt Source_fct_pt
Pointer to volumetric source function.
double interpolated_dudx_spherical_nst(const Vector< double > &s, const unsigned &i, const unsigned &j) const
Return matrix entry dudx(i,j) of the FE interpolated velocity derivative at local coordinate s.
virtual double dshape_and_dtest_eulerian_at_knot_spherical_nst(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Compute the shape functions and derivatives w.r.t. global coords at ipt-th integration point Return J...
double pressure_integral() const
Integral of pressure over element.
virtual void pshape_spherical_nst(const Vector< double > &s, Shape &psi) const =0
Compute the pressure shape functions at local coordinate s.
void compute_error_e(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, FiniteElement::SteadyExactSolutionFctPt exact_soln_dr_pt, FiniteElement::SteadyExactSolutionFctPt exact_soln_dtheta_pt, double &u_error, double &u_norm, double &p_error, double &p_norm)
Validate against exact solution. Solution is provided direct from exact_soln function....
void output(std::ostream &outfile)
Output function: x,y,[z],u,v,[w],p in tecplot format. Default number of plot points.
void full_output(std::ostream &outfile)
Full output function: x,y,[z],u,v,[w],p,du/dt,dv/dt,[dw/dt],dissipation in tecplot format....
virtual unsigned u_index_spherical_nst(const unsigned &i) const
Return the index at which the i-th unknown velocity component is stored. The default value,...
double * Density_Ratio_pt
Pointer to the density ratio (relative to the density used in the definition of the Reynolds number)
const double & re() const
Reynolds number.
SphericalNavierStokesSourceFctPt source_fct_pt() const
Access function for the source-function pointer. Const version.
virtual double dshape_and_dtest_eulerian_spherical_nst(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Compute the shape functions and derivatives w.r.t. global coords at local coordinate s....
static double Default_Physical_Ratio_Value
Static default value for the physical ratios (all are initialised to one)
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,...
double * ReSt_pt
Pointer to global Reynolds number x Strouhal number (=Womersley)
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)
static Vector< double > Gamma
Vector to decide whether the stress-divergence form is used or not.
double *& viscosity_ratio_pt()
Pointer to Viscosity Ratio.
const double & re_invro() const
Global Reynolds number multiplied by inverse Rossby number.
double u_spherical_nst(const unsigned &n, const unsigned &i) const
Velocity i at local node n. Uses suitably interpolated value for hanging nodes. The use of u_index_sp...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Add the elemental contribution to the jacobian matrix. and the residuals vector. Note that this funct...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Compute the element's residual Vector.
virtual double p_spherical_nst(const unsigned &n_p) const =0
Pressure at local pressure "node" n_p Uses suitably interpolated value for hanging nodes.
double *& re_invfr_pt()
Pointer to global inverse Froude number.
const double & density_ratio() const
Density ratio for element: Element's density relative to the viscosity used in the definition of the ...
Vector< double > * G_pt
Pointer to global gravity Vector.
virtual void get_body_force_spherical_nst(const double &time, const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &result)
Calculate the body force at a given time and local and/or Eulerian position. This function is virtual...
static int Pressure_not_stored_at_node
Static "magic" number that indicates that the pressure is not stored at a node.
double du_dt_spherical_nst(const unsigned &n, const unsigned &i) const
i-th component of du/dt at local node n. Uses suitably interpolated value for hanging nodes.
double * ReInvRo_pt
Pointer to global Reynolds number x inverse Rossby number (used when in a rotating frame)
double *& re_pt()
Pointer to Reynolds number.
void output(FILE *file_pt)
C-style output function: x,y,[z],u,v,[w],p in tecplot format. Default number of plot points.
virtual int p_nodal_index_spherical_nst() const
Return the index at which the pressure is stored if it is stored at the nodes. If not stored at the n...
double * Re_pt
Pointer to global Reynolds number.
double * Viscosity_Ratio_pt
Pointer to the viscosity ratio (relative to the viscosity used in the definition of the Reynolds numb...
Vector< double > actual_dr(const Vector< double > &x)
const double & re_invfr() const
Global inverse Froude number.
virtual void pshape_spherical_nst(const Vector< double > &s, Shape &psi, Shape &test) const =0
Compute the pressure shape and test functions at local coordinate s.
void disable_ALE()
Disable ALE, i.e. assert the mesh is not moving – you do this at your own risk!
SphericalNavierStokesBodyForceFctPt body_force_fct_pt() const
Access function for the body-force pointer. Const version.
double cot(const double &th) const
Include a cot function to simplify the NS momentum and jacobian expressions.
double *& re_invro_pt()
Pointer to global inverse Froude number.
double interpolated_p_spherical_nst(const Vector< double > &s) const
Return FE interpolated pressure at local coordinate s.
void compute_shear_stress(std::ostream &outfile)
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact solution specified via function pointer at a given number of plot points....
double kin_energy() const
Get integral of kinetic energy over element.
void output_vorticity(std::ostream &outfile, const unsigned &nplot)
Output function: x,y,[z], [omega_x,omega_y,[and/or omega_z]] in tecplot format. nplot points in each ...
double dissipation() const
Return integral of dissipation over element.
static double Default_Physical_Constant_Value
Static default value for the physical constants (all initialised to zero)
void get_traction(const Vector< double > &s, const Vector< double > &N, Vector< double > &traction)
Compute traction (on the viscous scale) exerted onto the fluid at local coordinate s....
double(* SphericalNavierStokesSourceFctPt)(const double &time, const Vector< double > &x)
Function pointer to source function fct(t,x) x is a Vector!
Vector< double > actual_dth(const Vector< double > &x)
const Vector< double > & g() const
Vector of gravitational components.
double * ReInvFr_pt
Pointer to global Reynolds number x inverse Froude number (= Bond number / Capillary number)
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
Add the element's contribution to its residuals vector, jacobian matrix and mass matrix.
virtual void fill_in_generic_residual_contribution_spherical_nst(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Compute the residuals for the Navier–Stokes equations; flag=1(or 0): do (or don't) compute the Jacobi...
static Vector< double > Default_Gravity_vector
Static default value for the gravity vector.
SphericalNavierStokesBodyForceFctPt Body_force_fct_pt
Pointer to body force function.
double *& re_st_pt()
Pointer to product of Reynolds and Strouhal number (=Womersley number)
double u_spherical_nst(const unsigned &t, const unsigned &n, const unsigned &i) const
Velocity i at local node n at timestep t (t=0: present; t>0: previous). Uses suitably interpolated va...
virtual int p_local_eqn(const unsigned &n) const =0
Access function for the local equation number information for the pressure. p_local_eqn[n] = local eq...
SphericalNavierStokesEquations()
Constructor: NULL the body force and source function and make sure the ALE terms are included by defa...
virtual unsigned npres_spherical_nst() const =0
Function to return number of pressure degrees of freedom.
void extract_velocity(std::ostream &outfile)
virtual void fix_pressure(const unsigned &p_dof, const double &p_value)=0
Pin p_dof-th pressure dof and set it to value specified by p_value.
double interpolated_u_spherical_nst(const Vector< double > &s, const unsigned &i) const
Return FE interpolated velocity u[i] at local coordinate s.
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when time-derivatives are computed....
void output_veloc(std::ostream &outfile, const unsigned &nplot, const unsigned &t)
Output function: x,y,[z],u,v,[w] in tecplot format. nplot points in each coordinate direction at time...
void get_load(const Vector< double > &s, const Vector< double > &N, Vector< double > &load)
This implements a pure virtual function defined in the FSIFluidElement class. The function computes t...
void(* SphericalNavierStokesBodyForceFctPt)(const double &time, const Vector< double > &x, Vector< double > &body_force)
Function pointer to body force function fct(t,x,f(x)) x is a Vector!
void enable_ALE()
(Re-)enable ALE, i.e. take possible mesh motion into account when evaluating the time-derivative....
const double & re_st() const
Product of Reynolds and Strouhal number (=Womersley number)
void interpolated_u_spherical_nst(const Vector< double > &s, Vector< double > &veloc) const
Compute vector of FE interpolated velocity u at local coordinate s.
Vector< double > *& g_pt()
Pointer to Vector of gravitational components.
SphericalNavierStokesBodyForceFctPt & body_force_fct_pt()
Access function for the body-force pointer.
virtual double get_source_spherical_nst(double time, const unsigned &ipt, const Vector< double > &x)
Calculate the source fct at given time and Eulerian position.
const double & viscosity_ratio() const
Viscosity ratio for element: Element's viscosity relative to the viscosity used in the definition of ...
Vector< double > actual(const Vector< double > &x)
double *& density_ratio_pt()
Pointer to Density ratio.
////////////////////////////////////////////////////////////////////// //////////////////////////////...
unsigned ntstorage() const
Return the number of doubles required to represent history (one for steady)
virtual double weight(const unsigned &i, const unsigned &j) const
Access function for j-th weight for the i-th derivative.
std::string type() const
Return string that indicates the type of the timestepper (e.g. "BDF", "Newmark", etc....
void shape< 2 >(const double &s, double *Psi)
1D shape functions specialised to linear order (2 Nodes)
//////////////////////////////////////////////////////////////////// ////////////////////////////////...