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",
121 OOMPH_CURRENT_FUNCTION,
122 OOMPH_EXCEPTION_LOCATION);
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";
150 error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
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";
163 error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
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];
250 if (flux > 0.0) flux = 0.0;
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) =
570 (body_force_pls[j] - body_force[j]) / eps_fd;
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
832 Shape& test)
const = 0;
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;
1057 OOMPH_CURRENT_FUNCTION,
1058 OOMPH_EXCEPTION_LOCATION);
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;
1086 error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
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);
1119 file_out << exact_soln[
i] << std::endl;
1142 std::stringstream error_stream;
1143 error_stream <<
"These Navier Stokes elements only store " << DIM + 1
1145 <<
"but i is currently " <<
i << std::endl;
1147 error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
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();
1694 DShape dpsifdx(n_node, DIM);
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();
1723 data.resize(2 *
dim + 1);
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())
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);
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(
3177 return this->
dim() + 1;
3186 if (fld == this->
dim())
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);
3294 template<
class ELEMENT>
3307 template<
class ELEMENT>
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
A class that represents a collection of data; each Data object may contain many different individual ...
long & eqn_number(const unsigned &i)
Return the equation number of the i-th stored variable.
void pin(const unsigned &i)
Pin the i-th stored variable.
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
void set_value(const unsigned &i, const double &value_)
Set the i-th stored data value to specified value. The only reason that we require an explicit set fu...
unsigned nvalue() const
Return number of values stored in data object (incl pinned ones).
double value(const unsigned &i) const
Return i-th stored value. This function is not virtual so that it can be inlined. This means that if ...
unsigned ntstorage() const
Return total number of doubles stored per value to record time history of each value (one for steady ...
/////////////////////////////////////////////////////////////////////////
FaceElements are elements that coincide with the faces of higher-dimensional "bulk" elements....
int & face_index()
Index of the face (a number that uniquely identifies the face in the element)
void outer_unit_normal(const Vector< double > &s, Vector< double > &unit_normal) const
Compute outer unit normal at the specified local coordinate.
double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s. Overloaded to get information from bulk...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
A general Finite Element class.
virtual unsigned nplot_points_paraview(const unsigned &nplot) const
Return the number of actual plot points for paraview plot with parameter nplot. Broken virtual; can b...
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
double nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n. Produces suitably interpolated values for hanging nodes...
virtual void output(std::ostream &outfile)
Output the element data — typically the values at the nodes in a format suitable for post-processing.
double size() const
Calculate the size of the element (length, area, volume,...) in Eulerian computational coordinates....
virtual void shape(const Vector< double > &s, Shape &psi) const =0
Calculate the geometric shape functions at local coordinate s. This function must be overloaded for e...
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
unsigned nnode() const
Return the number of nodes.
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Get cector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Compute the geometric shape functions and also first derivatives w.r.t. global coordinates at local c...
virtual void build_face_element(const int &face_index, FaceElement *face_element_pt)
Function for building a lower dimensional FaceElement on the specified face of the FiniteElement....
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as .
bool has_hanging_nodes() const
Return boolean to indicate if any of the element's nodes are geometrically hanging.
Helper class for elements that impose Robin boundary conditions on pressure advection diffusion probl...
FpPressureAdvDiffRobinBCElementBase()
Constructor.
virtual ~FpPressureAdvDiffRobinBCElementBase()
Empty virtual destructor.
virtual void fill_in_generic_residual_contribution_fp_press_adv_diff_robin_bc(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)=0
This function returns the residuals for the traction function. flag=1 (or 0): do (or don't) compute t...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
void output(std::ostream &outfile)
Overload the output function.
virtual void fill_in_generic_residual_contribution_fp_press_adv_diff_robin_bc(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
This function returns the residuals for the traction function. flag=1 (or 0): do (or don't) compute t...
FpPressureAdvDiffRobinBCElement(FiniteElement *const &element_pt, const int &face_index, const bool &called_from_refineable_constructor=false)
Constructor, which takes a "bulk" element and the value of the index and its limit....
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
This function returns the residuals and the jacobian.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
This function returns just the residuals.
void output(std::ostream &outfile, const unsigned &nplot)
Output function: x,y,[z],u,v,[w],p in tecplot format.
~FpPressureAdvDiffRobinBCElement()
Empty destructor.
static double Default_fd_jacobian_step
Double used for the default finite difference step in elemental jacobian calculations.
unsigned add_internal_data(Data *const &data_pt, const bool &fd=true)
Add a (pointer to an) internal data object to the element and return the index required to obtain it ...
Data *& internal_data_pt(const unsigned &i)
Return a pointer to i-th internal data object.
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
unsigned ninternal_data() const
Return the number of internal data objects.
int internal_local_eqn(const unsigned &i, const unsigned &j) const
Return the local equation number corresponding to the j-th value stored at the i-th internal data.
TimeStepper *& time_stepper_pt()
Access function for pointer to time stepper: Null if object is not time-dependent.
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.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
unsigned nscalar_paraview() const
Number of scalars/fields output by this element. Reimplements broken virtual function in base class.
void(* NavierStokesBodyForceFctPt)(const double &time, const Vector< double > &x, Vector< double > &body_force)
Function pointer to body force function fct(t,x,f(x)) x is a Vector!
double * Viscosity_Ratio_pt
Pointer to the viscosity ratio (relative to the viscosity used in the definition of the Reynolds numb...
void disable_ALE()
Disable ALE, i.e. assert the mesh is not moving – you do this at your own risk!
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
Return FE interpolated velocity u[i] at local coordinate s at time level t (t=0: present; t>0: histor...
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.
double kin_energy() const
Get integral of kinetic energy over element.
virtual unsigned u_index_nst(const unsigned &i) const
Return the index at which the i-th unknown velocity component is stored. The default value,...
virtual void build_fp_press_adv_diff_robin_bc_element(const unsigned &face_index)=0
Build FaceElements that apply the Robin boundary condition to the pressure advection diffusion proble...
double * Re_pt
Pointer to global Reynolds number.
double(* NavierStokesSourceFctPt)(const double &time, const Vector< double > &x)
Function pointer to source function fct(t,x) x is a Vector!
NavierStokesSourceFctPt & source_fct_pt()
Access function for the source-function pointer.
static int Pressure_not_stored_at_node
Static "magic" number that indicates that the pressure is not stored at a node.
double *& density_ratio_pt()
Pointer to Density ratio.
NavierStokesBodyForceFctPt body_force_fct_pt() const
Access function for the body-force pointer. Const version.
void fill_in_contribution_to_dresiduals_dparameter(double *const ¶meter_pt, Vector< double > &dres_dparam)
Compute the element's residual Vector.
void fill_in_contribution_to_hessian_vector_products(Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
Compute the hessian tensor vector products required to perform continuation of bifurcations analytica...
void get_pressure_and_velocity_mass_matrix_diagonal(Vector< double > &press_mass_diag, Vector< double > &veloc_mass_diag, const unsigned &which_one=0)
Compute the diagonal of the velocity/pressure mass matrices. If which one=0, both are computed,...
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
Compute the pressure shape and test functions and derivatives w.r.t. global coords at local coordinat...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Compute the element's residual Vector.
void interpolated_u_nst(const Vector< double > &s, Vector< double > &veloc) const
Compute vector of FE interpolated velocity u at local coordinate s.
virtual int p_nodal_index_nst() const
Return the index at which the pressure is stored if it is stored at the nodes. If not stored at the n...
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Validate against exact solution at given time Solution is provided via function pointer....
double interpolated_dudx_nst(const Vector< double > &s, const unsigned &i, const unsigned &j) const
Return FE interpolated derivatives of velocity component u[i] w.r.t spatial global coordinate directi...
void fill_in_pressure_advection_diffusion_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute the residuals and Jacobian for the associated pressure advection diffusion problem....
virtual double dshape_and_dtest_eulerian_at_knot_nst(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Compute the shape functions and derivatives w.r.t. global coords at ipt-th integration point Return J...
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
Pointer to source function pressure advection diffusion equation (only to be used during validation)
virtual void fill_in_generic_pressure_advection_diffusion_contribution_nst(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Compute the residuals for the associated pressure advection diffusion problem. Used by the Fp precond...
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
Add the element's contribution to its residuals vector, jacobian matrix and mass matrix.
virtual void pshape_nst(const Vector< double > &s, Shape &psi, Shape &test) const =0
Compute the pressure shape and test functions at local coordinate s.
virtual void get_dresidual_dnodal_coordinates(RankThreeTensor< double > &dresidual_dnodal_coordinates)
Compute derivatives of elemental residual vector with respect to nodal coordinates....
virtual void get_source_gradient_nst(const double &time, const unsigned &ipt, const Vector< double > &x, Vector< double > &gradient)
Get gradient of source term at (Eulerian) position x. This function is virtual to allow overloading i...
void scalar_value_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
Write values of the i-th scalar field at the plot points. Needs to be implemented for each new specif...
double interpolated_u_nst(const Vector< double > &s, const unsigned &i) const
Return FE interpolated velocity u[i] at local coordinate s.
double *& viscosity_ratio_pt()
Pointer to Viscosity Ratio.
double u_nst(const unsigned &n, const unsigned &i) const
Velocity i at local node n. Uses suitably interpolated value for hanging nodes. The use of u_index_ns...
NavierStokesPressureAdvDiffSourceFctPt & source_fct_for_pressure_adv_diff()
Access function for the source-function pointer for pressure advection diffusion (used for validation...
void scalar_value_fct_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt) const
Write values of the i-th scalar field at the plot points. Needs to be implemented for each new specif...
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)
NavierStokesPressureAdvDiffSourceFctPt source_fct_for_pressure_adv_diff() const
Access function for the source-function pointer for pressure advection diffusion (used for validation...
double d_kin_energy_dt() const
Get integral of time derivative of kinetic energy over element.
const double & re_st() const
Product of Reynolds and Strouhal number (=Womersley number)
virtual void get_body_force_nst(const double &time, const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &result)
Calculate the body force at a given time and local and/or Eulerian position. This function is virtual...
double interpolated_p_nst(const unsigned &t, const Vector< double > &s) const
Return FE interpolated pressure at local coordinate s at time level t.
void enable_ALE()
(Re-)enable ALE, i.e. take possible mesh motion into account when evaluating the time-derivative....
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)
Get gradient of body force term at (Eulerian) position x. This function is virtual to allow overloadi...
void compute_norm(double &norm)
Compute norm of solution: square of the L2 norm of the velocities.
double *& re_st_pt()
Pointer to product of Reynolds and Strouhal number (=Womersley number)
void fill_in_contribution_to_djacobian_dparameter(double *const ¶meter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
Compute the element's residual Vector and the jacobian matrix Virtual function can be overloaded by h...
NavierStokesBodyForceFctPt & body_force_fct_pt()
Access function for the body-force pointer.
const double & re_invfr() const
Global inverse Froude number.
int & pinned_fp_pressure_eqn()
Global eqn number of pressure dof that's pinned in pressure adv diff problem.
int Pinned_fp_pressure_eqn
Global eqn number of pressure dof that's pinned in pressure advection diffusion problem (defaults to ...
const double & viscosity_ratio() const
Viscosity ratio for element: Element's viscosity relative to the viscosity used in the definition of ...
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.
virtual void dinterpolated_u_nst_ddata(const Vector< double > &s, const unsigned &i, Vector< double > &du_ddata, Vector< unsigned > &global_eqn_number)
Compute the derivatives of the i-th component of velocity at point s with respect to all data that ca...
Vector< double > * G_pt
Pointer to global gravity Vector.
static double Default_Physical_Constant_Value
Static default value for the physical constants (all initialised to zero)
void output_pressure_advection_diffusion_robin_elements(std::ostream &outfile)
Output the FaceElements that apply the Robin boundary condition to the pressure advection diffusion p...
void full_output(std::ostream &outfile)
Full output function: x,y,[z],u,v,[w],p,du/dt,dv/dt,[dw/dt],dissipation in tecplot format....
void delete_pressure_advection_diffusion_robin_elements()
Delete the FaceElements that apply the Robin boundary condition to the pressure advection diffusion p...
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)
Compute the derivatives of the residuals for the Navier–Stokes equations with respect to a parameter ...
virtual double interpolated_p_nst(const Vector< double > &s) const
Return FE interpolated pressure at local coordinate s.
double * ReInvFr_pt
Pointer to global Reynolds number x inverse Froude number (= Bond number / Capillary number)
Vector< double > *& g_pt()
Pointer to Vector of gravitational components.
void output_veloc(std::ostream &outfile, const unsigned &nplot, const unsigned &t)
Output function: x,y,[z],u,v,[w] in tecplot format. nplot points in each coordinate direction at time...
static double Default_Physical_Ratio_Value
Static default value for the physical ratios (all are initialised to one)
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute the element's residual Vector and the jacobian matrix Virtual function can be overloaded by h...
void fill_in_pressure_advection_diffusion_residuals(Vector< double > &residuals)
Compute the residuals for the associated pressure advection diffusion problem. Used by the Fp precond...
void get_traction(const Vector< double > &s, const Vector< double > &N, Vector< double > &traction)
Compute traction (on the viscous scale) exerted onto the fluid at local coordinate s....
void pin_all_non_pressure_dofs(std::map< Data *, std::vector< int >> &eqn_number_backup)
Pin all non-pressure dofs and backup eqn numbers.
virtual double get_source_nst(const double &time, const unsigned &ipt, const Vector< double > &x)
Calculate the source fct at given time and Eulerian position.
double *& re_invfr_pt()
Pointer to global inverse Froude number.
double u_nst(const unsigned &t, const unsigned &n, const unsigned &i) const
Velocity i at local node n at timestep t (t=0: present; t>0: previous). Uses suitably interpolated va...
const double & re() const
Reynolds number.
virtual void fill_in_generic_residual_contribution_nst(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Compute the residuals for the Navier–Stokes equations. Flag=1 (or 0): do (or don't) compute the Jacob...
void output(FILE *file_pt)
C-style output function: x,y,[z],u,v,[w],p in tecplot format. Default number of plot points.
double * Density_Ratio_pt
Pointer to the density ratio (relative to the density used in the definition of the Reynolds number)
double *& re_pt()
Pointer to Reynolds number.
double pressure_integral() const
Integral of pressure over element.
static Vector< double > Default_Gravity_vector
Static default value for the gravity vector.
unsigned n_u_nst() const
Return the number of velocity components Used in FluidInterfaceElements.
double du_dt_nst(const unsigned &n, const unsigned &i) const
i-th component of du/dt at local node n. Uses suitably interpolated value for hanging nodes.
NavierStokesBodyForceFctPt Body_force_fct_pt
Pointer to body force function.
virtual double p_nst(const unsigned &n_p) const =0
Pressure at local pressure "node" n_p Uses suitably interpolated value for hanging nodes.
const double & density_ratio() const
Density ratio for element: Element's density relative to the viscosity used in the definition of the ...
void output(std::ostream &outfile)
Output function: x,y,[z],u,v,[w],p in tecplot format. Default number of plot points.
void output_vorticity(std::ostream &outfile, const unsigned &nplot)
Output function: x,y,[z], [omega_x,omega_y,[and/or omega_z]] in tecplot format. nplot points in each ...
double * ReSt_pt
Pointer to global Reynolds number x Strouhal number (=Womersley)
double(* NavierStokesPressureAdvDiffSourceFctPt)(const Vector< double > &x)
Function pointer to source function fct(x) for the pressure advection diffusion equation (only used d...
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when time-derivatives are computed....
NavierStokesEquations()
Constructor: NULL the body force and source function and make sure the ALE terms are included by defa...
virtual double dshape_and_dtest_eulerian_nst(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Compute the shape functions and derivatives w.r.t. global coords at local coordinate s....
const Vector< double > & g() const
Vector of gravitational components.
void get_load(const Vector< double > &s, const Vector< double > &N, Vector< double > &load)
This implements a pure virtual function defined in the FSIFluidElement class. The function computes t...
double dissipation() const
Return integral of dissipation over element.
Vector< FpPressureAdvDiffRobinBCElementBase * > Pressure_advection_diffusion_robin_element_pt
Storage for FaceElements that apply Robin BC for pressure adv diff equation used in Fp preconditioner...
std::string scalar_name_paraview(const unsigned &i) const
Name of the i-th scalar field. Default implementation returns V1 for the first one,...
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
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of ...
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)
Add the element's contribution to its residuals vector, jacobian matrix and mass matrix.
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact solution specified via function pointer at a given number of plot points....
void point_output_data(const Vector< double > &s, Vector< double > &data)
Output solution in data vector at local cordinates s: x,y [,z], u,v,[w], p.
static Vector< double > Gamma
Vector to decide whether the stress-divergence form is used or not.
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
TimeStepper *& position_time_stepper_pt()
Return a pointer to the position timestepper.
An OomphLibError object which should be thrown when an run-time error is encountered....
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
Crouzeix Raviart upgraded to become projectable.
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Specify the values associated with field fld. The information is returned in a vector of pairs which ...
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
ProjectableCrouzeixRaviartElement()
Constructor [this was only required explicitly from gcc 4.5.2 onwards...].
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Return interpolated field fld at local coordinate s, at time level t (t=0: present; t>0: history valu...
unsigned nfields_for_projection()
Number of fields to be projected: dim+1, corresponding to velocity components and pressure.
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Return Jacobian of mapping and shape functions of field fld at local coordinate s.
unsigned nhistory_values_for_coordinate_projection()
Number of positional history values. (Note: count includes current value!)
unsigned nhistory_values_for_projection(const unsigned &fld)
Number of history values to be stored for fld-th field. Whatever the timestepper has set up for the v...
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld.
Wrapper class for projectable elements. Adds "projectability" to the underlying ELEMENT.
///////////////////////////////////////////////////////////////// ///////////////////////////////////...
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Specify the values associated with field fld. The information is returned in a vector of pairs which ...
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Return Jacobian of mapping and shape functions of field fld at local coordinate s.
unsigned nhistory_values_for_projection(const unsigned &fld)
Number of history values to be stored for fld-th field. Whatever the timestepper has set up for the v...
unsigned nhistory_values_for_coordinate_projection()
Number of positional history values (Note: count includes current value!)
unsigned nfields_for_projection()
Number of fields to be projected: dim+1, corresponding to velocity components and pressure.
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld.
ProjectableTaylorHoodElement()
Constructor [this was only required explicitly from gcc 4.5.2 onwards...].
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Return interpolated field fld at local coordinate s, at time level t (t=0: present; t>0: history valu...
/////////////////////////////////////////////////////////////////////////// /////////////////////////...
virtual unsigned required_nvalue(const unsigned &n) const
Number of values (pinned or dofs) required at local node n.
unsigned npres_nst() const
Return number of pressure values.
void output(std::ostream &outfile)
Redirect output to NavierStokesEquations output.
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.
double dpshape_and_dptest_eulerian_nst(const Vector< double > &s, Shape &ppsi, DShape &dppsidx, Shape &ptest, DShape &dptestdx) const
Pressure shape and test functions and their derivs w.r.t. to global coords at local coordinate s (tak...
double dshape_and_dtest_eulerian_at_knot_nst(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Velocity shape and test functions and their derivs w.r.t. to global coords at ipt-th integation point...
double p_nst(const unsigned &i) const
Return the i-th pressure value (Discontinous pressure interpolation – no need to cater for hanging no...
unsigned P_nst_internal_index
Internal index that indicates at which internal data the pressure is stored.
void identify_load_data(std::set< std::pair< Data *, unsigned >> &paired_load_data)
Add to the set paired_load_data pairs containing.
void build_fp_press_adv_diff_robin_bc_element(const unsigned &face_index)
Build FaceElements that apply the Robin boundary condition to the pressure advection diffusion proble...
void output(std::ostream &outfile, const unsigned &nplot)
Redirect output to NavierStokesEquations output.
static const unsigned Initial_Nvalue[]
Static array of ints to hold required number of variables at nodes.
void output(FILE *file_pt, const unsigned &nplot)
Redirect output to NavierStokesEquations output.
QCrouzeixRaviartElement()
Constructor, there are DIM+1 internal values (for the pressure)
void full_output(std::ostream &outfile)
Full output function: x,y,[z],u,v,[w],p,du/dt,dv/dt,[dw/dt],dissipation in tecplot format....
double 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
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of ...
void full_output(std::ostream &outfile, const unsigned &nplot)
Full output function: x,y,[z],u,v,[w],p,du/dt,dv/dt,[dw/dt],dissipation in tecplot format....
void identify_pressure_data(std::set< std::pair< Data *, unsigned >> &paired_pressure_data)
Add to the set paired_pressure_data pairs containing.
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contain...
double p_nst(const unsigned &t, const unsigned &i) const
Return the i-th pressure value (Discontinous pressure interpolation – no need to cater for hanging no...
void output(FILE *file_pt)
Redirect output to NavierStokesEquations output.
unsigned ndof_types() const
The number of "DOF types" that degrees of freedom in this element are sub-divided into: Velocity and ...
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.
double dshape_and_dtest_eulerian_nst(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Velocity shape and test functions and their derivs w.r.t. to global coords at local coordinate s (tak...
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
/////////////////////////////////////////////////////////////////////////
void output(FILE *file_pt, const unsigned &nplot)
Redirect output to NavierStokesEquations output.
int p_local_eqn(const unsigned &n) const
Return the local equation numbers for the pressure values.
unsigned ndof_types() const
Returns the number of "DOF types" that degrees of freedom in this element are sub-divided into: Veloc...
double dshape_and_dtest_eulerian_at_knot_nst(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Velocity shape and test functions and their derivs w.r.t. to global coords at local coordinate s (tak...
double dpshape_and_dptest_eulerian_nst(const Vector< double > &s, Shape &ppsi, DShape &dppsidx, Shape &ptest, DShape &dptestdx) const
Pressure shape and test functions and their derivs w.r.t. to global coords at local coordinate s (tak...
double dshape_and_dtest_eulerian_nst(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Velocity shape and test functions and their derivs w.r.t. to global coords at local coordinate s (tak...
double p_nst(const unsigned &n_p) const
Access function for the pressure values at local pressure node n_p (const version)
void identify_pressure_data(std::set< std::pair< Data *, unsigned >> &paired_pressure_data)
Add to the set paired_pressure_data pairs containing.
void output(std::ostream &outfile)
Redirect output to NavierStokesEquations output.
void identify_load_data(std::set< std::pair< Data *, unsigned >> &paired_load_data)
Add to the set paired_load_data pairs containing.
void output(FILE *file_pt)
Redirect output to NavierStokesEquations output.
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
Number of values (pinned or dofs) required at node n. Can be overwritten for hanging node version.
virtual int p_nodal_index_nst() const
Set the value at which the pressure is stored in the nodes.
static const unsigned Pconv[]
Static array of ints to hold conversion from pressure node numbers to actual node numbers.
double p_nst(const unsigned &t, const unsigned &n_p) const
Access function for the pressure values at local pressure node n_p (const version)
void output(std::ostream &outfile, const unsigned &nplot)
Redirect output to NavierStokesEquations output.
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contain...
void build_fp_press_adv_diff_robin_bc_element(const unsigned &face_index)
Build FaceElements that apply the Robin boundary condition to the pressure advection diffusion proble...
unsigned npres_nst() const
Return number of pressure values.
QTaylorHoodElement()
Constructor, no internal data points.
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.
static const unsigned Initial_Nvalue[]
Static array of ints to hold number of variables at node.
bool is_pressure_node(const unsigned &n) const
Deduce whether or not the provided node is a pressure node.
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
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of ...
////////////////////////////////////////////////////////////////// //////////////////////////////////...
////////////////////////////////////////////////////////////////// //////////////////////////////////...
RefineableElements are FiniteElements that may be subdivided into children to provide a better local ...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
A Class for nodes that deform elastically (i.e. position is an unknown in the problem)....
Data *const & variable_position_pt() const
Pointer to variable_position data (const version)
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
TemplateFreeNavierStokesEquationsBase()
Constructor (empty)
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)
virtual int p_local_eqn(const unsigned &n) const =0
Access function for the local equation number information for the pressure. p_local_eqn[n] = local eq...
virtual int p_nodal_index_nst() const =0
Return the index at which the pressure is stored if it is stored at the nodes. If not stored at the n...
virtual void fill_in_pressure_advection_diffusion_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)=0
Compute the residuals and Jacobian for the associated pressure advection diffusion problem....
virtual int & pinned_fp_pressure_eqn()=0
Global eqn number of pressure dof that's pinned in pressure adv diff problem.
virtual void build_fp_press_adv_diff_robin_bc_element(const unsigned &face_index)=0
Build FaceElements that apply the Robin boundary condition to the pressure advection diffusion proble...
virtual void delete_pressure_advection_diffusion_robin_elements()=0
Delete the FaceElements that apply the Robin boundary condition to the pressure advection diffusion p...
virtual void fill_in_pressure_advection_diffusion_residuals(Vector< double > &residuals)=0
Compute the residuals for the associated pressure advection diffusion problem. Used by the Fp precond...
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
Compute the diagonal of the velocity/pressure mass matrices. If which one=0, both are computed,...
////////////////////////////////////////////////////////////////////// //////////////////////////////...
unsigned ntstorage() const
Return the number of doubles required to represent history (one for steady)
virtual double weight(const unsigned &i, const unsigned &j) const
Access function for j-th weight for the i-th derivative.
bool is_steady() const
Flag to indicate if a timestepper has been made steady (possibly temporarily to switch off time-depen...
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
void dshape< 2 >(const double &s, double *DPsi)
Derivatives of 1D shape functions specialised to linear order (2 Nodes)
void shape< 2 >(const double &s, double *Psi)
1D shape functions specialised to linear order (2 Nodes)
std::string to_string(T object, unsigned float_precision=8)
Conversion function that should work for anything with operator<< defined (at least all basic types).
//////////////////////////////////////////////////////////////////// ////////////////////////////////...