32 #ifndef OOMPH_AXISYMMETRIC_NAVIER_STOKES_ELEMENTS_HEADER 33 #define OOMPH_AXISYMMETRIC_NAVIER_STOKES_ELEMENTS_HEADER 37 #include <oomph-lib-config.h> 41 #include "../generic/Qelements.h" 42 #include "../generic/fsi.h" 43 #include "../generic/projection.h" 326 Shape &test)
const=0;
339 for(
unsigned i=0;
i<3;
i++) {result[
i] = 0.0;}
344 (*Body_force_fct_pt)(time,x,result);
370 for(
unsigned i=0;
i<2;
i++)
374 for(
unsigned j=0;j<3;j++)
376 d_body_force_dx(j,
i)=(body_force_pls[j]-body_force[j])/eps_fd;
417 double source_pls = 0.0;
419 for(
unsigned i=0;
i<2;
i++)
423 gradient[
i] = (source_pls-source)/eps_fd;
458 double*
const ¶meter_pt,
477 ALE_is_disabled(false)
593 const unsigned n_time = time_stepper_pt->
ntstorage();
596 for(
unsigned t=0;
t<n_time;
t++)
609 ALE_is_disabled=
true;
618 ALE_is_disabled=
false;
623 virtual double p_axi_nst(
const unsigned &n_p)
const=0;
660 const unsigned& which_one=0);
673 const unsigned& nplot)
const 680 for (
unsigned iplot=0;iplot<num_plot_points;iplot++)
699 std::stringstream error_stream;
701 <<
"Axisymmetric Navier-Stokes Elements only store 4 fields " 705 OOMPH_CURRENT_FUNCTION,
706 OOMPH_EXCEPTION_LOCATION);
729 std::stringstream error_stream;
731 <<
"Axisymmetric Navier-Stokes Elements only store 4 fields " 735 OOMPH_CURRENT_FUNCTION,
736 OOMPH_EXCEPTION_LOCATION);
746 for(
unsigned i=0;
i<2;
i++)
752 for(
unsigned i=0;
i<3;
i++)
772 void output(std::ostream &outfile,
const unsigned &nplot);
785 void output(FILE* file_pt,
const unsigned &nplot);
791 const unsigned &nplot,
798 const unsigned &nplot,
805 const unsigned &nplot,
816 double& error,
double& norm);
856 residuals,jacobian,mass_matrix,2);
864 dresidual_dnodal_coordinates);
892 double*
const ¶meter_pt,
899 parameter_pt,dres_dparam,djac_dparam,dmass_matrix_dparam,2);
908 unsigned n_node =
nnode();
914 for (
unsigned i=0;
i<3;
i++)
921 for(
unsigned l=0;l<n_node;l++)
930 const unsigned &
i)
const 933 unsigned n_node =
nnode();
943 double interpolated_u = 0.0;
945 for(
unsigned l=0;l<n_node;l++)
947 interpolated_u +=
nodal_value(l,u_nodal_index)*psi[l];
950 return(interpolated_u);
958 const unsigned &
i)
const 961 unsigned n_node =
nnode();
971 double interpolated_u = 0.0;
973 for(
unsigned l=0;l<n_node;l++)
975 interpolated_u +=
nodal_value(t,l,u_nodal_index)*psi[l];
978 return(interpolated_u);
994 unsigned n_node =
nnode();
1005 for(
unsigned l=0;l<n_node;l++)
1009 if(global_eqn >= 0) {++n_u_dof;}
1013 du_ddata.resize(n_u_dof,0.0);
1014 global_eqn_number.resize(n_u_dof,0);
1019 for(
unsigned l=0;l<n_node;l++)
1023 if (global_eqn >= 0)
1026 global_eqn_number[count] = global_eqn;
1028 du_ddata[count] = psi[l];
1047 double interpolated_p = 0.0;
1049 for(
unsigned l=0;l<n_pres;l++)
1054 return(interpolated_p);
1061 const unsigned &j)
const 1064 const unsigned n_node =
nnode();
1069 DShape dpsifds(n_node,2);
1078 double interpolated_duds = 0.0;
1081 for(
unsigned l=0;l<n_node;l++)
1083 interpolated_duds +=
nodal_value(l,u_nodal_index)*dpsifds(l,j);
1086 return(interpolated_duds);
1093 const unsigned &j)
const 1096 const unsigned n_node =
nnode();
1101 DShape dpsifdx(n_node,2);
1110 double interpolated_dudx = 0.0;
1113 for(
unsigned l=0;l<n_node;l++)
1115 interpolated_dudx +=
nodal_value(l,u_nodal_index)*dpsifdx(l,j);
1118 return(interpolated_dudx);
1124 const unsigned &
i)
const 1127 const unsigned n_node =
nnode();
1136 double interpolated_dudt = 0.0;
1139 for(
unsigned l=0;l<n_node;l++)
1144 return(interpolated_dudt);
1154 const unsigned &k)
const 1157 const unsigned n_node =
nnode();
1162 DShape dpsifds(n_node,2);
1189 inverse_jacobian,dpsifds,d_dpsifdx_dX);
1195 double interpolated_d_dudx_dX = 0.0;
1198 for(
unsigned l=0;l<n_node;l++)
1200 interpolated_d_dudx_dX +=
1201 nodal_value(l,u_nodal_index)*d_dpsifdx_dX(p,q,l,k);
1204 return(interpolated_d_dudx_dX);
1211 double result = 0.0;
1224 const unsigned n_dim =
dim();
1228 for(
unsigned ipt=0;ipt<n_intpt;ipt++)
1268 static const unsigned Initial_Nvalue[];
1297 const unsigned &ipt,
1359 void output(std::ostream &outfile,
const unsigned &n_plot)
1368 void output(FILE* file_pt,
const unsigned &n_plot)
1385 std::list<std::pair<unsigned long,unsigned> >& dof_lookup_list)
const;
1406 for(
unsigned i=0;
i<9;
i++)
1409 dtestdx(
i,0) = dpsidx(
i,0);
1410 dtestdx(
i,1) = dpsidx(
i,1);
1430 for(
unsigned i=0;
i<9;
i++)
1433 dtestdx(
i,0) = dpsidx(
i,0);
1434 dtestdx(
i,1) = dpsidx(
i,1);
1458 djacobian_dX,d_dpsidx_dX);
1462 for(
unsigned i=0;
i<9;
i++)
1466 for(
unsigned k=0;k<2;k++)
1468 dtestdx(
i,k) = dpsidx(
i,k);
1470 for(
unsigned p=0;p<2;p++)
1472 for(
unsigned q=0;q<9;q++)
1474 d_dtestdx_dX(p,q,
i,k) = d_dpsidx_dX(p,q,
i,k);
1503 for(
unsigned i=0;
i<3;
i++) test[
i] = psi[
i];
1546 static const unsigned Initial_Nvalue[];
1552 static const unsigned Pconv[];
1575 const unsigned &ipt,
1600 {
return Initial_Nvalue[n];}
1634 void output(std::ostream &outfile,
const unsigned &n_plot)
1642 void output(FILE* file_pt,
const unsigned &n_plot)
1659 std::list<std::pair<unsigned long, unsigned> >& dof_lookup_list)
const;
1681 for(
unsigned i=0;
i<9;
i++)
1684 dtestdx(
i,0) = dpsidx(
i,0);
1685 dtestdx(
i,1) = dpsidx(
i,1);
1706 for(
unsigned i=0;
i<9;
i++)
1709 dtestdx(
i,0) = dpsidx(
i,0);
1710 dtestdx(
i,1) = dpsidx(
i,1);
1734 djacobian_dX,d_dpsidx_dX);
1738 for(
unsigned i=0;
i<9;
i++)
1742 for(
unsigned k=0;k<2;k++)
1744 dtestdx(
i,k) = dpsidx(
i,k);
1746 for(
unsigned p=0;p<2;p++)
1748 for(
unsigned q=0;q<9;q++)
1750 d_dtestdx_dX(p,q,
i,k) = d_dpsidx_dX(p,q,
i,k);
1768 double psi1[2], psi2[2];
1775 for(
unsigned i=0;
i<2;
i++)
1777 for(
unsigned j=0;j<2;j++)
1780 psi[2*
i + j] = psi2[
i]*psi1[j];
1794 for(
unsigned i=0;
i<4;
i++) test[
i] = psi[
i];
1822 template<
class TAYLOR_HOOD_ELEMENT>
1844 unsigned nnod=this->
nnode();
1845 for (
unsigned j=0;j<nnod;j++)
1848 data_values.push_back(std::make_pair(this->
node_pt(j),fld));
1855 unsigned Pconv_size=3;
1856 for (
unsigned j=0;j<Pconv_size;j++)
1859 unsigned vertex_index=this->Pconv[j];
1860 data_values.push_back(std::make_pair(this->
node_pt(vertex_index),fld));
1904 unsigned n_dim=this->
dim();
1905 unsigned n_node=this->
nnode();
1912 Shape psif(n_node),testf(n_node);
1913 DShape dpsifdx(n_node,n_dim), dtestfdx(n_node,n_dim);
1922 Shape testf(n_node);
1923 DShape dpsifdx(n_node,n_dim), dtestfdx(n_node,n_dim);
1937 const unsigned &fld,
1940 unsigned n_node=this->
nnode();
1960 double interpolated_u = 0.0;
1963 for(
unsigned l=0;l<n_node;l++)
1965 interpolated_u += this->
nodal_value(t,l,u_nodal_index)*psi[l];
1967 return interpolated_u;
1982 return this->
nnode();
2009 template<
class ELEMENT>
2022 template<
class ELEMENT>
2034 template<
class CROUZEIX_RAVIART_ELEMENT>
2056 const unsigned n_node=this->
nnode();
2057 for (
unsigned n=0;n<n_node;n++)
2060 data_values.push_back(std::make_pair(this->
node_pt(n),fld));
2069 for(
unsigned j=0;j<n_press;j++)
2071 data_values.push_back(
2116 unsigned n_dim=this->
dim();
2117 unsigned n_node=this->
nnode();
2124 Shape psif(n_node),testf(n_node);
2125 DShape dpsifdx(n_node,n_dim), dtestfdx(n_node,n_dim);
2134 Shape testf(n_node);
2135 DShape dpsifdx(n_node,n_dim), dtestfdx(n_node,n_dim);
2149 const unsigned &fld,
2177 return this->
nnode();
2204 template<
class ELEMENT>
2217 template<
class ELEMENT>
virtual void pshape_axi_nst(const Vector< double > &s, Shape &psi) const =0
Compute the pressure shape functions at local coordinate s.
const double & re_st() const
Product of Reynolds and Strouhal number (=Womersley number)
double *& re_invro_pt()
Pointer to global inverse Froude number.
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Specify the values associated with field fld. The information is returned in a vector of pairs which ...
double interpolated_duds_axi_nst(const Vector< double > &s, const unsigned &i, const unsigned &j) const
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 ...
double *& re_invfr_pt()
Pointer to global inverse Froude number.
unsigned ndof_types() const
The number of "DOF types" that degrees of freedom in this element are sub-divided into: Velocity (3 c...
virtual double dshape_and_dtest_eulerian_at_knot_axi_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 p_axi_nst(const unsigned &i) const
Return the pressure values at internal dof i_internal (Discontinous pressure interpolation – no need...
unsigned P_axi_nst_internal_index
static int Pressure_not_stored_at_node
Static "magic" number that indicates that the pressure is not stored at a node.
void(*&)(const double &time, const Vector< double > &x, Vector< double > &f) axi_nst_body_force_fct_pt()
Access function for the body-force pointer.
virtual unsigned required_nvalue(const unsigned &n) const
Number of values that must be stored at local node n by the element. The default is 0...
virtual void dinterpolated_u_axi_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...
unsigned nfields_for_projection()
Number of fields to be projected: dim+1, corresponding to velocity components and pressure...
unsigned nhistory_values_for_projection(const unsigned &fld)
Number of history values to be stored for fld-th field. Whatever the timestepper has set up for the v...
double get_source_fct(const double &time, const unsigned &ipt, const Vector< double > &x)
Calculate the source fct at given time and Eulerian position.
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Get cector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
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-...
const double & viscosity_ratio() const
Viscosity ratio for element: Element's viscosity relative to the viscosity used in the definition of ...
static Vector< double > Gamma
Vector to decide whether the stress-divergence form is used or not.
double interpolated_u_axi_nst(const unsigned &t, const Vector< double > &s, const unsigned &i) const
Return FE interpolated velocity u[i] at local coordinate s.
virtual void fill_in_generic_residual_contribution_axi_nst(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Compute the residuals for the Navier–Stokes equations; flag=2 or 1 or 0: compute the Jacobian and/or...
virtual void fill_in_generic_dresidual_contribution_axi_nst(double *const ¶meter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam, unsigned flag)
Compute the derivative of residuals for the Navier–Stokes equations; with respect to a parameeter fl...
void pshape_axi_nst(const Vector< double > &s, Shape &psi) const
Pressure shape functions at local coordinate s.
virtual void get_viscosity_ratio_axisym_nst(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, double &visc_ratio)
Calculate the viscosity ratio relative to the viscosity used in the definition of the Reynolds number...
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld.
virtual unsigned nplot_points_paraview(const unsigned &nplot) const
Return the number of actual plot points for paraview plot with parameter nplot. Broken virtual; can b...
void fill_in_contribution_to_djacobian_dparameter(double *const ¶meter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
Compute the element's residual Vector and the jacobian matrix Virtual function can be overloaded by h...
unsigned n_u_nst() const
Return the number of velocity components for use in general FluidInterface clas.
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Specify the values associated with field fld. The information is returned in a vector of pairs which ...
const double Pi
50 digits from maple
static double Default_Physical_Ratio_Value
Static default value for the physical ratios (all are initialised to one)
Crouzeix Raviart upgraded to become projectable.
virtual void get_source_fct_gradient(const double &time, const unsigned &ipt, const Vector< double > &x, Vector< double > &gradient)
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld.
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. Function prints as many components as are returned in solution Vector.
A general Finite Element class.
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Return Jacobian of mapping and shape functions of field fld at local coordinate s.
unsigned nhistory_values_for_coordinate_projection()
Number of positional history values (includes current value!)
virtual double dshape_and_dtest_eulerian_axi_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...
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as ...
void output(std::ostream &outfile)
Redirect output to NavierStokesEquations output.
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when the time-derivatives are computed...
unsigned nscalar_paraview() const
Number of scalars/fields output by this element. Reimplements broken virtual function in base class...
TimeStepper *& position_time_stepper_pt()
Return a pointer to the position timestepper.
virtual void dJ_eulerian_dnodal_coordinates(const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
A template-free interface that calculates the derivative of the jacobian of a mapping with respect to...
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...
const double & re_invro() const
Global Reynolds number multiplied by inverse Rossby number.
void strain_rate(const Vector< double > &s, DenseMatrix< double > &strain_rate) const
Strain-rate tensor: where (in that order)
void output(std::ostream &outfile, const unsigned &n_plot)
Redirect output to NavierStokesEquations output.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Compute the element's residual Vector.
long & eqn_number(const unsigned &i)
Return the equation number of the i-th stored variable.
const double & density_ratio() const
Density ratio for element: Element's density relative to the viscosity used in the definition of the ...
void fill_in_contribution_to_dresiduals_dparameter(double *const ¶meter_pt, Vector< double > &dres_dparam)
Compute the element's residual Vector.
double(* Source_fct_pt)(const double &time, const Vector< double > &x)
Pointer to volumetric source function.
unsigned npres_axi_nst() const
Return number of pressure values.
double kin_energy() const
Get integral of kinetic energy over element.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
double * Viscosity_Ratio_pt
Pointer to the viscosity ratio (relative to the viscosity used in the definition of the Reynolds numb...
void pin(const unsigned &i)
Pin the i-th stored variable.
void disable_ALE()
Disable ALE, i.e. assert the mesh is not moving – you do this at your own risk!
void output(std::ostream &outfile)
Redirect output to NavierStokesEquations output.
double interpolated_d_dudx_dX_axi_nst(const Vector< double > &s, const unsigned &p, const unsigned &q, const unsigned &i, const unsigned &k) const
Return FE interpolated derivatives w.r.t. nodal coordinates X_{pq} of the spatial derivatives of the ...
double pressure_integral() const
Integral of pressure over element.
const double & re() const
Reynolds number.
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
double * ReInvRo_pt
Pointer to global Reynolds number x inverse Rossby number (used when in a rotating frame) ...
virtual double local_to_eulerian_mapping(const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
Calculate the mapping from local to Eulerian coordinates, given the derivatives of the shape function...
unsigned nhistory_values_for_projection(const unsigned &fld)
Number of history values to be stored for fld-th field. Whatever the timestepper has set up for the v...
TimeStepper *& time_stepper_pt()
Access function for pointer to time stepper: Null if object is not time-dependent.
unsigned ndof_types() const
Returns the number of "DOF types" that degrees of freedom in this element are sub-divided into: Veloc...
unsigned nhistory_values_for_coordinate_projection()
Number of positional history values (includes current value!)
Vector< double > *& g_pt()
Pointer to Vector of gravitational components.
AxisymmetricQCrouzeixRaviartElement()
Constructor, there are three internal values (for the pressure)
std::string scalar_name_paraview(const unsigned &i) const
Name of the i-th scalar field. Default implementation returns V1 for the first one, V2 for the second etc. Can (should!) be overloaded with more meaningful names in specific elements.
virtual void get_body_force_axi_nst(const double &time, const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &result)
Calculate the body force fct at a given time and Eulerian position.
unsigned nodal_dimension() const
Return the required Eulerian dimension of the nodes in this element.
double *& re_st_pt()
Pointer to product of Reynolds and Strouhal number (=Womersley number)
virtual void get_body_force_gradient_axi_nst(const double &time, const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, DenseMatrix< double > &d_body_force_dx)
Get gradient of body force term at (Eulerian) position x. Computed via function pointer (if set) or b...
void(* Body_force_fct_pt)(const double &time, const Vector< double > &x, Vector< double > &result)
Pointer to body force function.
void fix_pressure(const unsigned &n_p, const double &pvalue)
Fix the pressure at local pressure node n_p.
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...
virtual unsigned u_index_axi_nst(const unsigned &i) const
Return the index at which the i-th unknown velocity component.
double interpolated_dudx_axi_nst(const Vector< double > &s, const unsigned &i, const unsigned &j) const
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as ...
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
const Vector< double > & g() const
Vector of gravitational components.
virtual unsigned npres_axi_nst() const =0
Function to return number of pressure degrees of freedom.
double interpolated_u_axi_nst(const Vector< double > &s, const unsigned &i) const
Return FE interpolated velocity u[i] at local coordinate s.
void enable_ALE()
(Re-)enable ALE, i.e. take possible mesh motion into account when evaluating the time-derivative. Note: By default, ALE is enabled, at the expense of possibly creating unnecessary work in problems where the mesh is, in fact, stationary.
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 scalar_value_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
Write values of the i-th scalar field at the plot points. Needs to be implemented for each new specif...
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
virtual void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned > > &dof_lookup_list) const
Create a list of pairs for the unknowns that this element is "in charge of" – ignore any unknowns as...
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. Plot at a given number of plot points and compute L2 error and L2 norm of velocity solution over element.
int p_local_eqn(const unsigned &n) const
Overload the access function for the pressure's local equation numbers.
void output(std::ostream &outfile, const unsigned &n_plot)
Redirect output to NavierStokesEquations output.
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
double * Re_pt
Pointer to global Reynolds number.
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...
Data *& internal_data_pt(const unsigned &i)
Return a pointer to i-th internal data object.
unsigned nfields_for_projection()
Number of fields to be projected: dim+1, corresponding to velocity components and pressure...
virtual void get_dresidual_dnodal_coordinates(RankThreeTensor< double > &dresidual_dnodal_coordinates)
Compute derivatives of elemental residual vector with respect to nodal coordinates. This function computes these terms analytically and overwrites the default implementation in the FiniteElement base class. dresidual_dnodal_coordinates(l,i,j) = d res(l) / dX_{ij}.
void pshape_axi_nst(const Vector< double > &s, Shape &psi) const
Pressure shape functions at local coordinate s.
A class that represents a collection of data; each Data object may contain many different individual ...
virtual void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Function to compute the geometric shape functions and derivatives w.r.t. local coordinates at local c...
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
double p_axi_nst(const unsigned &n_p) const
Access function for the pressure values at local pressure node n_p (const version) ...
double(*&)(const double &time, const Vector< double > &x) source_fct_pt()
Access function for the source-function pointer.
double du_dt_axi_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...
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
virtual int p_nodal_index_axi_nst() const
Which nodal value represents the pressure?
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...
unsigned ntstorage() const
Return the number of doubles required to represent history (one for steady)
void output(FILE *file_pt)
Redirect output to NavierStokesEquations output.
void interpolated_u_axi_nst(const Vector< double > &s, Vector< double > &veloc) const
Compute vector of FE interpolated velocity u at local coordinate s.
unsigned ntstorage() const
Return total number of doubles stored per value to record time history of each value (one for steady ...
static Vector< double > Default_Gravity_vector
Static default value for the gravity vector.
double *& density_ratio_pt()
Pointer to Density ratio.
double *& viscosity_ratio_pt()
Pointer to Viscosity Ratio.
double *& re_pt()
Pointer to Reynolds number.
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
void output(std::ostream &outfile)
Output function: x,y,[z],u,v,[w],p in tecplot format. Default number of plot points.
double compute_physical_size() const
Compute the volume of the element.
double * ReInvFr_pt
Pointer to global Reynolds number x inverse Froude number (= Bond number / Capillary number) ...
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...
Axisymmetric Taylor Hood upgraded to become projectable.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
void fix_pressure(const unsigned &p_dof, const double &pvalue)
Function to fix the internal pressure dof idof_internal.
double dshape_and_dtest_eulerian_at_knot_axi_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...
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Return interpolated field fld at local coordinate s, at time level t (t=0: present; t>0: history valu...
double dissipation() const
Return integral of dissipation over element.
int p_local_eqn(const unsigned &n) const
Overload the access function for the pressure's local.
double interpolated_p_axi_nst(const Vector< double > &s) const
Return FE interpolated pressure at local coordinate s.
void traction(const Vector< double > &s, const Vector< double > &N, Vector< double > &traction) const
Compute traction (on the viscous scale) at local coordinate s for outer unit normal N...
double dshape_and_dtest_eulerian_axi_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 point_output_data(const Vector< double > &s, Vector< double > &data)
Output solution in data vector at local cordinates s: r,z,u_r,u_z,u_phi,p.
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 ...
Vector< double > * G_pt
Pointer to global gravity Vector.
double * ReSt_pt
Pointer to global Reynolds number x Strouhal number (=Womersley)
bool is_steady() const
Flag to indicate if a timestepper has been made steady (possibly temporarily to switch off time-depen...
void output(FILE *file_pt, const unsigned &n_plot)
Redirect output to NavierStokesEquations output.
const double & re_invfr() const
Global inverse Froude number.
void output(FILE *file_pt, const unsigned &n_plot)
Redirect output to NavierStokesEquations output.
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
void output(FILE *file_pt)
Redirect output to NavierStokesEquations output.
unsigned npres_axi_nst() const
Return number of pressure values.
void fill_in_contribution_to_djacobian_and_dmass_matrix_dparameter(double *const ¶meter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam)
virtual int p_nodal_index_axi_nst() const
Which nodal value represents the pressure?
virtual double p_axi_nst(const unsigned &n_p) const =0
Pressure at local pressure "node" n_p Uses suitably interpolated value for hanging nodes...
virtual void d_dshape_eulerian_dnodal_coordinates(const double &det_jacobian, const DenseMatrix< double > &jacobian, const DenseMatrix< double > &djacobian_dX, const DenseMatrix< double > &inverse_jacobian, const DShape &dpsids, RankFourTensor< double > &d_dpsidx_dX) const
A template-free interface that calculates the derivative w.r.t. the nodal coordinates of the derivat...
unsigned nnode() const
Return the number of nodes.
Base class for time-stepping schemes. Timestepper provides an approximation of the temporal derivativ...
void fill_in_contribution_to_hessian_vector_products(Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
Compute the hessian tensor vector products required to perform continuation of bifurcations analytica...
virtual double weight(const unsigned &i, const unsigned &j) const
Access function for j-th weight for the i-th derivative.
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...
static double Default_fd_jacobian_step
Double used for the default finite difference step in elemental jacobian calculations.
double dshape_and_dtest_eulerian_axi_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...
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Return Jacobian of mapping and shape functions of field fld at local coordinate s.
virtual double J_eulerian(const Vector< double > &s) const
Return the Jacobian of mapping from local to global coordinates at local position s...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute the element's residual Vector and the jacobian matrix Virtual function can be overloaded by h...
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...
double * Density_Ratio_pt
Pointer to the density ratio (relative to the density used in the definition of the Reynolds number) ...
double interpolated_dudt_axi_nst(const Vector< double > &s, const unsigned &i) const
AxisymmetricQTaylorHoodElement()
Constructor, no internal data points.
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Return interpolated field fld at local coordinate s, at time level t (t=0: present; t>0: history valu...
double dshape_and_dtest_eulerian_at_knot_axi_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...
std::string to_string(T object, unsigned float_precision=8)
Conversion function that should work for anything with operator<< defined (at least all basic types)...
void output(FILE *file_pt)
Output function: x,y,[z],u,v,[w],p in tecplot format. Default number of plot points.
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...
AxisymmetricNavierStokesEquations()
Constructor: NULL the body force and source function.
static double Default_Physical_Constant_Value
Static default value for the physical constants (all initialised to zero)
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...
void shape< 2 >(const double &s, double *Psi)
1D shape functions specialised to linear order (2 Nodes)
unsigned u_index_nst(const unsigned &i) const
Return the index at which the i-th unknown velocity component is stored with a common interface for u...