28 #ifndef OOMPH_AXISYMMETRIC_NAVIER_STOKES_ELEMENTS_HEADER
29 #define OOMPH_AXISYMMETRIC_NAVIER_STOKES_ELEMENTS_HEADER
33 #include <oomph-lib-config.h>
37 #include "../generic/Qelements.h"
38 #include "../generic/fsi.h"
39 #include "../generic/projection.h"
187 DShape& dtestdx)
const = 0;
197 DShape& dtestdx)
const = 0;
232 for (
unsigned i = 0;
i < 3;
i++)
240 (*Body_force_fct_pt)(time,
x, result);
266 for (
unsigned i = 0;
i < 2;
i++)
270 for (
unsigned j = 0;
j < 3;
j++)
272 d_body_force_dx(
j,
i) =
319 double source_pls = 0.0;
321 for (
unsigned i = 0;
i < 2;
i++)
325 gradient[
i] = (source_pls -
source) / eps_fd;
360 double*
const& parameter_pt,
398 const double&
re()
const
548 for (
unsigned t = 0;
t < n_time;
t++)
576 virtual double p_axi_nst(
const unsigned& n_p)
const = 0;
615 const unsigned& which_one = 0);
628 const unsigned& nplot)
const
635 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
653 std::stringstream error_stream;
655 <<
"Axisymmetric Navier-Stokes Elements only store 4 fields "
682 std::stringstream error_stream;
684 <<
"Axisymmetric Navier-Stokes Elements only store 4 fields "
697 for (
unsigned i = 0;
i < 2;
i++)
703 for (
unsigned i = 0;
i < 3;
i++)
723 void output(std::ostream& outfile,
const unsigned& nplot);
736 void output(FILE* file_pt,
const unsigned& nplot);
742 const unsigned& nplot,
749 const unsigned& nplot,
756 const unsigned& nplot,
810 residuals, jacobian, mass_matrix, 2);
837 double*
const& parameter_pt,
853 double*
const& parameter_pt,
860 parameter_pt, dres_dparam, djac_dparam, dmass_matrix_dparam, 2);
869 unsigned n_node =
nnode();
875 for (
unsigned i = 0;
i < 3;
i++)
882 for (
unsigned l = 0; l < n_node; l++)
891 const unsigned&
i)
const
894 unsigned n_node =
nnode();
904 double interpolated_u = 0.0;
906 for (
unsigned l = 0; l < n_node; l++)
908 interpolated_u +=
nodal_value(l, u_nodal_index) * psi[l];
911 return (interpolated_u);
919 const unsigned&
i)
const
922 unsigned n_node =
nnode();
932 double interpolated_u = 0.0;
934 for (
unsigned l = 0; l < n_node; l++)
936 interpolated_u +=
nodal_value(
t, l, u_nodal_index) * psi[l];
939 return (interpolated_u);
955 unsigned n_node =
nnode();
965 unsigned n_u_dof = 0;
966 for (
unsigned l = 0; l < n_node; l++)
977 du_ddata.resize(n_u_dof, 0.0);
978 global_eqn_number.resize(n_u_dof, 0);
983 for (
unsigned l = 0; l < n_node; l++)
990 global_eqn_number[count] = global_eqn;
992 du_ddata[count] = psi[l];
1011 double interpolated_p = 0.0;
1013 for (
unsigned l = 0; l < n_pres; l++)
1015 interpolated_p +=
p_axi_nst(l) * psi[l];
1018 return (interpolated_p);
1025 const unsigned&
j)
const
1028 const unsigned n_node =
nnode();
1033 DShape dpsifds(n_node, 2);
1042 double interpolated_duds = 0.0;
1045 for (
unsigned l = 0; l < n_node; l++)
1047 interpolated_duds +=
nodal_value(l, u_nodal_index) * dpsifds(l,
j);
1050 return (interpolated_duds);
1057 const unsigned&
j)
const
1060 const unsigned n_node =
nnode();
1065 DShape dpsifdx(n_node, 2);
1074 double interpolated_dudx = 0.0;
1077 for (
unsigned l = 0; l < n_node; l++)
1079 interpolated_dudx +=
nodal_value(l, u_nodal_index) * dpsifdx(l,
j);
1082 return (interpolated_dudx);
1088 const unsigned&
i)
const
1091 const unsigned n_node =
nnode();
1100 double interpolated_dudt = 0.0;
1103 for (
unsigned l = 0; l < n_node; l++)
1108 return (interpolated_dudt);
1118 const unsigned&
k)
const
1121 const unsigned n_node =
nnode();
1126 DShape dpsifds(n_node, 2);
1153 det, jacobian, djacobian_dX, inverse_jacobian, dpsifds, d_dpsifdx_dX);
1159 double interpolated_d_dudx_dX = 0.0;
1162 for (
unsigned l = 0; l < n_node; l++)
1164 interpolated_d_dudx_dX +=
1168 return (interpolated_d_dudx_dX);
1175 double result = 0.0;
1188 const unsigned n_dim =
dim();
1192 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
1195 for (
unsigned i = 0;
i < n_dim;
i++)
1201 for (
unsigned i = 0;
i < n_dim_eulerian;
i++)
1213 result +=
x[0] *
w *
J;
1257 const unsigned& ipt,
1267 const unsigned& ipt,
1339 void output(std::ostream& outfile,
const unsigned& n_plot)
1352 void output(FILE* file_pt,
const unsigned& n_plot)
1371 std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list)
const;
1392 for (
unsigned i = 0;
i < 9;
i++)
1395 dtestdx(
i, 0) = dpsidx(
i, 0);
1396 dtestdx(
i, 1) = dpsidx(
i, 1);
1418 for (
unsigned i = 0;
i < 9;
i++)
1421 dtestdx(
i, 0) = dpsidx(
i, 0);
1422 dtestdx(
i, 1) = dpsidx(
i, 1);
1438 const unsigned& ipt,
1449 ipt, psi, dpsidx, djacobian_dX, d_dpsidx_dX);
1453 for (
unsigned i = 0;
i < 9;
i++)
1457 for (
unsigned k = 0;
k < 2;
k++)
1459 dtestdx(
i,
k) = dpsidx(
i,
k);
1461 for (
unsigned p = 0;
p < 2;
p++)
1463 for (
unsigned q = 0;
q < 9;
q++)
1465 d_dtestdx_dX(
p,
q,
i,
k) = d_dpsidx_dX(
p,
q,
i,
k);
1493 for (
unsigned i = 0;
i < 3;
i++)
test[
i] = psi[
i];
1555 const unsigned& ipt,
1565 const unsigned& ipt,
1641 void output(std::ostream& outfile,
const unsigned& n_plot)
1653 void output(FILE* file_pt,
const unsigned& n_plot)
1673 std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list)
const;
1694 for (
unsigned i = 0;
i < 9;
i++)
1697 dtestdx(
i, 0) = dpsidx(
i, 0);
1698 dtestdx(
i, 1) = dpsidx(
i, 1);
1720 for (
unsigned i = 0;
i < 9;
i++)
1723 dtestdx(
i, 0) = dpsidx(
i, 0);
1724 dtestdx(
i, 1) = dpsidx(
i, 1);
1740 const unsigned& ipt,
1751 ipt, psi, dpsidx, djacobian_dX, d_dpsidx_dX);
1755 for (
unsigned i = 0;
i < 9;
i++)
1759 for (
unsigned k = 0;
k < 2;
k++)
1761 dtestdx(
i,
k) = dpsidx(
i,
k);
1763 for (
unsigned p = 0;
p < 2;
p++)
1765 for (
unsigned q = 0;
q < 9;
q++)
1767 d_dtestdx_dX(
p,
q,
i,
k) = d_dpsidx_dX(
p,
q,
i,
k);
1784 double psi1[2], psi2[2];
1791 for (
unsigned i = 0;
i < 2;
i++)
1793 for (
unsigned j = 0;
j < 2;
j++)
1796 psi[2 *
i +
j] = psi2[
i] * psi1[
j];
1810 for (
unsigned i = 0;
i < 4;
i++)
test[
i] = psi[
i];
1840 template<
class TAYLOR_HOOD_ELEMENT>
1860 unsigned nnod = this->
nnode();
1861 for (
unsigned j = 0;
j < nnod;
j++)
1864 data_values.push_back(std::make_pair(this->
node_pt(
j), fld));
1871 unsigned Pconv_size = 3;
1872 for (
unsigned j = 0;
j < Pconv_size;
j++)
1875 unsigned vertex_index = this->Pconv[
j];
1876 data_values.push_back(
1877 std::make_pair(this->
node_pt(vertex_index), fld));
1920 unsigned n_dim = this->
dim();
1921 unsigned n_node = this->
nnode();
1926 this->pshape_axi_nst(
s, psi);
1928 Shape psif(n_node), testf(n_node);
1929 DShape dpsifdx(n_node, n_dim), dtestfdx(n_node, n_dim);
1932 double J = this->dshape_and_dtest_eulerian_axi_nst(
1933 s, psif, dpsifdx, testf, dtestfdx);
1938 Shape testf(n_node);
1939 DShape dpsifdx(n_node, n_dim), dtestfdx(n_node, n_dim);
1942 double J = this->dshape_and_dtest_eulerian_axi_nst(
1943 s, psi, dpsifdx, testf, dtestfdx);
1952 const unsigned& fld,
1955 unsigned n_node = this->
nnode();
1960 return this->interpolated_p_axi_nst(
s);
1966 unsigned u_nodal_index = this->u_index_axi_nst(fld);
1972 this->
shape(s, psi);
1975 double interpolated_u = 0.0;
1978 for (
unsigned l = 0; l < n_node; l++)
1980 interpolated_u += this->
nodal_value(t, l, u_nodal_index) * psi[l];
1982 return interpolated_u;
1992 return this->npres_axi_nst();
1996 return this->
nnode();
2006 return this->p_local_eqn(
j);
2010 const unsigned u_nodal_index = this->u_index_axi_nst(fld);
2021 template<
class ELEMENT>
2034 template<
class ELEMENT>
2047 template<
class CROUZEIX_RAVIART_ELEMENT>
2067 const unsigned n_node = this->
nnode();
2068 for (
unsigned n = 0;
n < n_node;
n++)
2071 data_values.push_back(std::make_pair(this->
node_pt(
n), fld));
2078 const unsigned n_press = this->npres_axi_nst();
2080 for (
unsigned j = 0;
j < n_press;
j++)
2082 data_values.push_back(std::make_pair(
2126 unsigned n_dim = this->
dim();
2127 unsigned n_node = this->
nnode();
2132 this->pshape_axi_nst(
s, psi);
2134 Shape psif(n_node), testf(n_node);
2135 DShape dpsifdx(n_node, n_dim), dtestfdx(n_node, n_dim);
2138 double J = this->dshape_and_dtest_eulerian_axi_nst(
2139 s, psif, dpsifdx, testf, dtestfdx);
2144 Shape testf(n_node);
2145 DShape dpsifdx(n_node, n_dim), dtestfdx(n_node, n_dim);
2148 double J = this->dshape_and_dtest_eulerian_axi_nst(
2149 s, psi, dpsifdx, testf, dtestfdx);
2158 const unsigned& fld,
2167 return this->interpolated_p_axi_nst(
s);
2172 return this->interpolated_u_axi_nst(
t,
s, fld);
2182 return this->npres_axi_nst();
2186 return this->
nnode();
2196 return this->p_local_eqn(
j);
2200 const unsigned u_nodal_index = this->u_index_axi_nst(fld);
2225 std::set<std::pair<Data*, unsigned>>& paired_load_data)
2231 unsigned u_index[
DIM];
2232 for (
unsigned i = 0;
i <
DIM;
i++)
2238 unsigned n_node = this->
nnode();
2239 for (
unsigned n = 0;
n < n_node;
n++)
2243 for (
unsigned i = 0;
i <
DIM;
i++)
2245 paired_load_data.insert(std::make_pair(this->
node_pt(
n), u_index[
i]));
2262 std::set<std::pair<Data*, unsigned>>& paired_pressure_data)
2269 for (
unsigned l = 0; l < n_pres; l++)
2273 paired_pressure_data.insert(
2319 template<
class ELEMENT>
2332 template<
class ELEMENT>
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
int data[]
Definition: Map_placement_new.cpp:1
RowVector3d w
Definition: Matrix_resize_int.cpp:3
void load(Archive &ar, ParticleHandler &handl)
Definition: Particles.h:21
float * p
Definition: Tutorial_Map_using.cpp:9
Definition: axisym_navier_stokes_elements.h:115
double dissipation() const
Return integral of dissipation over element.
Definition: axisym_navier_stokes_elements.cc:594
void point_output_data(const Vector< double > &s, Vector< double > &data)
Definition: axisym_navier_stokes_elements.h:694
static Vector< double > Default_Gravity_vector
Static default value for the gravity vector.
Definition: axisym_navier_stokes_elements.h:130
double pressure_integral() const
Integral of pressure over element.
Definition: axisym_navier_stokes_elements.cc:869
void traction(const Vector< double > &s, const Vector< double > &N, Vector< double > &traction) const
Definition: axisym_navier_stokes_elements.cc:649
unsigned nscalar_paraview() const
Definition: axisym_navier_stokes_elements.h:619
unsigned u_index_nst(const unsigned &i) const
Definition: axisym_navier_stokes_elements.h:516
double(* Source_fct_pt)(const double &time, const Vector< double > &x)
Pointer to volumetric source function.
Definition: axisym_navier_stokes_elements.h:168
virtual void pshape_axi_nst(const Vector< double > &s, Shape &psi, Shape &test) const =0
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)
Definition: axisym_navier_stokes_elements.cc:2510
Vector< double > * G_pt
Pointer to global gravity Vector.
Definition: axisym_navier_stokes_elements.h:160
double interpolated_duds_axi_nst(const Vector< double > &s, const unsigned &i, const unsigned &j) const
Definition: axisym_navier_stokes_elements.h:1023
double interpolated_dudx_axi_nst(const Vector< double > &s, const unsigned &i, const unsigned &j) const
Definition: axisym_navier_stokes_elements.h:1055
const double & density_ratio() const
Definition: axisym_navier_stokes_elements.h:459
void fill_in_contribution_to_hessian_vector_products(Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
Definition: axisym_navier_stokes_elements.cc:3239
void scalar_value_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
Definition: axisym_navier_stokes_elements.h:626
virtual double dshape_and_dtest_eulerian_at_knot_axi_nst(const unsigned &ipt, Shape &psi, DShape &dpsidx, RankFourTensor< double > &d_dpsidx_dX, Shape &test, DShape &dtestdx, RankFourTensor< double > &d_dtestdx_dX, DenseMatrix< double > &djacobian_dX) const =0
double * ReInvFr_pt
Definition: axisym_navier_stokes_elements.h:153
const double & re_invro() const
Global Reynolds number multiplied by inverse Rossby number.
Definition: axisym_navier_stokes_elements.h:434
double compute_physical_size() const
Compute the volume of the element.
Definition: axisym_navier_stokes_elements.h:1172
void output_veloc(std::ostream &outfile, const unsigned &nplot, const unsigned &t)
Definition: axisym_navier_stokes_elements.cc:422
virtual double p_axi_nst(const unsigned &n_p) const =0
double interpolated_u_axi_nst(const Vector< double > &s, const unsigned &i) const
Return FE interpolated velocity u[i] at local coordinate s.
Definition: axisym_navier_stokes_elements.h:890
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Compute the element's residual Vector.
Definition: axisym_navier_stokes_elements.h:780
double *& re_invro_pt()
Pointer to global inverse Froude number.
Definition: axisym_navier_stokes_elements.h:440
std::string scalar_name_paraview(const unsigned &i) const
Definition: axisym_navier_stokes_elements.h:667
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: axisym_navier_stokes_elements.cc:310
virtual int p_nodal_index_axi_nst() const
Which nodal value represents the pressure?
Definition: axisym_navier_stokes_elements.h:579
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.
Definition: axisym_navier_stokes_elements.h:917
const double & viscosity_ratio() const
Definition: axisym_navier_stokes_elements.h:472
const Vector< double > & g() const
Vector of gravitational components.
Definition: axisym_navier_stokes_elements.h:446
const double & re() const
Reynolds number.
Definition: axisym_navier_stokes_elements.h:398
double *& density_ratio_pt()
Pointer to Density ratio.
Definition: axisym_navier_stokes_elements.h:465
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: axisym_navier_stokes_elements.h:793
void output(FILE *file_pt)
Definition: axisym_navier_stokes_elements.h:728
double * Re_pt
Pointer to global Reynolds number.
Definition: axisym_navier_stokes_elements.h:146
virtual void pshape_axi_nst(const Vector< double > &s, Shape &psi) const =0
Compute the pressure shape functions at local coordinate s.
AxisymmetricNavierStokesEquations()
Constructor: NULL the body force and source function.
Definition: axisym_navier_stokes_elements.h:375
double *& re_pt()
Pointer to Reynolds number.
Definition: axisym_navier_stokes_elements.h:410
virtual void get_viscosity_ratio_axisym_nst(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, double &visc_ratio)
Definition: axisym_navier_stokes_elements.h:339
double interpolated_d_dudx_dX_axi_nst(const Vector< double > &s, const unsigned &p, const unsigned &q, const unsigned &i, const unsigned &k) const
Definition: axisym_navier_stokes_elements.h:1114
void interpolated_u_axi_nst(const Vector< double > &s, Vector< double > &veloc) const
Compute vector of FE interpolated velocity u at local coordinate s.
Definition: axisym_navier_stokes_elements.h:865
double *& re_st_pt()
Pointer to product of Reynolds and Strouhal number (=Womersley number)
Definition: axisym_navier_stokes_elements.h:416
double * ReSt_pt
Pointer to global Reynolds number x Strouhal number (=Womersley)
Definition: axisym_navier_stokes_elements.h:149
void disable_ALE()
Definition: axisym_navier_stokes_elements.h:560
virtual double dshape_and_dtest_eulerian_at_knot_axi_nst(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
void(*&)(const double &time, const Vector< double > &x, Vector< double > &f) axi_nst_body_force_fct_pt()
Access function for the body-force pointer.
Definition: axisym_navier_stokes_elements.h:484
double get_source_fct(const double &time, const unsigned &ipt, const Vector< double > &x)
Calculate the source fct at given time and Eulerian position.
Definition: axisym_navier_stokes_elements.h:286
double * Density_Ratio_pt
Definition: axisym_navier_stokes_elements.h:141
static int Pressure_not_stored_at_node
Definition: axisym_navier_stokes_elements.h:119
double du_dt_axi_nst(const unsigned &n, const unsigned &i) const
Definition: axisym_navier_stokes_elements.h:531
static double Default_Physical_Ratio_Value
Definition: axisym_navier_stokes_elements.h:127
double kin_energy() const
Get integral of kinetic energy over element.
Definition: axisym_navier_stokes_elements.cc:821
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Definition: axisym_navier_stokes_elements.cc:155
double *& re_invfr_pt()
Pointer to global inverse Froude number.
Definition: axisym_navier_stokes_elements.h:428
void(* Body_force_fct_pt)(const double &time, const Vector< double > &x, Vector< double > &result)
Pointer to body force function.
Definition: axisym_navier_stokes_elements.h:163
virtual unsigned npres_axi_nst() const =0
Function to return number of pressure degrees of freedom.
double *& viscosity_ratio_pt()
Pointer to Viscosity Ratio.
Definition: axisym_navier_stokes_elements.h:478
virtual double dshape_and_dtest_eulerian_axi_nst(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
static double Default_Physical_Constant_Value
Navier–Stokes equations static data.
Definition: axisym_navier_stokes_elements.h:123
double interpolated_dudt_axi_nst(const Vector< double > &s, const unsigned &i) const
Definition: axisym_navier_stokes_elements.h:1087
double(*&)(const double &time, const Vector< double > &x) source_fct_pt()
Access function for the source-function pointer.
Definition: axisym_navier_stokes_elements.h:492
bool ALE_is_disabled
Definition: axisym_navier_stokes_elements.h:173
virtual void fill_in_generic_residual_contribution_axi_nst(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: axisym_navier_stokes_elements.cc:914
virtual void get_source_fct_gradient(const double &time, const unsigned &ipt, const Vector< double > &x, Vector< double > &gradient)
Definition: axisym_navier_stokes_elements.h:305
unsigned n_u_nst() const
Definition: axisym_navier_stokes_elements.h:523
void fill_in_contribution_to_dresiduals_dparameter(double *const ¶meter_pt, Vector< double > &dres_dparam)
Compute the element's residual Vector.
Definition: axisym_navier_stokes_elements.h:821
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)
Definition: axisym_navier_stokes_elements.h:852
const double & re_st() const
Product of Reynolds and Strouhal number (=Womersley number)
Definition: axisym_navier_stokes_elements.h:404
const double & re_invfr() const
Global inverse Froude number.
Definition: axisym_navier_stokes_elements.h:422
double interpolated_p_axi_nst(const Vector< double > &s) const
Return FE interpolated pressure at local coordinate s.
Definition: axisym_navier_stokes_elements.h:1001
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)
Definition: axisym_navier_stokes_elements.h:247
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
Definition: axisym_navier_stokes_elements.h:803
double * ReInvRo_pt
Definition: axisym_navier_stokes_elements.h:157
void get_pressure_and_velocity_mass_matrix_diagonal(Vector< double > &press_mass_diag, Vector< double > &veloc_mass_diag, const unsigned &which_one=0)
Definition: axisym_navier_stokes_elements.cc:63
virtual void dinterpolated_u_axi_nst_ddata(const Vector< double > &s, const unsigned &i, Vector< double > &du_ddata, Vector< unsigned > &global_eqn_number)
Definition: axisym_navier_stokes_elements.h:948
void enable_ALE()
Definition: axisym_navier_stokes_elements.h:569
void strain_rate(const Vector< double > &s, DenseMatrix< double > &strain_rate) const
Strain-rate tensor: where (in that order)
Definition: axisym_navier_stokes_elements.cc:725
virtual int p_local_eqn(const unsigned &n) const =0
void fill_in_contribution_to_djacobian_dparameter(double *const ¶meter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
Definition: axisym_navier_stokes_elements.h:836
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.
Definition: axisym_navier_stokes_elements.h:222
double * Viscosity_Ratio_pt
Definition: axisym_navier_stokes_elements.h:137
virtual void get_dresidual_dnodal_coordinates(RankThreeTensor< double > &dresidual_dnodal_coordinates)
Definition: axisym_navier_stokes_elements.cc:1567
Vector< double > *& g_pt()
Pointer to Vector of gravitational components.
Definition: axisym_navier_stokes_elements.h:452
virtual unsigned u_index_axi_nst(const unsigned &i) const
Definition: axisym_navier_stokes_elements.h:506
void output(std::ostream &outfile)
Definition: axisym_navier_stokes_elements.h:715
static Vector< double > Gamma
Vector to decide whether the stress-divergence form is used or not.
Definition: axisym_navier_stokes_elements.h:393
Definition: axisym_navier_stokes_elements.h:1234
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Definition: axisym_navier_stokes_elements.cc:3639
void output(FILE *file_pt, const unsigned &n_plot)
Redirect output to NavierStokesEquations output.
Definition: axisym_navier_stokes_elements.h:1352
void get_traction(const Vector< double > &s, const Vector< double > &N, Vector< double > &traction) const
Compute traction at local coordinate s for outer unit normal N.
Definition: axisym_navier_stokes_elements.cc:3720
virtual unsigned required_nvalue(const unsigned &n) const
Number of values (pinned or dofs) required at local node n.
Definition: axisym_navier_stokes_elements.cc:3711
void fix_pressure(const unsigned &p_dof, const double &pvalue)
Function to fix the internal pressure dof idof_internal.
Definition: axisym_navier_stokes_elements.h:1314
void pshape_axi_nst(const Vector< double > &s, Shape &psi) const
Pressure shape functions at local coordinate s.
Definition: axisym_navier_stokes_elements.h:1478
unsigned ndof_types() const
Definition: axisym_navier_stokes_elements.h:1359
void output(FILE *file_pt)
Redirect output to NavierStokesEquations output.
Definition: axisym_navier_stokes_elements.h:1346
unsigned P_axi_nst_internal_index
Definition: axisym_navier_stokes_elements.h:1242
unsigned npres_axi_nst() const
Return number of pressure values.
Definition: axisym_navier_stokes_elements.h:1308
AxisymmetricQCrouzeixRaviartElement()
Constructor, there are three internal values (for the pressure)
Definition: axisym_navier_stokes_elements.h:1288
static const unsigned Initial_Nvalue[]
Static array of ints to hold required number of variables at nodes.
Definition: axisym_navier_stokes_elements.h:1237
double p_axi_nst(const unsigned &i) const
Definition: axisym_navier_stokes_elements.h:1302
double dshape_and_dtest_eulerian_axi_nst(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: axisym_navier_stokes_elements.h:1382
int p_local_eqn(const unsigned &n) const
Definition: axisym_navier_stokes_elements.h:1327
double dshape_and_dtest_eulerian_at_knot_axi_nst(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: axisym_navier_stokes_elements.h:1408
void output(std::ostream &outfile)
Redirect output to NavierStokesEquations output.
Definition: axisym_navier_stokes_elements.h:1333
void output(std::ostream &outfile, const unsigned &n_plot)
Redirect output to NavierStokesEquations output.
Definition: axisym_navier_stokes_elements.h:1339
Definition: axisym_navier_stokes_elements.h:1532
void fix_pressure(const unsigned &n_p, const double &pvalue)
Fix the pressure at local pressure node n_p.
Definition: axisym_navier_stokes_elements.h:1616
virtual int p_nodal_index_axi_nst() const
Which nodal value represents the pressure?
Definition: axisym_navier_stokes_elements.h:1597
static const unsigned Pconv[]
Definition: axisym_navier_stokes_elements.h:1540
double dshape_and_dtest_eulerian_axi_nst(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: axisym_navier_stokes_elements.h:1684
void output(FILE *file_pt, const unsigned &n_plot)
Redirect output to NavierStokesEquations output.
Definition: axisym_navier_stokes_elements.h:1653
int p_local_eqn(const unsigned &n) const
Definition: axisym_navier_stokes_elements.h:1629
unsigned npres_axi_nst() const
Return number of pressure values.
Definition: axisym_navier_stokes_elements.h:1610
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Definition: axisym_navier_stokes_elements.cc:3741
void output(std::ostream &outfile, const unsigned &n_plot)
Redirect output to NavierStokesEquations output.
Definition: axisym_navier_stokes_elements.h:1641
virtual unsigned required_nvalue(const unsigned &n) const
Definition: axisym_navier_stokes_elements.h:1591
void output(FILE *file_pt)
Redirect output to NavierStokesEquations output.
Definition: axisym_navier_stokes_elements.h:1647
unsigned ndof_types() const
Definition: axisym_navier_stokes_elements.h:1661
void output(std::ostream &outfile)
Redirect output to NavierStokesEquations output.
Definition: axisym_navier_stokes_elements.h:1635
void pshape_axi_nst(const Vector< double > &s, Shape &psi) const
Pressure shape functions at local coordinate s.
Definition: axisym_navier_stokes_elements.h:1780
double dshape_and_dtest_eulerian_at_knot_axi_nst(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: axisym_navier_stokes_elements.h:1710
void get_traction(const Vector< double > &s, const Vector< double > &N, Vector< double > &traction) const
Compute traction at local coordinate s for outer unit normal N.
Definition: axisym_navier_stokes_elements.cc:3794
AxisymmetricQTaylorHoodElement()
Constructor, no internal data points.
Definition: axisym_navier_stokes_elements.h:1584
static const unsigned Initial_Nvalue[]
Static array of ints to hold number of variables at node.
Definition: axisym_navier_stokes_elements.h:1535
double p_axi_nst(const unsigned &n_p) const
Definition: axisym_navier_stokes_elements.h:1604
long & eqn_number(const unsigned &i)
Return the equation number of the i-th stored variable.
Definition: nodes.h:367
void pin(const unsigned &i)
Pin the i-th stored variable.
Definition: nodes.h:385
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
void set_value(const unsigned &i, const double &value_)
Definition: nodes.h:271
double value(const unsigned &i) const
Definition: nodes.h:293
unsigned ntstorage() const
Definition: nodes.cc:879
Axisymmetric FSI Element.
Definition: axisym_navier_stokes_elements.h:2212
FSIAxisymmetricQTaylorHoodElement()
Constructor.
Definition: axisym_navier_stokes_elements.h:2215
void identify_load_data(std::set< std::pair< Data *, unsigned >> &paired_load_data)
Definition: axisym_navier_stokes_elements.h:2224
void get_load(const Vector< double > &s, const Vector< double > &N, Vector< double > &load)
Definition: axisym_navier_stokes_elements.h:2282
void identify_pressure_data(std::set< std::pair< Data *, unsigned >> &paired_pressure_data)
Definition: axisym_navier_stokes_elements.h:2261
FaceGeometry()
Definition: axisym_navier_stokes_elements.h:1505
FaceGeometry()
Definition: axisym_navier_stokes_elements.h:1821
FaceGeometry()
Definition: axisym_navier_stokes_elements.h:2299
FaceGeometry()
Definition: axisym_navier_stokes_elements.h:1517
FaceGeometry()
Definition: axisym_navier_stokes_elements.h:1833
FaceGeometry()
Definition: axisym_navier_stokes_elements.h:2311
FaceGeometry()
Definition: axisym_navier_stokes_elements.h:2338
FaceGeometry()
Definition: axisym_navier_stokes_elements.h:2040
FaceGeometry()
Definition: axisym_navier_stokes_elements.h:2324
FaceGeometry()
Definition: axisym_navier_stokes_elements.h:2026
Definition: elements.h:4998
Definition: elements.h:1313
virtual unsigned nplot_points_paraview(const unsigned &nplot) const
Definition: elements.h:2862
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
Definition: elements.cc:2749
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
Definition: elements.cc:3962
virtual void dJ_eulerian_dnodal_coordinates(const DenseMatrix< double > &jacobian, const DShape &dpsids, DenseMatrix< double > &djacobian_dX) const
Definition: elements.cc:2669
virtual void shape(const Vector< double > &s, Shape &psi) const =0
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432
unsigned dim() const
Definition: elements.h:2611
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Definition: elements.h:1759
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Definition: elements.h:3148
virtual double local_to_eulerian_mapping(const DShape &dpsids, DenseMatrix< double > &jacobian, DenseMatrix< double > &inverse_jacobian) const
Definition: elements.h:1508
virtual double dshape_eulerian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidx) const
Definition: elements.cc:3325
virtual void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Definition: elements.h:1981
unsigned nodal_dimension() const
Return the required Eulerian dimension of the nodes in this element.
Definition: elements.h:2484
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Definition: elements.h:1765
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Definition: elements.cc:3298
virtual double J_eulerian(const Vector< double > &s) const
Definition: elements.cc:4103
static double Default_fd_jacobian_step
Definition: elements.h:1198
Data *& internal_data_pt(const unsigned &i)
Return a pointer to i-th internal data object.
Definition: elements.h:622
int internal_local_eqn(const unsigned &i, const unsigned &j) const
Definition: elements.h:267
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227
unsigned add_internal_data(Data *const &data_pt, const bool &fd=true)
Definition: elements.cc:62
TimeStepper *& time_stepper_pt()
Definition: geom_objects.h:192
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
Definition: matrices.h:74
Definition: elements.h:5231
TimeStepper *& position_time_stepper_pt()
Return a pointer to the position timestepper.
Definition: nodes.h:1022
Definition: oomph_definitions.h:222
Definition: elements.h:3439
Crouzeix Raviart upgraded to become projectable.
Definition: axisym_navier_stokes_elements.h:2050
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld.
Definition: axisym_navier_stokes_elements.h:2178
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Definition: axisym_navier_stokes_elements.h:2157
unsigned nhistory_values_for_coordinate_projection()
Number of positional history values (includes current value!)
Definition: axisym_navier_stokes_elements.h:2115
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Definition: axisym_navier_stokes_elements.h:2058
unsigned nhistory_values_for_projection(const unsigned &fld)
Definition: axisym_navier_stokes_elements.h:2101
unsigned nfields_for_projection()
Definition: axisym_navier_stokes_elements.h:2093
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Definition: axisym_navier_stokes_elements.h:2122
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
Definition: axisym_navier_stokes_elements.h:2192
Axisymmetric Taylor Hood upgraded to become projectable.
Definition: axisym_navier_stokes_elements.h:1843
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Definition: axisym_navier_stokes_elements.h:1916
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld.
Definition: axisym_navier_stokes_elements.h:1988
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
Definition: axisym_navier_stokes_elements.h:2002
unsigned nfields_for_projection()
Definition: axisym_navier_stokes_elements.h:1887
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Definition: axisym_navier_stokes_elements.h:1951
unsigned nhistory_values_for_coordinate_projection()
Number of positional history values (includes current value!)
Definition: axisym_navier_stokes_elements.h:1909
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Definition: axisym_navier_stokes_elements.h:1851
unsigned nhistory_values_for_projection(const unsigned &fld)
Definition: axisym_navier_stokes_elements.h:1895
Definition: projection.h:183
Definition: Qelements.h:459
A Rank 4 Tensor class.
Definition: matrices.h:1701
A Rank 3 Tensor class.
Definition: matrices.h:1370
Definition: timesteppers.h:231
unsigned ntstorage() const
Definition: timesteppers.h:601
virtual double weight(const unsigned &i, const unsigned &j) const
Access function for j-th weight for the i-th derivative.
Definition: timesteppers.h:594
bool is_steady() const
Definition: timesteppers.h:389
@ N
Definition: constructor.cpp:22
static int f(const TensorMap< Tensor< int, 3 > > &tensor)
Definition: cxx11_tensor_map.cpp:237
RealScalar s
Definition: level1_cplx_impl.h:130
char char char int int * k
Definition: level2_impl.h:374
#define DIM
Definition: linearised_navier_stokes_elements.h:44
EIGEN_DEVICE_FUNC const Scalar & q
Definition: SpecialFunctionsImpl.h:2019
void body_force(const double &time, const Vector< double > &x, Vector< double > &result)
Definition: axisym_linear_elasticity/cylinder/cylinder.cc:96
void source(const Vector< double > &x, Vector< double > &f)
Source function.
Definition: unstructured_two_d_circle.cc:46
int error
Definition: calibrate.py:297
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286
const double Pi
50 digits from maple
Definition: oomph_utilities.h:157
void shape< 2 >(const double &s, double *Psi)
1D shape functions specialised to linear order (2 Nodes)
Definition: shape.h:608
std::string to_string(T object, unsigned float_precision=8)
Definition: oomph_utilities.h:189
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
list x
Definition: plotDoE.py:28
t
Definition: plotPSD.py:36
Definition: indexed_view.cpp:20
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
void product(const MatrixType &m)
Definition: product.h:42
void set(Container &c, Position position, const Value &value)
Definition: stdlist_overload.cpp:36
const char Y
Definition: test/EulerAngles.cpp:32
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2