28 #ifndef OOMPH_GENERALISED_NEWTONIAN_AXISYMMETRIC_NAVIER_STOKES_ELEMENTS_HEADER
29 #define OOMPH_GENERALISED_NEWTONIAN_AXISYMMETRIC_NAVIER_STOKES_ELEMENTS_HEADER
33 #include <oomph-lib-config.h>
39 #include "../generic/Qelements.h"
40 #include "../generic/fsi.h"
41 #include "../generic/projection.h"
42 #include "../generic/generalised_newtonian_constitutive_models.h"
202 DShape& dtestdx)
const = 0;
212 DShape& dtestdx)
const = 0;
247 for (
unsigned i = 0;
i < 3;
i++)
255 (*Body_force_fct_pt)(time,
x, result);
281 for (
unsigned i = 0;
i < 2;
i++)
285 for (
unsigned j = 0;
j < 3;
j++)
287 d_body_force_dx(
j,
i) =
334 double source_pls = 0.0;
336 for (
unsigned i = 0;
i < 2;
i++)
340 gradient[
i] = (source_pls -
source) / eps_fd;
378 double*
const& parameter_pt,
420 const double&
re()
const
587 for (
unsigned t = 0;
t < n_time;
t++)
615 virtual double p_axi_nst(
const unsigned& n_p)
const = 0;
630 double& max_invariant,
631 double& min_viscosity,
632 double& max_viscosity)
const;
689 const unsigned& which_one = 0);
702 const unsigned& nplot)
const
709 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
727 std::stringstream error_stream;
729 <<
"Axisymmetric Navier-Stokes Elements only store 4 fields "
756 std::stringstream error_stream;
758 <<
"Axisymmetric Navier-Stokes Elements only store 4 fields "
771 for (
unsigned i = 0;
i < 2;
i++)
777 for (
unsigned i = 0;
i < 3;
i++)
797 void output(std::ostream& outfile,
const unsigned& nplot);
810 void output(FILE* file_pt,
const unsigned& nplot);
816 const unsigned& nplot,
823 const unsigned& nplot,
830 const unsigned& nplot,
889 residuals, jacobian, mass_matrix, 2);
916 double*
const& parameter_pt,
932 double*
const& parameter_pt,
939 parameter_pt, dres_dparam, djac_dparam, dmass_matrix_dparam, 2);
948 unsigned n_node =
nnode();
954 for (
unsigned i = 0;
i < 3;
i++)
961 for (
unsigned l = 0; l < n_node; l++)
970 const unsigned&
i)
const
973 unsigned n_node =
nnode();
983 double interpolated_u = 0.0;
985 for (
unsigned l = 0; l < n_node; l++)
987 interpolated_u +=
nodal_value(l, u_nodal_index) * psi[l];
990 return (interpolated_u);
998 const unsigned&
i)
const
1001 unsigned n_node =
nnode();
1011 double interpolated_u = 0.0;
1013 for (
unsigned l = 0; l < n_node; l++)
1015 interpolated_u +=
nodal_value(
t, l, u_nodal_index) * psi[l];
1018 return (interpolated_u);
1034 unsigned n_node =
nnode();
1044 unsigned n_u_dof = 0;
1045 for (
unsigned l = 0; l < n_node; l++)
1049 if (global_eqn >= 0)
1056 du_ddata.resize(n_u_dof, 0.0);
1057 global_eqn_number.resize(n_u_dof, 0);
1062 for (
unsigned l = 0; l < n_node; l++)
1066 if (global_eqn >= 0)
1069 global_eqn_number[count] = global_eqn;
1071 du_ddata[count] = psi[l];
1090 double interpolated_p = 0.0;
1092 for (
unsigned l = 0; l < n_pres; l++)
1094 interpolated_p +=
p_axi_nst(l) * psi[l];
1097 return (interpolated_p);
1104 const unsigned&
j)
const
1107 const unsigned n_node =
nnode();
1112 DShape dpsifds(n_node, 2);
1121 double interpolated_duds = 0.0;
1124 for (
unsigned l = 0; l < n_node; l++)
1126 interpolated_duds +=
nodal_value(l, u_nodal_index) * dpsifds(l,
j);
1129 return (interpolated_duds);
1136 const unsigned&
j)
const
1139 const unsigned n_node =
nnode();
1144 DShape dpsifdx(n_node, 2);
1153 double interpolated_dudx = 0.0;
1156 for (
unsigned l = 0; l < n_node; l++)
1158 interpolated_dudx +=
nodal_value(l, u_nodal_index) * dpsifdx(l,
j);
1161 return (interpolated_dudx);
1167 const unsigned&
i)
const
1170 const unsigned n_node =
nnode();
1179 double interpolated_dudt = 0.0;
1182 for (
unsigned l = 0; l < n_node; l++)
1187 return (interpolated_dudt);
1197 const unsigned&
k)
const
1200 const unsigned n_node =
nnode();
1205 DShape dpsifds(n_node, 2);
1232 det, jacobian, djacobian_dX, inverse_jacobian, dpsifds, d_dpsifdx_dX);
1238 double interpolated_d_dudx_dX = 0.0;
1241 for (
unsigned l = 0; l < n_node; l++)
1243 interpolated_d_dudx_dX +=
1247 return (interpolated_d_dudx_dX);
1254 double result = 0.0;
1267 const unsigned n_dim =
dim();
1271 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
1274 for (
unsigned i = 0;
i < n_dim;
i++)
1280 for (
unsigned i = 0;
i < n_dim_eulerian;
i++)
1292 result +=
x[0] *
w *
J;
1336 const unsigned& ipt,
1346 const unsigned& ipt,
1419 void output(std::ostream& outfile,
const unsigned& n_plot)
1433 void output(FILE* file_pt,
const unsigned& n_plot)
1453 std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list)
const;
1474 for (
unsigned i = 0;
i < 9;
i++)
1477 dtestdx(
i, 0) = dpsidx(
i, 0);
1478 dtestdx(
i, 1) = dpsidx(
i, 1);
1500 for (
unsigned i = 0;
i < 9;
i++)
1503 dtestdx(
i, 0) = dpsidx(
i, 0);
1504 dtestdx(
i, 1) = dpsidx(
i, 1);
1520 const unsigned& ipt,
1531 ipt, psi, dpsidx, djacobian_dX, d_dpsidx_dX);
1535 for (
unsigned i = 0;
i < 9;
i++)
1539 for (
unsigned k = 0;
k < 2;
k++)
1541 dtestdx(
i,
k) = dpsidx(
i,
k);
1543 for (
unsigned p = 0;
p < 2;
p++)
1545 for (
unsigned q = 0;
q < 9;
q++)
1547 d_dtestdx_dX(
p,
q,
i,
k) = d_dpsidx_dX(
p,
q,
i,
k);
1575 for (
unsigned i = 0;
i < 3;
i++)
test[
i] = psi[
i];
1639 const unsigned& ipt,
1649 const unsigned& ipt,
1726 void output(std::ostream& outfile,
const unsigned& n_plot)
1739 void output(FILE* file_pt,
const unsigned& n_plot)
1759 std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list)
const;
1780 for (
unsigned i = 0;
i < 9;
i++)
1783 dtestdx(
i, 0) = dpsidx(
i, 0);
1784 dtestdx(
i, 1) = dpsidx(
i, 1);
1806 for (
unsigned i = 0;
i < 9;
i++)
1809 dtestdx(
i, 0) = dpsidx(
i, 0);
1810 dtestdx(
i, 1) = dpsidx(
i, 1);
1826 const unsigned& ipt,
1837 ipt, psi, dpsidx, djacobian_dX, d_dpsidx_dX);
1841 for (
unsigned i = 0;
i < 9;
i++)
1845 for (
unsigned k = 0;
k < 2;
k++)
1847 dtestdx(
i,
k) = dpsidx(
i,
k);
1849 for (
unsigned p = 0;
p < 2;
p++)
1851 for (
unsigned q = 0;
q < 9;
q++)
1853 d_dtestdx_dX(
p,
q,
i,
k) = d_dpsidx_dX(
p,
q,
i,
k);
1870 double psi1[2], psi2[2];
1877 for (
unsigned i = 0;
i < 2;
i++)
1879 for (
unsigned j = 0;
j < 2;
j++)
1882 psi[2 *
i +
j] = psi2[
i] * psi1[
j];
1896 for (
unsigned i = 0;
i < 4;
i++)
test[
i] = psi[
i];
1928 template<
class TAYLOR_HOOD_ELEMENT>
1948 unsigned nnod = this->
nnode();
1949 for (
unsigned j = 0;
j < nnod;
j++)
1952 data_values.push_back(std::make_pair(this->
node_pt(
j), fld));
1959 unsigned Pconv_size = 3;
1960 for (
unsigned j = 0;
j < Pconv_size;
j++)
1963 unsigned vertex_index = this->Pconv[
j];
1964 data_values.push_back(
1965 std::make_pair(this->
node_pt(vertex_index), fld));
2008 unsigned n_dim = this->
dim();
2009 unsigned n_node = this->
nnode();
2014 this->pshape_axi_nst(
s, psi);
2016 Shape psif(n_node), testf(n_node);
2017 DShape dpsifdx(n_node, n_dim), dtestfdx(n_node, n_dim);
2020 double J = this->dshape_and_dtest_eulerian_axi_nst(
2021 s, psif, dpsifdx, testf, dtestfdx);
2026 Shape testf(n_node);
2027 DShape dpsifdx(n_node, n_dim), dtestfdx(n_node, n_dim);
2030 double J = this->dshape_and_dtest_eulerian_axi_nst(
2031 s, psi, dpsifdx, testf, dtestfdx);
2040 const unsigned& fld,
2043 unsigned n_node = this->
nnode();
2048 return this->interpolated_p_axi_nst(
s);
2054 unsigned u_nodal_index = this->u_index_axi_nst(fld);
2060 this->
shape(s, psi);
2063 double interpolated_u = 0.0;
2066 for (
unsigned l = 0; l < n_node; l++)
2068 interpolated_u += this->
nodal_value(t, l, u_nodal_index) * psi[l];
2070 return interpolated_u;
2080 return this->npres_axi_nst();
2084 return this->
nnode();
2094 return this->p_local_eqn(
j);
2098 const unsigned u_nodal_index = this->u_index_axi_nst(fld);
2109 template<
class ELEMENT>
2123 template<
class ELEMENT>
2136 template<
class CROUZEIX_RAVIART_ELEMENT>
2156 const unsigned n_node = this->
nnode();
2157 for (
unsigned n = 0;
n < n_node;
n++)
2160 data_values.push_back(std::make_pair(this->
node_pt(
n), fld));
2167 const unsigned n_press = this->npres_axi_nst();
2169 for (
unsigned j = 0;
j < n_press;
j++)
2171 data_values.push_back(std::make_pair(
2215 unsigned n_dim = this->
dim();
2216 unsigned n_node = this->
nnode();
2221 this->pshape_axi_nst(
s, psi);
2223 Shape psif(n_node), testf(n_node);
2224 DShape dpsifdx(n_node, n_dim), dtestfdx(n_node, n_dim);
2227 double J = this->dshape_and_dtest_eulerian_axi_nst(
2228 s, psif, dpsifdx, testf, dtestfdx);
2233 Shape testf(n_node);
2234 DShape dpsifdx(n_node, n_dim), dtestfdx(n_node, n_dim);
2237 double J = this->dshape_and_dtest_eulerian_axi_nst(
2238 s, psi, dpsifdx, testf, dtestfdx);
2247 const unsigned& fld,
2256 return this->interpolated_p_axi_nst(
s);
2261 return this->interpolated_u_axi_nst(
t,
s, fld);
2271 return this->npres_axi_nst();
2275 return this->
nnode();
2285 return this->p_local_eqn(
j);
2289 const unsigned u_nodal_index = this->u_index_axi_nst(fld);
2300 template<
class ELEMENT>
2314 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
float * p
Definition: Tutorial_Map_using.cpp:9
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
FaceGeometry()
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1601
FaceGeometry()
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1920
FaceGeometry()
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:2320
FaceGeometry()
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:2129
FaceGeometry()
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1588
FaceGeometry()
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1907
FaceGeometry()
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:2306
FaceGeometry()
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:2115
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
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:119
unsigned u_index_nst(const unsigned &i) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:556
virtual void get_source_fct_gradient(const double &time, const unsigned &ipt, const Vector< double > &x, Vector< double > &gradient)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:320
double *& re_st_pt()
Pointer to product of Reynolds and Strouhal number (=Womersley number)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:438
void fill_in_contribution_to_djacobian_dparameter(double *const ¶meter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:915
void get_pressure_and_velocity_mass_matrix_diagonal(Vector< double > &press_mass_diag, Vector< double > &veloc_mass_diag, const unsigned &which_one=0)
Definition: generalised_newtonian_axisym_navier_stokes_elements.cc:71
static bool Pre_multiply_by_viscosity_ratio
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:139
virtual void extrapolated_strain_rate(const unsigned &ipt, DenseMatrix< double > &strain_rate) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:665
virtual void dinterpolated_u_axi_nst_ddata(const Vector< double > &s, const unsigned &i, Vector< double > &du_ddata, Vector< unsigned > &global_eqn_number)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1027
const double & re_invro() const
Global Reynolds number multiplied by inverse Rossby number.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:456
virtual unsigned u_index_axi_nst(const unsigned &i) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:546
virtual int p_local_eqn(const unsigned &n) const =0
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: generalised_newtonian_axisym_navier_stokes_elements.h:931
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: generalised_newtonian_axisym_navier_stokes_elements.cc:3369
double interpolated_u_axi_nst(const Vector< double > &s, const unsigned &i) const
Return FE interpolated velocity u[i] at local coordinate s.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:969
void fill_in_contribution_to_hessian_vector_products(Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
Definition: generalised_newtonian_axisym_navier_stokes_elements.cc:4098
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: generalised_newtonian_axisym_navier_stokes_elements.h:506
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: generalised_newtonian_axisym_navier_stokes_elements.h:1193
void output(std::ostream &outfile)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:789
double dissipation() const
Return integral of dissipation over element.
Definition: generalised_newtonian_axisym_navier_stokes_elements.cc:601
virtual void pshape_axi_nst(const Vector< double > &s, Shape &psi, Shape &test) const =0
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Definition: generalised_newtonian_axisym_navier_stokes_elements.cc:163
unsigned n_u_nst() const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:563
double * ReInvRo_pt
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:165
const double & viscosity_ratio() const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:494
void output(FILE *file_pt)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:802
double *& viscosity_ratio_pt()
Pointer to Viscosity Ratio.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:500
Vector< double > *& g_pt()
Pointer to Vector of gravitational components.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:474
double kin_energy() const
Get integral of kinetic energy over element.
Definition: generalised_newtonian_axisym_navier_stokes_elements.cc:1025
virtual void extrapolated_strain_rate(const Vector< double > &s, DenseMatrix< double > &strain_rate) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.cc:977
double pressure_integral() const
Integral of pressure over element.
Definition: generalised_newtonian_axisym_navier_stokes_elements.cc:1075
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)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:426
static Vector< double > Gamma
Vector to decide whether the stress-divergence form is used or not.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:415
virtual unsigned npres_axi_nst() const =0
Function to return number of pressure degrees of freedom.
void enable_ALE()
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:608
void max_and_min_invariant_and_viscosity(double &min_invariant, double &max_invariant, double &min_viscosity, double &max_viscosity) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.cc:1119
void use_current_strainrate_to_compute_second_invariant()
Switch to fully implict evaluation of non-Newtonian const eqn.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:526
std::string scalar_name_paraview(const unsigned &i) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:741
const double & re_invfr() const
Global inverse Froude number.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:444
double compute_physical_size() const
Compute the volume of the element.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1251
bool Use_extrapolated_strainrate_to_compute_second_invariant
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:183
double *& density_ratio_pt()
Pointer to Density ratio.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:487
GeneralisedNewtonianConstitutiveEquation< 3 > * Constitutive_eqn_pt
Pointer to the generalised Newtonian constitutive equation.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:179
double * Viscosity_Ratio_pt
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:145
static Vector< double > Default_Gravity_vector
Static default value for the gravity vector.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:134
double * ReInvFr_pt
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:161
static int Pressure_not_stored_at_node
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:123
double(*&)(const double &time, const Vector< double > &x) source_fct_pt()
Access function for the source-function pointer.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:514
void output_veloc(std::ostream &outfile, const unsigned &nplot, const unsigned &t)
Definition: generalised_newtonian_axisym_navier_stokes_elements.cc:430
double * Re_pt
Pointer to global Reynolds number.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:154
GeneralisedNewtonianAxisymmetricNavierStokesEquations()
Constructor: NULL the body force and source function.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:393
const double & re() const
Reynolds number.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:420
virtual void get_dresidual_dnodal_coordinates(RankThreeTensor< double > &dresidual_dnodal_coordinates)
Definition: generalised_newtonian_axisym_navier_stokes_elements.cc:2421
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
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: generalised_newtonian_axisym_navier_stokes_elements.cc:318
GeneralisedNewtonianConstitutiveEquation< 3 > *& constitutive_eqn_pt()
Access function for the constitutive equation pointer.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:520
void fill_in_contribution_to_dresiduals_dparameter(double *const ¶meter_pt, Vector< double > &dres_dparam)
Compute the element's residual Vector.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:900
double du_dt_axi_nst(const unsigned &n, const unsigned &i) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:570
void point_output_data(const Vector< double > &s, Vector< double > &data)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:768
double * Density_Ratio_pt
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:149
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: generalised_newtonian_axisym_navier_stokes_elements.h:301
void traction(const Vector< double > &s, const Vector< double > &N, Vector< double > &traction) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.cc:657
virtual double dshape_and_dtest_eulerian_axi_nst(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Compute the element's residual Vector.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:854
virtual void get_viscosity_ratio_axisym_nst(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, double &visc_ratio)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:357
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:867
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: generalised_newtonian_axisym_navier_stokes_elements.h:237
void(* Body_force_fct_pt)(const double &time, const Vector< double > &x, Vector< double > &result)
Pointer to body force function.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:171
double *& re_invro_pt()
Pointer to global inverse Froude number.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:462
double interpolated_duds_axi_nst(const Vector< double > &s, const unsigned &i, const unsigned &j) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1102
double * ReSt_pt
Pointer to global Reynolds number x Strouhal number (=Womersley)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:157
virtual int p_nodal_index_axi_nst() const
Which nodal value represents the pressure?
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:618
unsigned nscalar_paraview() const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:693
double(* Source_fct_pt)(const double &time, const Vector< double > &x)
Pointer to volumetric source function.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:176
virtual double p_axi_nst(const unsigned &n_p) const =0
static double Default_Physical_Constant_Value
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:127
static double Default_Physical_Ratio_Value
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:131
double *& re_invfr_pt()
Pointer to global inverse Froude number.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:450
double interpolated_p_axi_nst(const Vector< double > &s) const
Return FE interpolated pressure at local coordinate s.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1080
const double & density_ratio() const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:481
double interpolated_dudx_axi_nst(const Vector< double > &s, const unsigned &i, const unsigned &j) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1134
double interpolated_dudt_axi_nst(const Vector< double > &s, const unsigned &i) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1166
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: generalised_newtonian_axisym_navier_stokes_elements.h:996
virtual void fill_in_generic_residual_contribution_axi_nst(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: generalised_newtonian_axisym_navier_stokes_elements.cc:1167
void scalar_value_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:700
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:882
void use_extrapolated_strainrate_to_compute_second_invariant()
Use extrapolation for non-Newtonian const eqn.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:532
Vector< double > * G_pt
Pointer to global gravity Vector.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:168
virtual double dshape_and_dtest_eulerian_at_knot_axi_nst(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
double *& re_pt()
Pointer to Reynolds number.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:432
const Vector< double > & g() const
Vector of gravitational components.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:468
bool ALE_is_disabled
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:188
void disable_ALE()
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:599
void strain_rate(const Vector< double > &s, DenseMatrix< double > &strain_rate) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.cc:749
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: generalised_newtonian_axisym_navier_stokes_elements.h:262
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: generalised_newtonian_axisym_navier_stokes_elements.h:944
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1313
void pshape_axi_nst(const Vector< double > &s, Shape &psi) const
Pressure shape functions at local coordinate s.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1561
double p_axi_nst(const unsigned &i) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1382
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: generalised_newtonian_axisym_navier_stokes_elements.cc:4581
double dshape_and_dtest_eulerian_at_knot_axi_nst(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1490
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.cc:4499
static const unsigned Initial_Nvalue[]
Static array of ints to hold required number of variables at nodes.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1316
double dshape_and_dtest_eulerian_axi_nst(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1464
int p_local_eqn(const unsigned &n) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1407
void output(std::ostream &outfile, const unsigned &n_plot)
Redirect output to NavierStokesEquations output.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1419
virtual unsigned required_nvalue(const unsigned &n) const
Number of values (pinned or dofs) required at local node n.
Definition: generalised_newtonian_axisym_navier_stokes_elements.cc:4573
void output(FILE *file_pt)
Redirect output to NavierStokesEquations output.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1427
GeneralisedNewtonianAxisymmetricQCrouzeixRaviartElement()
Constructor, there are three internal values (for the pressure)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1367
void fix_pressure(const unsigned &p_dof, const double &pvalue)
Function to fix the internal pressure dof idof_internal.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1394
void output(std::ostream &outfile)
Redirect output to NavierStokesEquations output.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1413
void output(FILE *file_pt, const unsigned &n_plot)
Redirect output to NavierStokesEquations output.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1433
unsigned npres_axi_nst() const
Return number of pressure values.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1388
unsigned ndof_types() const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1441
unsigned P_axi_nst_internal_index
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1321
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1616
virtual unsigned required_nvalue(const unsigned &n) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1676
void output(FILE *file_pt, const unsigned &n_plot)
Redirect output to NavierStokesEquations output.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1739
GeneralisedNewtonianAxisymmetricQTaylorHoodElement()
Constructor, no internal data points.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1668
void output(std::ostream &outfile)
Redirect output to NavierStokesEquations output.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1720
static const unsigned Pconv[]
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1624
unsigned ndof_types() const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1747
virtual int p_nodal_index_axi_nst() const
Which nodal value represents the pressure?
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1682
unsigned npres_axi_nst() const
Return number of pressure values.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1695
void output(std::ostream &outfile, const unsigned &n_plot)
Redirect output to NavierStokesEquations output.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1726
void pshape_axi_nst(const Vector< double > &s, Shape &psi) const
Pressure shape functions at local coordinate s.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1867
int p_local_eqn(const unsigned &n) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1714
double dshape_and_dtest_eulerian_at_knot_axi_nst(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1796
static const unsigned Initial_Nvalue[]
Static array of ints to hold number of variables at node.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1619
void output(FILE *file_pt)
Redirect output to NavierStokesEquations output.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1733
double p_axi_nst(const unsigned &n_p) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1689
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: generalised_newtonian_axisym_navier_stokes_elements.cc:4659
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.cc:4604
double dshape_and_dtest_eulerian_axi_nst(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1770
void fix_pressure(const unsigned &n_p, const double &pvalue)
Fix the pressure at local pressure node n_p.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1701
Crouzeix Raviart upgraded to become projectable.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:2139
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:2267
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:2246
unsigned nhistory_values_for_coordinate_projection()
Number of positional history values (includes current value!)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:2204
unsigned nhistory_values_for_projection(const unsigned &fld)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:2190
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:2281
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:2147
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:2211
unsigned nfields_for_projection()
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:2182
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1931
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:2004
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1939
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:2090
unsigned nfields_for_projection()
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1975
unsigned nhistory_values_for_projection(const unsigned &fld)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1983
unsigned nhistory_values_for_coordinate_projection()
Number of positional history values (includes current value!)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:1997
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld.
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:2076
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Definition: generalised_newtonian_axisym_navier_stokes_elements.h:2039
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
Definition: generalised_newtonian_constitutive_models.h:72
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
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
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
const char Y
Definition: test/EulerAngles.cpp:32
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2