28 #ifndef OOMPH_NAVIER_STOKES_ELEMENTS_HEADER
29 #define OOMPH_NAVIER_STOKES_ELEMENTS_HEADER
33 #include <oomph-lib-config.h>
38 #include "../generic/Qelements.h"
39 #include "../generic/fsi.h"
40 #include "../generic/projection.h"
84 template<
class ELEMENT>
97 const bool& called_from_refineable_constructor =
false)
107 if (!called_from_refineable_constructor)
110 if (element_pt->
dim() == 3)
119 throw OomphLibError(
"This flux element will not work correctly "
120 "if nodes are hanging\n",
144 std::ostringstream error_message;
146 <<
"fill_in_contribution_to_residuals() must not be called directly.\n"
147 <<
"since it uses the local equation numbering of the bulk element\n"
148 <<
"which calls the relevant helper function directly.\n";
157 std::ostringstream error_message;
159 <<
"fill_in_contribution_to_jacobian() must not be called directly.\n"
160 <<
"since it uses the local equation numbering of the bulk element\n"
161 <<
"which calls the relevant helper function directly.\n";
173 void output(std::ostream& outfile,
const unsigned& nplot)
188 template<
class ELEMENT>
194 unsigned my_dim = this->dim();
205 unsigned n_intpt = this->integral_pt()->nweight();
209 int local_unknown = 0;
212 ELEMENT* bulk_el_pt =
dynamic_cast<ELEMENT*
>(this->bulk_element_pt());
215 unsigned n_pres = bulk_el_pt->npres_nst();
218 double re = bulk_el_pt->re();
221 Shape psip(n_pres), testp(n_pres);
224 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
227 double w = this->integral_pt()->weight(ipt);
230 for (
unsigned i = 0;
i < my_dim;
i++)
231 s[
i] = this->integral_pt()->knot(ipt,
i);
234 s_bulk = this->local_coordinate_in_bulk(
s);
237 this->outer_unit_normal(ipt, unit_normal);
240 bulk_el_pt->interpolated_u_nst(s_bulk, veloc);
244 for (
unsigned i = 0;
i < my_dim + 1;
i++)
246 flux += veloc[
i] * unit_normal[
i];
253 double interpolated_press = bulk_el_pt->interpolated_p_nst(s_bulk);
256 bulk_el_pt->pshape_nst(s_bulk, psip, testp);
259 double J = this->J_eulerian(
s);
266 for (
unsigned l = 0; l < n_pres; l++)
268 local_eqn = bulk_el_pt->p_local_eqn(l);
273 residuals[local_eqn] -= re *
flux * interpolated_press * testp[l] *
W;
279 for (
unsigned l2 = 0; l2 < n_pres; l2++)
281 local_unknown = bulk_el_pt->p_local_eqn(l2);
284 if (local_unknown >= 0)
286 jacobian(local_eqn, local_unknown) -=
287 re *
flux * psip[l2] * testp[l] *
W;
343 std::map<
Data*, std::vector<int>>& eqn_number_backup) = 0;
349 const unsigned& face_index) = 0;
364 const unsigned& which_one = 0) = 0;
391 template<
unsigned DIM>
488 DShape& dtestdx)
const = 0;
498 DShape& dtestdx)
const = 0;
527 for (
unsigned i = 0;
i <
DIM;
i++)
535 (*Body_force_fct_pt)(time,
x, result);
563 for (
unsigned i = 0;
i <
DIM;
i++)
567 for (
unsigned j = 0;
j <
DIM;
j++)
569 d_body_force_dx(
j,
i) =
621 double source_pls = 0.0;
623 for (
unsigned i = 0;
i <
DIM;
i++)
627 gradient[
i] = (source_pls -
source) / eps_fd;
660 double*
const& parameter_pt,
703 const double&
re()
const
841 DShape& dptestdx)
const = 0;
848 double u_nst(
const unsigned&
n,
const unsigned&
i)
const
855 double u_nst(
const unsigned&
t,
const unsigned&
n,
const unsigned&
i)
const
892 const unsigned u_nodal_index = this->
u_index_nst(i);
897 for (
unsigned t = 0;
t < n_time;
t++)
925 virtual double p_nst(
const unsigned& n_p)
const = 0;
928 virtual double p_nst(
const unsigned&
t,
const unsigned& n_p)
const = 0;
931 virtual void fix_pressure(
const unsigned& p_dof,
const double& p_value) = 0;
1010 const unsigned& which_one = 0);
1023 const unsigned& nplot)
const
1031 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
1052 std::stringstream error_stream;
1053 error_stream <<
"These Navier Stokes elements only store " <<
DIM + 1
1055 <<
"but i is currently " <<
i << std::endl;
1068 std::ofstream& file_out,
1070 const unsigned& nplot,
1078 std::stringstream error_stream;
1081 error_stream <<
"These Navier Stokes elements only store " <<
DIM + 1
1082 <<
" fields, but i is currently " <<
i << std::endl;
1100 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
1106 for (
unsigned j = 0;
j <
DIM;
j++)
1116 (*exact_soln_pt)(time, spatial_coordinates,
exact_soln);
1142 std::stringstream error_stream;
1143 error_stream <<
"These Navier Stokes elements only store " <<
DIM + 1
1145 <<
"but i is currently " <<
i << std::endl;
1163 void output(std::ostream& outfile,
const unsigned& nplot);
1175 void output(FILE* file_pt,
const unsigned& nplot);
1189 void full_output(std::ostream& outfile,
const unsigned& nplot);
1196 const unsigned& nplot,
1208 const unsigned& nplot,
1215 const unsigned& nplot,
1290 residuals, jacobian, mass_matrix, 2);
1310 double*
const& parameter_pt,
1326 double*
const& parameter_pt,
1333 parameter_pt, dres_dparam, djac_dparam, dmass_matrix_dparam, 2);
1352 residuals, jacobian, 1);
1358 std::map<
Data*, std::vector<int>>& eqn_number_backup)
1363 for (
unsigned j = 0;
j < nint;
j++)
1366 if (eqn_number_backup[data_pt].
size() == 0)
1368 unsigned nvalue = data_pt->
nvalue();
1369 eqn_number_backup[data_pt].resize(nvalue);
1370 for (
unsigned i = 0;
i < nvalue;
i++)
1373 eqn_number_backup[data_pt][
i] = data_pt->
eqn_number(
i);
1382 unsigned nnod = this->
nnode();
1383 for (
unsigned j = 0;
j < nnod;
j++)
1386 if (eqn_number_backup[nod_pt].
size() == 0)
1388 unsigned nvalue = nod_pt->
nvalue();
1389 eqn_number_backup[nod_pt].resize(nvalue);
1390 for (
unsigned i = 0;
i < nvalue;
i++)
1417 if (solid_nod_pt != 0)
1420 if (eqn_number_backup[solid_posn_data_pt].
size() == 0)
1422 unsigned nvalue = solid_posn_data_pt->
nvalue();
1423 eqn_number_backup[solid_posn_data_pt].resize(nvalue);
1424 for (
unsigned i = 0;
i < nvalue;
i++)
1427 eqn_number_backup[solid_posn_data_pt][
i] =
1431 solid_posn_data_pt->
pin(
i);
1444 const unsigned& face_index) = 0;
1450 std::ostream& outfile)
1453 for (
unsigned e = 0;
e < nel;
e++)
1457 outfile <<
"ZONE" << std::endl;
1462 for (
unsigned ipt = 0; ipt <
n; ipt++)
1464 for (
unsigned i = 0;
i <
DIM - 1;
i++)
1470 for (
unsigned i = 0;
i <
DIM;
i++)
1472 outfile <<
x[
i] <<
" ";
1474 for (
unsigned i = 0;
i <
DIM;
i++)
1476 outfile << unit_normal[
i] <<
" ";
1478 outfile << std::endl;
1489 for (
unsigned e = 0;
e < nel;
e++)
1509 unsigned n_node =
nnode();
1515 for (
unsigned i = 0;
i <
DIM;
i++)
1522 for (
unsigned l = 0; l < n_node; l++)
1533 unsigned n_node =
nnode();
1543 double interpolated_u = 0.0;
1545 for (
unsigned l = 0; l < n_node; l++)
1547 interpolated_u +=
nodal_value(l, u_nodal_index) * psi[l];
1550 return (interpolated_u);
1557 const unsigned&
i)
const
1560 unsigned n_node =
nnode();
1572 double interpolated_u = 0.0;
1574 for (
unsigned l = 0; l < n_node; l++)
1576 interpolated_u +=
nodal_value(
t, l, u_nodal_index) * psi[l];
1579 return (interpolated_u);
1593 unsigned n_node =
nnode();
1603 unsigned n_u_dof = 0;
1604 for (
unsigned l = 0; l < n_node; l++)
1608 if (global_eqn >= 0)
1615 du_ddata.resize(n_u_dof, 0.0);
1616 global_eqn_number.resize(n_u_dof, 0);
1621 for (
unsigned l = 0; l < n_node; l++)
1625 if (global_eqn >= 0)
1628 global_eqn_number[count] = global_eqn;
1630 du_ddata[count] = psi[l];
1649 double interpolated_p = 0.0;
1651 for (
unsigned l = 0; l < n_pres; l++)
1653 interpolated_p +=
p_nst(l) * psi[l];
1656 return (interpolated_p);
1671 double interpolated_p = 0.0;
1673 for (
unsigned l = 0; l < n_pres; l++)
1675 interpolated_p +=
p_nst(
t, l) * psi[l];
1678 return (interpolated_p);
1686 const unsigned&
j)
const
1689 const unsigned n_node =
nnode();
1703 double interpolated_dudx = 0.0;
1706 for (
unsigned l = 0; l < n_node; l++)
1708 interpolated_dudx +=
nodal_value(l, u_nodal_index) * dpsifdx(l,
j);
1711 return (interpolated_dudx);
1720 unsigned dim =
s.size();
1726 for (
unsigned i = 0;
i <
dim;
i++)
1746 template<
unsigned DIM>
1782 const unsigned& ipt,
1824 double p_nst(
const unsigned&
t,
const unsigned&
i)
const
1876 std::set<std::pair<Data*, unsigned>>& paired_load_data);
1887 std::set<std::pair<Data*, unsigned>>& paired_pressure_data);
1897 void output(std::ostream& outfile,
const unsigned& nplot)
1910 void output(FILE* file_pt,
const unsigned& nplot)
1947 std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list)
const;
1957 template<
unsigned DIM>
1966 double J = this->dshape_eulerian(
s, psi, dpsidx);
1979 template<
unsigned DIM>
1981 DIM>::dshape_and_dtest_eulerian_at_knot_nst(
const unsigned& ipt,
1988 double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
2009 const unsigned& ipt,
2019 const double J = this->dshape_eulerian_at_knot(
2020 ipt, psi, dpsidx, djacobian_dX, d_dpsidx_dX);
2024 for (
unsigned i = 0;
i < 9;
i++)
2028 for (
unsigned k = 0;
k < 2;
k++)
2030 dtestdx(
i,
k) = dpsidx(
i,
k);
2032 for (
unsigned p = 0;
p < 2;
p++)
2034 for (
unsigned q = 0;
q < 9;
q++)
2036 d_dtestdx_dX(
p,
q,
i,
k) = d_dpsidx_dX(
p,
q,
i,
k);
2059 const unsigned& ipt,
2069 const double J = this->dshape_eulerian_at_knot(
2070 ipt, psi, dpsidx, djacobian_dX, d_dpsidx_dX);
2074 for (
unsigned i = 0;
i < 27;
i++)
2078 for (
unsigned k = 0;
k < 3;
k++)
2080 dtestdx(
i,
k) = dpsidx(
i,
k);
2082 for (
unsigned p = 0;
p < 3;
p++)
2084 for (
unsigned q = 0;
q < 27;
q++)
2086 d_dtestdx_dX(
p,
q,
i,
k) = d_dpsidx_dX(
p,
q,
i,
k);
2129 dppsidx(0, 0) = 0.0;
2130 dppsidx(1, 0) = 1.0;
2131 dppsidx(2, 0) = 0.0;
2133 dppsidx(0, 1) = 0.0;
2134 dppsidx(1, 1) = 0.0;
2135 dppsidx(2, 1) = 1.0;
2141 dshape_local(
s, psi, dpsi);
2147 const double det = local_to_eulerian_mapping(dpsi, inverse_jacobian);
2151 transform_derivatives(inverse_jacobian, dppsidx);
2165 template<
unsigned DIM>
2171 this->pshape_nst(
s, psi);
2211 dppsidx(0, 0) = 0.0;
2212 dppsidx(1, 0) = 1.0;
2213 dppsidx(2, 0) = 0.0;
2214 dppsidx(3, 0) = 0.0;
2216 dppsidx(0, 1) = 0.0;
2217 dppsidx(1, 1) = 0.0;
2218 dppsidx(2, 1) = 1.0;
2219 dppsidx(3, 1) = 0.0;
2221 dppsidx(0, 2) = 0.0;
2222 dppsidx(1, 2) = 0.0;
2223 dppsidx(2, 2) = 0.0;
2224 dppsidx(3, 2) = 1.0;
2230 dshape_local(
s, psi, dpsi);
2236 const double det = local_to_eulerian_mapping(dpsi, inverse_jacobian);
2240 transform_derivatives(inverse_jacobian, dppsidx);
2305 template<
unsigned DIM>
2340 const unsigned& ipt,
2372 return static_cast<int>(
DIM);
2390 double p_nst(
const unsigned&
t,
const unsigned& n_p)
const
2407 return static_cast<unsigned>(
pow(2.0,
static_cast<int>(
DIM)));
2418 return std::find(this->Pconv, this->Pconv + n_p,
n) != this->Pconv + n_p;
2450 std::set<std::pair<Data*, unsigned>>& paired_load_data);
2462 std::set<std::pair<Data*, unsigned>>& paired_pressure_data);
2472 void output(std::ostream& outfile,
const unsigned& nplot)
2484 void output(FILE* file_pt,
const unsigned& nplot)
2504 std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list)
const;
2514 template<
unsigned DIM>
2523 double J = this->dshape_eulerian(
s, psi, dpsidx);
2539 template<
unsigned DIM>
2541 const unsigned& ipt,
2548 double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
2573 double psi1[2], psi2[2];
2574 double dpsi1[2], dpsi2[2];
2584 for (
unsigned i = 0;
i < 2;
i++)
2586 for (
unsigned j = 0;
j < 2;
j++)
2589 ppsi[2 *
i +
j] = psi2[
i] * psi1[
j];
2590 dppsidx(2 *
i +
j, 0) = psi2[
i] * dpsi1[
j];
2591 dppsidx(2 *
i +
j, 1) = dpsi2[
i] * psi1[
j];
2599 dshape_local(
s, psi, dpsi);
2605 const double det = local_to_eulerian_mapping(dpsi, inverse_jacobian);
2609 transform_derivatives(inverse_jacobian, dppsidx);
2631 const unsigned& ipt,
2641 const double J = this->dshape_eulerian_at_knot(
2642 ipt, psi, dpsidx, djacobian_dX, d_dpsidx_dX);
2646 for (
unsigned i = 0;
i < 9;
i++)
2650 for (
unsigned k = 0;
k < 2;
k++)
2652 dtestdx(
i,
k) = dpsidx(
i,
k);
2654 for (
unsigned p = 0;
p < 2;
p++)
2656 for (
unsigned q = 0;
q < 9;
q++)
2658 d_dtestdx_dX(
p,
q,
i,
k) = d_dpsidx_dX(
p,
q,
i,
k);
2680 const unsigned& ipt,
2690 const double J = this->dshape_eulerian_at_knot(
2691 ipt, psi, dpsidx, djacobian_dX, d_dpsidx_dX);
2695 for (
unsigned i = 0;
i < 27;
i++)
2699 for (
unsigned k = 0;
k < 3;
k++)
2701 dtestdx(
i,
k) = dpsidx(
i,
k);
2703 for (
unsigned p = 0;
p < 3;
p++)
2705 for (
unsigned q = 0;
q < 27;
q++)
2707 d_dtestdx_dX(
p,
q,
i,
k) = d_dpsidx_dX(
p,
q,
i,
k);
2727 double psi1[2], psi2[2];
2734 for (
unsigned i = 0;
i < 2;
i++)
2736 for (
unsigned j = 0;
j < 2;
j++)
2739 psi[2 *
i +
j] = psi2[
i] * psi1[
j];
2759 double psi1[2], psi2[2], psi3[2];
2760 double dpsi1[2], dpsi2[2], dpsi3[2];
2772 for (
unsigned i = 0;
i < 2;
i++)
2774 for (
unsigned j = 0;
j < 2;
j++)
2776 for (
unsigned k = 0;
k < 2;
k++)
2779 ppsi[4 *
i + 2 *
j +
k] = psi3[
i] * psi2[
j] * psi1[
k];
2780 dppsidx(4 *
i + 2 *
j +
k, 0) = psi3[
i] * psi2[
j] * dpsi1[
k];
2781 dppsidx(4 *
i + 2 *
j +
k, 1) = psi3[
i] * dpsi2[
j] * psi1[
k];
2782 dppsidx(4 *
i + 2 *
j +
k, 2) = dpsi3[
i] * psi2[
j] * psi1[
k];
2791 dshape_local(
s, psi, dpsi);
2797 const double det = local_to_eulerian_mapping(dpsi, inverse_jacobian);
2801 transform_derivatives(inverse_jacobian, dppsidx);
2820 double psi1[2], psi2[2], psi3[2];
2829 for (
unsigned i = 0;
i < 2;
i++)
2831 for (
unsigned j = 0;
j < 2;
j++)
2833 for (
unsigned k = 0;
k < 2;
k++)
2836 psi[4 *
i + 2 *
j +
k] = psi3[
i] * psi2[
j] * psi1[
k];
2846 template<
unsigned DIM>
2852 this->pshape_nst(
s, psi);
2911 template<
class TAYLOR_HOOD_ELEMENT>
2933 if (fld < this->dim())
2936 unsigned nnod = this->nnode();
2937 for (
unsigned j = 0;
j < nnod;
j++)
2940 data_values.push_back(std::make_pair(this->node_pt(
j), fld));
2947 unsigned Pconv_size = this->dim() + 1;
2948 for (
unsigned j = 0;
j < Pconv_size;
j++)
2951 unsigned vertex_index = this->Pconv[
j];
2952 data_values.push_back(
2953 std::make_pair(this->node_pt(vertex_index), fld));
2965 return this->dim() + 1;
2974 if (fld == this->dim())
2977 return this->node_pt(0)->ntstorage();
2981 return this->node_pt(0)->ntstorage();
2989 return this->node_pt(0)->position_time_stepper_pt()->ntstorage();
2998 unsigned n_dim = this->dim();
2999 unsigned n_node = this->nnode();
3004 this->pshape_nst(
s, psi);
3006 Shape psif(n_node), testf(n_node);
3007 DShape dpsifdx(n_node, n_dim), dtestfdx(n_node, n_dim);
3010 double J = this->dshape_and_dtest_eulerian_nst(
3011 s, psif, dpsifdx, testf, dtestfdx);
3016 Shape testf(n_node);
3017 DShape dpsifdx(n_node, n_dim), dtestfdx(n_node, n_dim);
3021 this->dshape_and_dtest_eulerian_nst(
s, psi, dpsifdx, testf, dtestfdx);
3030 const unsigned& fld,
3033 unsigned n_dim = this->dim();
3034 unsigned n_node = this->nnode();
3039 return this->interpolated_p_nst(
t,
s);
3045 unsigned u_nodal_index = this->u_index_nst(fld);
3051 this->
shape(s, psi);
3054 double interpolated_u = 0.0;
3057 for (
unsigned l = 0; l < n_node; l++)
3059 interpolated_u += this->nodal_value(
t, l, u_nodal_index) * psi[l];
3061 return interpolated_u;
3069 if (fld == this->dim())
3071 return this->npres_nst();
3075 return this->nnode();
3083 if (fld == this->dim())
3085 return this->p_local_eqn(
j);
3089 const unsigned u_nodal_index = this->u_index_nst(fld);
3090 return this->nodal_local_eqn(
j, u_nodal_index);
3100 template<
class ELEMENT>
3113 template<
class ELEMENT>
3125 template<
class CROUZEIX_RAVIART_ELEMENT>
3146 if (fld < this->dim())
3149 const unsigned n_node = this->nnode();
3150 for (
unsigned n = 0;
n < n_node;
n++)
3153 data_values.push_back(std::make_pair(this->node_pt(
n), fld));
3160 const unsigned n_press = this->npres_nst();
3162 for (
unsigned j = 0;
j < n_press;
j++)
3164 data_values.push_back(std::make_pair(
3165 this->internal_data_pt(this->P_nst_internal_index),
j));
3177 return this->dim() + 1;
3186 if (fld == this->dim())
3193 return this->node_pt(0)->ntstorage();
3201 return this->node_pt(0)->position_time_stepper_pt()->ntstorage();
3210 unsigned n_dim = this->dim();
3211 unsigned n_node = this->nnode();
3216 this->pshape_nst(
s, psi);
3218 Shape psif(n_node), testf(n_node);
3219 DShape dpsifdx(n_node, n_dim), dtestfdx(n_node, n_dim);
3222 double J = this->dshape_and_dtest_eulerian_nst(
3223 s, psif, dpsifdx, testf, dtestfdx);
3228 Shape testf(n_node);
3229 DShape dpsifdx(n_node, n_dim), dtestfdx(n_node, n_dim);
3233 this->dshape_and_dtest_eulerian_nst(
s, psi, dpsifdx, testf, dtestfdx);
3242 const unsigned& fld,
3245 unsigned n_dim = this->dim();
3250 return this->interpolated_p_nst(
s);
3255 return this->interpolated_u_nst(
t,
s, fld);
3263 if (fld == this->dim())
3265 return this->npres_nst();
3269 return this->nnode();
3277 if (fld == this->dim())
3279 return this->p_local_eqn(
j);
3283 const unsigned u_nodal_index = this->u_index_nst(fld);
3284 return this->nodal_local_eqn(
j, u_nodal_index);
3294 template<
class ELEMENT>
3307 template<
class ELEMENT>
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
Array< double, 1, 3 > e(1./3., 0.5, 2.)
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
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
unsigned nvalue() const
Return number of values stored in data object (incl pinned ones).
Definition: nodes.h:483
double value(const unsigned &i) const
Definition: nodes.h:293
Definition: elements.h:4338
int & face_index()
Definition: elements.h:4626
void outer_unit_normal(const Vector< double > &s, Vector< double > &unit_normal) const
Compute outer unit normal at the specified local coordinate.
Definition: elements.cc:6006
double interpolated_x(const Vector< double > &s, const unsigned &i) const
Definition: elements.h:4528
FaceGeometry()
Definition: navier_stokes_elements.h:3312
FaceGeometry()
Definition: navier_stokes_elements.h:3118
FaceGeometry()
Definition: navier_stokes_elements.h:2279
FaceGeometry()
Definition: navier_stokes_elements.h:2291
FaceGeometry()
Definition: navier_stokes_elements.h:2887
FaceGeometry()
Definition: navier_stokes_elements.h:2899
FaceGeometry()
Definition: navier_stokes_elements.h:3299
FaceGeometry()
Definition: navier_stokes_elements.h:3105
FaceGeometry()
Definition: navier_stokes_elements.h:2258
FaceGeometry()
Definition: navier_stokes_elements.h:2268
FaceGeometry()
Definition: navier_stokes_elements.h:2865
FaceGeometry()
Definition: navier_stokes_elements.h:2875
Definition: elements.h:4998
Definition: elements.h:1313
virtual unsigned nplot_points_paraview(const unsigned &nplot) const
Definition: elements.h:2862
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 void output(std::ostream &outfile)
Definition: elements.h:3050
double size() const
Definition: elements.cc:4290
virtual void build_face_element(const int &face_index, FaceElement *face_element_pt)
Definition: elements.cc:5132
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 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
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
bool has_hanging_nodes() const
Definition: elements.h:2470
Definition: navier_stokes_elements.h:53
FpPressureAdvDiffRobinBCElementBase()
Constructor.
Definition: navier_stokes_elements.h:56
virtual ~FpPressureAdvDiffRobinBCElementBase()
Empty virtual destructor.
Definition: navier_stokes_elements.h:59
virtual void fill_in_generic_residual_contribution_fp_press_adv_diff_robin_bc(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)=0
Definition: navier_stokes_elements.h:89
void output(std::ostream &outfile)
Overload the output function.
Definition: navier_stokes_elements.h:167
virtual void fill_in_generic_residual_contribution_fp_press_adv_diff_robin_bc(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Definition: navier_stokes_elements.h:190
FpPressureAdvDiffRobinBCElement(FiniteElement *const &element_pt, const int &face_index, const bool &called_from_refineable_constructor=false)
Definition: navier_stokes_elements.h:94
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
This function returns the residuals and the jacobian.
Definition: navier_stokes_elements.h:154
void fill_in_contribution_to_residuals(Vector< double > &residuals)
This function returns just the residuals.
Definition: navier_stokes_elements.h:142
void output(std::ostream &outfile, const unsigned &nplot)
Output function: x,y,[z],u,v,[w],p in tecplot format.
Definition: navier_stokes_elements.h:173
~FpPressureAdvDiffRobinBCElement()
Empty destructor.
Definition: navier_stokes_elements.h:132
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
unsigned ninternal_data() const
Return the number of internal data objects.
Definition: elements.h:823
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.
Definition: matrices.h:74
Definition: elements.h:5231
Definition: navier_stokes_elements.h:395
unsigned nscalar_paraview() const
Definition: navier_stokes_elements.h:1014
void(* NavierStokesBodyForceFctPt)(const double &time, const Vector< double > &x, Vector< double > &body_force)
Definition: navier_stokes_elements.h:399
double * Viscosity_Ratio_pt
Definition: navier_stokes_elements.h:436
void disable_ALE()
Definition: navier_stokes_elements.h:909
virtual double p_nst(const unsigned &t, const unsigned &n_p) const =0
Pressure at local pressure "node" n_p at time level t.
double interpolated_u_nst(const unsigned &t, const Vector< double > &s, const unsigned &i) const
Definition: navier_stokes_elements.h:1555
void get_vorticity(const Vector< double > &s, Vector< double > &vorticity) const
Compute the vorticity vector at local coordinate s.
NavierStokesSourceFctPt Source_fct_pt
Pointer to volumetric source function.
Definition: navier_stokes_elements.h:461
double kin_energy() const
Get integral of kinetic energy over element.
Definition: navier_stokes_elements.cc:1416
virtual unsigned u_index_nst(const unsigned &i) const
Definition: navier_stokes_elements.h:866
virtual void build_fp_press_adv_diff_robin_bc_element(const unsigned &face_index)=0
double * Re_pt
Pointer to global Reynolds number.
Definition: navier_stokes_elements.h:445
double(* NavierStokesSourceFctPt)(const double &time, const Vector< double > &x)
Definition: navier_stokes_elements.h:405
NavierStokesSourceFctPt & source_fct_pt()
Access function for the source-function pointer.
Definition: navier_stokes_elements.h:789
static int Pressure_not_stored_at_node
Definition: navier_stokes_elements.h:418
double *& density_ratio_pt()
Pointer to Density ratio.
Definition: navier_stokes_elements.h:747
NavierStokesBodyForceFctPt body_force_fct_pt() const
Access function for the body-force pointer. Const version.
Definition: navier_stokes_elements.h:783
void fill_in_contribution_to_dresiduals_dparameter(double *const ¶meter_pt, Vector< double > &dres_dparam)
Compute the element's residual Vector.
Definition: navier_stokes_elements.h:1294
void fill_in_contribution_to_hessian_vector_products(Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
Definition: navier_stokes_elements.cc:2105
void get_pressure_and_velocity_mass_matrix_diagonal(Vector< double > &press_mass_diag, Vector< double > &veloc_mass_diag, const unsigned &which_one=0)
Definition: navier_stokes_elements.cc:71
void get_vorticity(const Vector< double > &s, double &vorticity) const
Compute the scalar vorticity at local coordinate s (2D)
virtual double dpshape_and_dptest_eulerian_nst(const Vector< double > &s, Shape &ppsi, DShape &dppsidx, Shape &ptest, DShape &dptestdx) const =0
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Compute the element's residual Vector.
Definition: navier_stokes_elements.h:1260
void interpolated_u_nst(const Vector< double > &s, Vector< double > &veloc) const
Compute vector of FE interpolated velocity u at local coordinate s.
Definition: navier_stokes_elements.h:1505
virtual int p_nodal_index_nst() const
Definition: navier_stokes_elements.h:936
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Definition: navier_stokes_elements.cc:289
double interpolated_dudx_nst(const Vector< double > &s, const unsigned &i, const unsigned &j) const
Definition: navier_stokes_elements.h:1684
void fill_in_pressure_advection_diffusion_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: navier_stokes_elements.h:1348
virtual double dshape_and_dtest_eulerian_at_knot_nst(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
virtual void pshape_nst(const Vector< double > &s, Shape &psi) const =0
Compute the pressure shape functions at local coordinate s.
virtual unsigned npres_nst() const =0
Function to return number of pressure degrees of freedom.
NavierStokesPressureAdvDiffSourceFctPt Press_adv_diff_source_fct_pt
Definition: navier_stokes_elements.h:465
virtual void fill_in_generic_pressure_advection_diffusion_contribution_nst(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Definition: navier_stokes_elements.cc:1605
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
Definition: navier_stokes_elements.h:1283
virtual void pshape_nst(const Vector< double > &s, Shape &psi, Shape &test) const =0
virtual void get_dresidual_dnodal_coordinates(RankThreeTensor< double > &dresidual_dnodal_coordinates)
Definition: navier_stokes_elements.cc:2123
virtual void get_source_gradient_nst(const double &time, const unsigned &ipt, const Vector< double > &x, Vector< double > &gradient)
Definition: navier_stokes_elements.h:607
void scalar_value_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
Definition: navier_stokes_elements.h:1021
double interpolated_u_nst(const Vector< double > &s, const unsigned &i) const
Return FE interpolated velocity u[i] at local coordinate s.
Definition: navier_stokes_elements.h:1530
double *& viscosity_ratio_pt()
Pointer to Viscosity Ratio.
Definition: navier_stokes_elements.h:734
double u_nst(const unsigned &n, const unsigned &i) const
Definition: navier_stokes_elements.h:848
NavierStokesPressureAdvDiffSourceFctPt & source_fct_for_pressure_adv_diff()
Definition: navier_stokes_elements.h:802
void scalar_value_fct_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) const
Definition: navier_stokes_elements.h:1067
void strain_rate(const Vector< double > &s, DenseMatrix< double > &strain_rate) const
Strain-rate tensor: 1/2 (du_i/dx_j + du_j/dx_i)
Definition: navier_stokes_elements.cc:1185
NavierStokesPressureAdvDiffSourceFctPt source_fct_for_pressure_adv_diff() const
Definition: navier_stokes_elements.h:809
double d_kin_energy_dt() const
Get integral of time derivative of kinetic energy over element.
Definition: navier_stokes_elements.cc:1460
const double & re_st() const
Product of Reynolds and Strouhal number (=Womersley number)
Definition: navier_stokes_elements.h:709
virtual void get_body_force_nst(const double &time, const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &result)
Definition: navier_stokes_elements.h:517
double interpolated_p_nst(const unsigned &t, const Vector< double > &s) const
Return FE interpolated pressure at local coordinate s at time level t.
Definition: navier_stokes_elements.h:1661
void enable_ALE()
Definition: navier_stokes_elements.h:918
virtual void get_body_force_gradient_nst(const double &time, const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, DenseMatrix< double > &d_body_force_dx)
Definition: navier_stokes_elements.h:544
void compute_norm(double &norm)
Compute norm of solution: square of the L2 norm of the velocities.
Definition: navier_stokes_elements.cc:186
double *& re_st_pt()
Pointer to product of Reynolds and Strouhal number (=Womersley number)
Definition: navier_stokes_elements.h:721
void fill_in_contribution_to_djacobian_dparameter(double *const ¶meter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
Definition: navier_stokes_elements.h:1309
NavierStokesBodyForceFctPt & body_force_fct_pt()
Access function for the body-force pointer.
Definition: navier_stokes_elements.h:777
const double & re_invfr() const
Global inverse Froude number.
Definition: navier_stokes_elements.h:753
int & pinned_fp_pressure_eqn()
Definition: navier_stokes_elements.h:817
int Pinned_fp_pressure_eqn
Definition: navier_stokes_elements.h:479
const double & viscosity_ratio() const
Definition: navier_stokes_elements.h:728
virtual void fix_pressure(const unsigned &p_dof, const double &p_value)=0
Pin p_dof-th pressure dof and set it to value specified by p_value.
NavierStokesSourceFctPt source_fct_pt() const
Access function for the source-function pointer. Const version.
Definition: navier_stokes_elements.h:795
virtual void dinterpolated_u_nst_ddata(const Vector< double > &s, const unsigned &i, Vector< double > &du_ddata, Vector< unsigned > &global_eqn_number)
Definition: navier_stokes_elements.h:1587
Vector< double > * G_pt
Pointer to global gravity Vector.
Definition: navier_stokes_elements.h:455
static double Default_Physical_Constant_Value
Navier–Stokes equations static data.
Definition: navier_stokes_elements.h:422
void output_pressure_advection_diffusion_robin_elements(std::ostream &outfile)
Definition: navier_stokes_elements.h:1449
void full_output(std::ostream &outfile)
Definition: navier_stokes_elements.h:1180
void delete_pressure_advection_diffusion_robin_elements()
Definition: navier_stokes_elements.h:1486
virtual void fill_in_generic_dresidual_contribution_nst(double *const ¶meter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam, unsigned flag)
Definition: navier_stokes_elements.cc:2087
virtual double interpolated_p_nst(const Vector< double > &s) const
Return FE interpolated pressure at local coordinate s.
Definition: navier_stokes_elements.h:1639
double * ReInvFr_pt
Definition: navier_stokes_elements.h:452
Vector< double > *& g_pt()
Pointer to Vector of gravitational components.
Definition: navier_stokes_elements.h:771
void output_veloc(std::ostream &outfile, const unsigned &nplot, const unsigned &t)
Definition: navier_stokes_elements.cc:687
static double Default_Physical_Ratio_Value
Navier–Stokes equations static data.
Definition: navier_stokes_elements.h:426
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: navier_stokes_elements.h:1273
void fill_in_pressure_advection_diffusion_residuals(Vector< double > &residuals)
Definition: navier_stokes_elements.h:1339
void get_traction(const Vector< double > &s, const Vector< double > &N, Vector< double > &traction)
Definition: navier_stokes_elements.cc:1098
void pin_all_non_pressure_dofs(std::map< Data *, std::vector< int >> &eqn_number_backup)
Pin all non-pressure dofs and backup eqn numbers.
Definition: navier_stokes_elements.h:1357
virtual double get_source_nst(const double &time, const unsigned &ipt, const Vector< double > &x)
Definition: navier_stokes_elements.h:585
double *& re_invfr_pt()
Pointer to global inverse Froude number.
Definition: navier_stokes_elements.h:759
double u_nst(const unsigned &t, const unsigned &n, const unsigned &i) const
Definition: navier_stokes_elements.h:855
const double & re() const
Reynolds number.
Definition: navier_stokes_elements.h:703
virtual void fill_in_generic_residual_contribution_nst(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Definition: navier_stokes_elements.cc:1774
void output(FILE *file_pt)
Definition: navier_stokes_elements.h:1167
double * Density_Ratio_pt
Definition: navier_stokes_elements.h:440
double *& re_pt()
Pointer to Reynolds number.
Definition: navier_stokes_elements.h:715
double pressure_integral() const
Integral of pressure over element.
Definition: navier_stokes_elements.cc:1558
static Vector< double > Default_Gravity_vector
Static default value for the gravity vector.
Definition: navier_stokes_elements.h:429
unsigned n_u_nst() const
Definition: navier_stokes_elements.h:873
double du_dt_nst(const unsigned &n, const unsigned &i) const
Definition: navier_stokes_elements.h:880
NavierStokesBodyForceFctPt Body_force_fct_pt
Pointer to body force function.
Definition: navier_stokes_elements.h:458
virtual double p_nst(const unsigned &n_p) const =0
const double & density_ratio() const
Definition: navier_stokes_elements.h:741
void output(std::ostream &outfile)
Definition: navier_stokes_elements.h:1155
void output_vorticity(std::ostream &outfile, const unsigned &nplot)
Definition: navier_stokes_elements.cc:978
double * ReSt_pt
Pointer to global Reynolds number x Strouhal number (=Womersley)
Definition: navier_stokes_elements.h:448
double(* NavierStokesPressureAdvDiffSourceFctPt)(const Vector< double > &x)
Definition: navier_stokes_elements.h:412
bool ALE_is_disabled
Definition: navier_stokes_elements.h:470
NavierStokesEquations()
Definition: navier_stokes_elements.h:677
virtual double dshape_and_dtest_eulerian_nst(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
const Vector< double > & g() const
Vector of gravitational components.
Definition: navier_stokes_elements.h:765
void get_load(const Vector< double > &s, const Vector< double > &N, Vector< double > &load)
Definition: navier_stokes_elements.h:992
double dissipation() const
Return integral of dissipation over element.
Definition: navier_stokes_elements.cc:1046
Vector< FpPressureAdvDiffRobinBCElementBase * > Pressure_advection_diffusion_robin_element_pt
Definition: navier_stokes_elements.h:475
std::string scalar_name_paraview(const unsigned &i) const
Definition: navier_stokes_elements.h:1127
virtual double dshape_and_dtest_eulerian_at_knot_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 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: navier_stokes_elements.h:1325
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: navier_stokes_elements.cc:573
void point_output_data(const Vector< double > &s, Vector< double > &data)
Definition: navier_stokes_elements.h:1717
static Vector< double > Gamma
Vector to decide whether the stress-divergence form is used or not.
Definition: navier_stokes_elements.h:698
Definition: oomph_definitions.h:222
Definition: elements.h:3439
Crouzeix Raviart upgraded to become projectable.
Definition: navier_stokes_elements.h:3128
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Definition: navier_stokes_elements.h:3140
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
Definition: navier_stokes_elements.h:3275
ProjectableCrouzeixRaviartElement()
Definition: navier_stokes_elements.h:3132
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Definition: navier_stokes_elements.h:3241
unsigned nfields_for_projection()
Definition: navier_stokes_elements.h:3175
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Definition: navier_stokes_elements.h:3206
unsigned nhistory_values_for_coordinate_projection()
Definition: navier_stokes_elements.h:3199
unsigned nhistory_values_for_projection(const unsigned &fld)
Definition: navier_stokes_elements.h:3184
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld.
Definition: navier_stokes_elements.h:3261
Definition: projection.h:183
Taylor Hood upgraded to become projectable.
Definition: navier_stokes_elements.h:2914
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Definition: navier_stokes_elements.h:2927
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Definition: navier_stokes_elements.h:2994
unsigned nhistory_values_for_projection(const unsigned &fld)
Definition: navier_stokes_elements.h:2972
unsigned nhistory_values_for_coordinate_projection()
Definition: navier_stokes_elements.h:2987
unsigned nfields_for_projection()
Definition: navier_stokes_elements.h:2963
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
Definition: navier_stokes_elements.h:3081
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld.
Definition: navier_stokes_elements.h:3067
ProjectableTaylorHoodElement()
Definition: navier_stokes_elements.h:2918
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Definition: navier_stokes_elements.h:3029
Definition: navier_stokes_elements.h:1749
virtual unsigned required_nvalue(const unsigned &n) const
Number of values (pinned or dofs) required at local node n.
Definition: navier_stokes_elements.cc:2581
unsigned npres_nst() const
Return number of pressure values.
Definition: navier_stokes_elements.h:1830
void output(std::ostream &outfile)
Redirect output to NavierStokesEquations output.
Definition: navier_stokes_elements.h:1891
void pshape_nst(const Vector< double > &s, Shape &psi) const
Pressure shape functions at local coordinate s.
int p_local_eqn(const unsigned &n) const
Return the local equation numbers for the pressure values.
Definition: navier_stokes_elements.h:1845
double dpshape_and_dptest_eulerian_nst(const Vector< double > &s, Shape &ppsi, DShape &dppsidx, Shape &ptest, DShape &dptestdx) const
double dshape_and_dtest_eulerian_at_knot_nst(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: navier_stokes_elements.h:1981
double p_nst(const unsigned &i) const
Definition: navier_stokes_elements.h:1816
unsigned P_nst_internal_index
Definition: navier_stokes_elements.h:1757
void identify_load_data(std::set< std::pair< Data *, unsigned >> &paired_load_data)
Definition: navier_stokes_elements.cc:2598
void build_fp_press_adv_diff_robin_bc_element(const unsigned &face_index)
Definition: navier_stokes_elements.h:1861
void output(std::ostream &outfile, const unsigned &nplot)
Redirect output to NavierStokesEquations output.
Definition: navier_stokes_elements.h:1897
static const unsigned Initial_Nvalue[]
Static array of ints to hold required number of variables at nodes.
Definition: navier_stokes_elements.h:1752
void output(FILE *file_pt, const unsigned &nplot)
Redirect output to NavierStokesEquations output.
Definition: navier_stokes_elements.h:1910
QCrouzeixRaviartElement()
Constructor, there are DIM+1 internal values (for the pressure)
Definition: navier_stokes_elements.h:1794
void full_output(std::ostream &outfile)
Definition: navier_stokes_elements.h:1919
double dshape_and_dtest_eulerian_at_knot_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
void full_output(std::ostream &outfile, const unsigned &nplot)
Definition: navier_stokes_elements.h:1927
void identify_pressure_data(std::set< std::pair< Data *, unsigned >> &paired_pressure_data)
Definition: navier_stokes_elements.cc:2635
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Definition: navier_stokes_elements.cc:2662
double p_nst(const unsigned &t, const unsigned &i) const
Definition: navier_stokes_elements.h:1824
void output(FILE *file_pt)
Redirect output to NavierStokesEquations output.
Definition: navier_stokes_elements.h:1904
unsigned ndof_types() const
Definition: navier_stokes_elements.h:1935
void fix_pressure(const unsigned &p_dof, const double &p_value)
Pin p_dof-th pressure dof and set it to value specified by p_value.
Definition: navier_stokes_elements.h:1851
double dshape_and_dtest_eulerian_nst(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: navier_stokes_elements.h:1958
Definition: Qelements.h:459
Definition: navier_stokes_elements.h:2308
void output(FILE *file_pt, const unsigned &nplot)
Redirect output to NavierStokesEquations output.
Definition: navier_stokes_elements.h:2484
int p_local_eqn(const unsigned &n) const
Return the local equation numbers for the pressure values.
Definition: navier_stokes_elements.h:2376
unsigned ndof_types() const
Definition: navier_stokes_elements.h:2492
double dshape_and_dtest_eulerian_at_knot_nst(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: navier_stokes_elements.h:2540
double dpshape_and_dptest_eulerian_nst(const Vector< double > &s, Shape &ppsi, DShape &dppsidx, Shape &ptest, DShape &dptestdx) const
double dshape_and_dtest_eulerian_nst(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: navier_stokes_elements.h:2515
double p_nst(const unsigned &n_p) const
Definition: navier_stokes_elements.h:2383
void pshape_nst(const Vector< double > &s, Shape &psi, Shape &test) const
Pressure shape and test functions at local coordinte s.
Definition: navier_stokes_elements.h:2847
void identify_pressure_data(std::set< std::pair< Data *, unsigned >> &paired_pressure_data)
Definition: navier_stokes_elements.cc:2798
void output(std::ostream &outfile)
Redirect output to NavierStokesEquations output.
Definition: navier_stokes_elements.h:2466
void identify_load_data(std::set< std::pair< Data *, unsigned >> &paired_load_data)
Definition: navier_stokes_elements.cc:2762
void output(FILE *file_pt)
Redirect output to NavierStokesEquations output.
Definition: navier_stokes_elements.h:2478
void pshape_nst(const Vector< double > &s, Shape &psi) const
Pressure shape functions at local coordinate s.
virtual unsigned required_nvalue(const unsigned &n) const
Definition: navier_stokes_elements.h:2355
virtual int p_nodal_index_nst() const
Set the value at which the pressure is stored in the nodes.
Definition: navier_stokes_elements.h:2370
static const unsigned Pconv[]
Definition: navier_stokes_elements.h:2316
double p_nst(const unsigned &t, const unsigned &n_p) const
Definition: navier_stokes_elements.h:2390
void output(std::ostream &outfile, const unsigned &nplot)
Redirect output to NavierStokesEquations output.
Definition: navier_stokes_elements.h:2472
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Definition: navier_stokes_elements.cc:2825
void build_fp_press_adv_diff_robin_bc_element(const unsigned &face_index)
Definition: navier_stokes_elements.h:2434
unsigned npres_nst() const
Return number of pressure values.
Definition: navier_stokes_elements.h:2405
QTaylorHoodElement()
Constructor, no internal data points.
Definition: navier_stokes_elements.h:2351
void fix_pressure(const unsigned &p_dof, const double &p_value)
Pin p_dof-th pressure dof and set it to value specified by p_value.
Definition: navier_stokes_elements.h:2423
static const unsigned Initial_Nvalue[]
Static array of ints to hold number of variables at node.
Definition: navier_stokes_elements.h:2311
bool is_pressure_node(const unsigned &n) const
Deduce whether or not the provided node is a pressure node.
Definition: navier_stokes_elements.h:2412
double dshape_and_dtest_eulerian_at_knot_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
A Rank 4 Tensor class.
Definition: matrices.h:1701
A Rank 3 Tensor class.
Definition: matrices.h:1370
Definition: refineable_elements.h:97
Data *const & variable_position_pt() const
Pointer to variable_position data (const version)
Definition: nodes.h:1765
Definition: navier_stokes_elements.h:308
TemplateFreeNavierStokesEquationsBase()
Constructor (empty)
Definition: navier_stokes_elements.h:311
virtual void pin_all_non_pressure_dofs(std::map< Data *, std::vector< int >> &eqn_number_backup)=0
Pin all non-pressure dofs and backup eqn numbers of all Data.
virtual ~TemplateFreeNavierStokesEquationsBase()
Virtual destructor (empty)
Definition: navier_stokes_elements.h:314
virtual int p_local_eqn(const unsigned &n) const =0
virtual int p_nodal_index_nst() const =0
virtual void fill_in_pressure_advection_diffusion_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)=0
virtual int & pinned_fp_pressure_eqn()=0
virtual void build_fp_press_adv_diff_robin_bc_element(const unsigned &face_index)=0
virtual void delete_pressure_advection_diffusion_robin_elements()=0
virtual void fill_in_pressure_advection_diffusion_residuals(Vector< double > &residuals)=0
virtual void get_pressure_and_velocity_mass_matrix_diagonal(Vector< double > &press_mass_diag, Vector< double > &veloc_mass_diag, const unsigned &which_one=0)=0
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
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_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 pow(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:625
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 flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
void exact_soln(const double &time, const Vector< double > &x, Vector< double > &soln)
Definition: unstructured_two_d_curved.cc:301
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
void shape(const double &s, double *Psi)
Definition: shape.h:564
void dshape< 2 >(const double &s, double *DPsi)
Derivatives of 1D shape functions specialised to linear order (2 Nodes)
Definition: shape.h:616
void shape< 2 >(const double &s, double *Psi)
1D shape functions specialised to linear order (2 Nodes)
Definition: shape.h:608
@ W
Definition: quadtree.h:63
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