30 #ifndef OOMPH_BOUSSINESQ_ELEMENTS_HEADER
31 #define OOMPH_BOUSSINESQ_ELEMENTS_HEADER
64 template<
unsigned DIM>
102 const double&
ra()
const
135 "This function hasn't been implemented for this element",
148 const unsigned& nplot)
const
151 "This function hasn't been implemented for this element",
175 void output(std::ostream& outfile,
const unsigned& nplot)
185 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
191 for (
unsigned i = 0;
i <
DIM;
i++)
197 for (
unsigned i = 0;
i <
DIM;
i++)
208 outfile << std::endl;
222 void output(FILE* file_pt,
const unsigned& n_plot)
229 const unsigned& Nplot,
239 const unsigned& Nplot,
312 for (
unsigned i = 0;
i <
DIM;
i++)
336 #ifdef USE_FD_JACOBIAN_FOR_BUOYANT_Q_ELEMENT
351 Vector<double>& residuals,
359 residuals, jacobian, mass_matrix);
365 #ifdef USE_OFF_DIAGONAL_FD_JACOBIAN_FOR_BUOYANT_Q_ELEMENT
369 void fill_in_off_diagonal_jacobian_blocks_by_fd(
374 unsigned u_nodal_nst[
DIM];
375 for (
unsigned i = 0;
i <
DIM;
i++)
384 unsigned n_dof = this->
ndof();
387 Vector<double> newres(n_dof);
390 int local_unknown = 0, local_eqn = 0;
396 unsigned n_node = this->
nnode();
402 for (
unsigned n = 0;
n < n_node;
n++)
405 for (
unsigned i = 0;
i <
DIM;
i++)
411 if (local_unknown >= 0)
417 double old_var = *value_pt;
420 *value_pt += fd_step;
426 for (
unsigned m = 0;
m < n_dof;
m++)
437 for (
unsigned m = 0;
m < n_node;
m++)
446 (newres[local_eqn] - residuals[local_eqn]) / fd_step;
447 jacobian(local_eqn, local_unknown) = sum;
464 if (local_unknown >= 0)
470 double old_var = *value_pt;
473 *value_pt += fd_step;
479 for (
unsigned m = 0;
m < n_dof;
m++)
489 for (
unsigned m = 0;
m < n_node;
m++)
492 for (
unsigned j = 0;
j <
DIM;
j++)
499 (newres[local_eqn] - residuals[local_eqn]) / fd_step;
500 jacobian(local_eqn, local_unknown) = sum;
527 residuals, jacobian);
531 fill_in_off_diagonal_jacobian_blocks_by_fd(residuals, jacobian);
538 Vector<double>& residuals,
554 fill_in_off_diagonal_jacobian_blocks_by_fd(residuals, jacobian);
574 unsigned u_nodal_nst[
DIM];
575 for (
unsigned i = 0;
i <
DIM;
i++)
584 const unsigned n_node = this->
nnode();
587 Shape psif(n_node), testf(n_node);
594 double Ra = this->
ra();
595 double Pe = this->
pe();
599 int local_eqn = 0, local_unknown = 0;
602 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
609 ipt, psif, dpsifdx, testf, dtestfdx);
619 for (
unsigned l = 0; l < n_node; l++)
624 for (
unsigned j = 0;
j <
DIM;
j++)
626 interpolated_du_adv_diff_dx[
j] += u_value * dpsifdx(l,
j);
633 for (
unsigned l = 0; l < n_node; l++)
640 for (
unsigned i = 0;
i <
DIM;
i++)
647 for (
unsigned l2 = 0; l2 < n_node; l2++)
652 if (local_unknown >= 0)
655 jacobian(local_eqn, local_unknown) +=
670 for (
unsigned l2 = 0; l2 < n_node; l2++)
673 for (
unsigned i2 = 0; i2 <
DIM; i2++)
678 if (local_unknown >= 0)
681 jacobian(local_eqn, local_unknown) -=
682 Pe * psif(l2) * interpolated_du_adv_diff_dx[i2] *
707 residuals, jacobian);
754 template<
unsigned int DIM>
756 :
public virtual QElement<DIM - 1, 3>
799 template<
unsigned DIM>
832 const double&
ra()
const
867 "This function hasn't been implemented for this element",
880 const unsigned& nplot)
const
883 "This function hasn't been implemented for this element",
905 void output(std::ostream& outfile,
const unsigned& nplot)
915 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
921 for (
unsigned i = 0;
i <
DIM;
i++)
927 for (
unsigned i = 0;
i <
DIM;
i++)
938 outfile << std::endl;
951 void output(FILE* file_pt,
const unsigned& n_plot)
958 const unsigned& Nplot,
968 const unsigned& Nplot,
1022 s, advection_values);
1025 for (
unsigned i = 0;
i <
DIM;
i++)
1027 values.push_back(nst_values[
i]);
1031 values.push_back(advection_values[0]);
1055 s, advection_values);
1058 for (
unsigned i = 0;
i <
DIM;
i++)
1060 values.push_back(nst_values[
i]);
1064 values.push_back(advection_values[0]);
1104 this->Ra_pt = cast_father_element_pt->
ra_pt();
1129 unsigned n_fluid_flux =
1136 unsigned n_temp_flux =
1144 for (
unsigned i = 0;
i < n_temp_flux;
i++)
1146 flux[n_fluid_flux +
i] = temp_flux[
i];
1163 unsigned n_fluid_flux =
1166 unsigned n_temp_flux =
1172 for (
unsigned i = 0;
i < n_fluid_flux;
i++)
1178 for (
unsigned i = 0;
i < n_temp_flux;
i++)
1180 flux_index[n_fluid_flux +
i] = 1;
1232 const unsigned& ipt,
1245 for (
unsigned i = 0;
i <
DIM;
i++)
1269 #ifdef USE_FD_JACOBIAN_FOR_REFINEABLE_BUOYANT_Q_ELEMENT
1275 residuals, jacobian);
1296 residuals, jacobian, mass_matrix);
1309 HangInfo *hang_info_pt = 0, *hang_info2_pt = 0;
1313 unsigned u_nodal_nst[
DIM];
1314 for (
unsigned i = 0;
i <
DIM;
i++)
1323 const unsigned n_node = this->
nnode();
1326 Shape psif(n_node), testf(n_node);
1333 double Ra = this->
ra();
1334 double Pe = this->
pe();
1338 int local_eqn = 0, local_unknown = 0;
1341 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
1348 ipt, psif, dpsifdx, testf, dtestfdx);
1358 for (
unsigned l = 0; l < n_node; l++)
1361 double u_value = this->
nodal_value(l, u_nodal_adv_diff);
1363 for (
unsigned j = 0;
j <
DIM;
j++)
1365 interpolated_du_adv_diff_dx[
j] += u_value * dpsifdx(l,
j);
1373 for (
unsigned l = 0; l < n_node; l++)
1378 unsigned n_master = 1;
1379 double hang_weight = 1.0;
1385 if (is_node_hanging)
1388 n_master = hang_info_pt->
nmaster();
1397 for (
unsigned m = 0;
m < n_master;
m++)
1400 if (is_node_hanging)
1416 for (
unsigned i = 0;
i <
DIM;
i++)
1419 if (is_node_hanging)
1435 unsigned n_master2 = 1;
1436 double hang_weight2 = 1.0;
1439 for (
unsigned l2 = 0; l2 < n_node; l2++)
1445 if (is_node2_hanging)
1448 n_master2 = hang_info2_pt->nmaster();
1457 for (
unsigned m2 = 0;
m2 < n_master2;
m2++)
1459 if (is_node2_hanging)
1463 hang_info2_pt->master_node_pt(
m2), u_nodal_adv_diff);
1465 hang_weight2 = hang_info2_pt->master_weight(
m2);
1476 if (local_unknown >= 0)
1479 jacobian(local_eqn, local_unknown) +=
1481 hang_weight * hang_weight2;
1493 if (is_node_hanging)
1510 unsigned n_master2 = 1;
1511 double hang_weight2 = 1.0;
1514 for (
unsigned l2 = 0; l2 < n_node; l2++)
1520 if (is_node2_hanging)
1523 n_master2 = hang_info2_pt->nmaster();
1532 for (
unsigned m2 = 0;
m2 < n_master2;
m2++)
1535 if (is_node2_hanging)
1538 hang_weight2 = hang_info2_pt->master_weight(
m2);
1548 for (
unsigned i2 = 0; i2 <
DIM; i2++)
1551 if (is_node2_hanging)
1555 hang_info2_pt->master_node_pt(
m2), u_nodal_nst[i2]);
1565 if (local_unknown >= 0)
1568 jacobian(local_eqn, local_unknown) -=
1569 Pe * psif(l2) * interpolated_du_adv_diff_dx[i2] *
1570 testf(l) *
W * hang_weight * hang_weight2;
Matrix< Scalar, Dynamic, Dynamic > DenseMatrix
Definition: BenchSparseUtil.h:23
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
RowVector3d w
Definition: Matrix_resize_int.cpp:3
Definition: 3d_ref_b_convect.cc:533
Definition: advection_diffusion_elements.h:52
AdvectionDiffusionEquations()
Definition: advection_diffusion_elements.h:68
double interpolated_u_adv_diff(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
Definition: advection_diffusion_elements.h:458
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: advection_diffusion_elements.h:435
void enable_ALE()
Definition: advection_diffusion_elements.h:135
void disable_ALE()
Definition: advection_diffusion_elements.h:125
const double & pe() const
Peclet number.
Definition: advection_diffusion_elements.h:318
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector (wrapper)
Definition: advection_diffusion_elements.h:421
Definition: boussinesq_elements.h:68
double * Ra_pt
Pointer to a private data member, the Rayleigh number.
Definition: boussinesq_elements.h:71
std::string scalar_name_paraview(const unsigned &i) const
Definition: boussinesq_elements.h:160
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Definition: boussinesq_elements.h:272
static double Default_Physical_Constant_Value
The static default value of the Rayleigh number.
Definition: boussinesq_elements.h:74
void scalar_value_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
Definition: boussinesq_elements.h:146
void output(std::ostream &outfile)
Overload the standard output function with the broken default.
Definition: boussinesq_elements.h:167
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function: Broken default.
Definition: boussinesq_elements.h:222
unsigned nscalar_paraview() const
Definition: boussinesq_elements.h:132
void disable_ALE()
Final override for disable ALE.
Definition: boussinesq_elements.h:114
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Definition: boussinesq_elements.h:258
const double & ra() const
Access function for the Rayleigh number (const version)
Definition: boussinesq_elements.h:102
unsigned required_nvalue(const unsigned &n) const
Definition: boussinesq_elements.h:95
double *& ra_pt()
Access function for the pointer to the Rayleigh number.
Definition: boussinesq_elements.h:108
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: boussinesq_elements.h:696
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
Definition: boussinesq_elements.h:716
void output(std::ostream &outfile, const unsigned &nplot)
Definition: boussinesq_elements.h:175
void unfix_pressure(const unsigned &p_dof)
Unpin p_dof-th pressure dof.
Definition: boussinesq_elements.h:87
BuoyantQCrouzeixRaviartElement()
Definition: boussinesq_elements.h:80
void output(FILE *file_pt)
C-style output function: Broken default.
Definition: boussinesq_elements.h:216
void fill_in_off_diagonal_jacobian_blocks_analytic(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: boussinesq_elements.h:564
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Definition: boussinesq_elements.h:323
unsigned u_index_adv_diff() const
Definition: boussinesq_elements.h:248
void output_fct(std::ostream &outfile, const unsigned &Nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output function for an exact solution: Broken default.
Definition: boussinesq_elements.h:228
void get_wind_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Definition: boussinesq_elements.h:283
void output_fct(std::ostream &outfile, const unsigned &Nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: boussinesq_elements.h:238
void get_body_force_nst(const double &time, const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &result)
Definition: boussinesq_elements.h:298
void enable_ALE()
Final override for enable ALE.
Definition: boussinesq_elements.h:122
void unpin(const unsigned &i)
Unpin the i-th stored variable.
Definition: nodes.h:391
double * value_pt(const unsigned &i) const
Definition: nodes.h:324
FaceGeometry()
Definition: boussinesq_elements.h:759
FaceGeometry()
Definition: boussinesq_elements.h:770
Definition: elements.h:4998
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: elements.h:1735
virtual void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output an exact solution over the element.
Definition: elements.h:3104
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
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Definition: elements.h:3161
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
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
virtual void compute_error(FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Calculate the norm of the error and that of the exact solution.
Definition: elements.h:3198
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Definition: elements.h:1759
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Definition: elements.h:3148
virtual unsigned nplot_points(const unsigned &nplot) const
Definition: elements.h:3186
double raw_nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2576
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Definition: elements.h:3174
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Definition: elements.h:1765
static double Default_fd_jacobian_step
Definition: elements.h:1198
unsigned ndof() const
Return the number of equations/dofs in the element.
Definition: elements.h:835
Data *& internal_data_pt(const unsigned &i)
Return a pointer to i-th internal data object.
Definition: elements.h:622
virtual void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
Definition: elements.cc:1322
double const & master_weight(const unsigned &i) const
Return weight for dofs on i-th master node.
Definition: nodes.h:808
Node *const & master_node_pt(const unsigned &i) const
Return a pointer to the i-th master node.
Definition: nodes.h:791
unsigned nmaster() const
Return the number of master nodes.
Definition: nodes.h:785
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
Definition: navier_stokes_elements.h:395
void disable_ALE()
Definition: navier_stokes_elements.h:909
virtual unsigned u_index_nst(const unsigned &i) const
Definition: navier_stokes_elements.h:866
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
void enable_ALE()
Definition: navier_stokes_elements.h:918
virtual double interpolated_p_nst(const Vector< double > &s) const
Return FE interpolated pressure at local coordinate s.
Definition: navier_stokes_elements.h:1639
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: navier_stokes_elements.h:1273
NavierStokesEquations()
Definition: navier_stokes_elements.h:677
const Vector< double > & g() const
Vector of gravitational components.
Definition: navier_stokes_elements.h:765
HangInfo *const & hanging_pt() const
Definition: nodes.h:1228
bool is_hanging() const
Test whether the node is geometrically hanging.
Definition: nodes.h:1285
Definition: oomph_definitions.h:222
Definition: elements.h:3439
Definition: advection_diffusion_elements.h:610
Definition: navier_stokes_elements.h:1749
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
unsigned P_nst_internal_index
Definition: navier_stokes_elements.h:1757
Definition: Qelements.h:459
Definition: refineable_advection_diffusion_elements.h:58
void further_build()
Further build: Copy source function pointer from father element.
Definition: refineable_advection_diffusion_elements.h:177
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
Definition: refineable_advection_diffusion_elements.h:77
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Definition: refineable_advection_diffusion_elements.h:84
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: refineable_advection_diffusion_elements.h:94
Definition: boussinesq_elements.h:803
void rebuild_from_sons(Mesh *&mesh_pt)
Definition: boussinesq_elements.h:1081
void fill_in_contribution_to_jacobian_and_mass_matrix(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
Definition: boussinesq_elements.h:1289
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Definition: boussinesq_elements.h:1205
void get_wind_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Definition: boussinesq_elements.h:1216
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: boussinesq_elements.h:1266
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the constituent elements' contribution to the residual vector.
Definition: boussinesq_elements.h:1252
RefineableBuoyantQCrouzeixRaviartElement()
Definition: boussinesq_elements.h:815
unsigned nrecovery_order()
The recovery order is that of the NavierStokes elements.
Definition: boussinesq_elements.h:1109
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Definition: boussinesq_elements.h:1126
unsigned nscalar_paraview() const
Definition: boussinesq_elements.h:864
static double Default_Physical_Constant_Value
The static default value of the Rayleigh number.
Definition: boussinesq_elements.h:809
void output(std::ostream &outfile, const unsigned &nplot)
Definition: boussinesq_elements.h:905
unsigned ncompound_fluxes()
Definition: boussinesq_elements.h:1153
Node * vertex_node_pt(const unsigned &j) const
Definition: boussinesq_elements.h:991
const double & ra() const
Access function for the Rayleigh number (const version)
Definition: boussinesq_elements.h:832
void further_build()
Definition: boussinesq_elements.h:1091
unsigned num_Z2_flux_terms()
Definition: boussinesq_elements.h:1116
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function: Broken default.
Definition: boussinesq_elements.h:951
double * Ra_pt
Pointer to a new physical variable, the Rayleigh number.
Definition: boussinesq_elements.h:806
unsigned required_nvalue(const unsigned &n) const
Definition: boussinesq_elements.h:825
void get_body_force_nst(const double &time, const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &result)
Definition: boussinesq_elements.h:1231
void fill_in_off_diagonal_jacobian_blocks_analytic(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: boussinesq_elements.h:1301
void enable_ALE()
Final override for enable ALE.
Definition: boussinesq_elements.h:853
void disable_ALE()
Final override for disable ALE.
Definition: boussinesq_elements.h:845
double *& ra_pt()
Access function for the pointer to the Rayleigh number.
Definition: boussinesq_elements.h:838
unsigned nvertex_node() const
Definition: boussinesq_elements.h:984
void output(FILE *file_pt)
C-style output function: Broken default.
Definition: boussinesq_elements.h:945
unsigned ncont_interpolated_values() const
Definition: boussinesq_elements.h:998
void scalar_value_paraview(std::ofstream &file_out, const unsigned &i, const unsigned &nplot) const
Definition: boussinesq_elements.h:878
void output_fct(std::ostream &outfile, const unsigned &Nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output function for an exact solution: Broken default.
Definition: boussinesq_elements.h:957
void output_fct(std::ostream &outfile, const unsigned &Nplot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: boussinesq_elements.h:967
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: boussinesq_elements.h:1007
void further_setup_hanging_nodes()
Definition: boussinesq_elements.h:1071
void output(std::ostream &outfile)
Overload the standard output function with the broken default.
Definition: boussinesq_elements.h:898
unsigned u_index_adv_diff() const
Definition: boussinesq_elements.h:977
void get_Z2_compound_flux_indices(Vector< unsigned > &flux_index)
Definition: boussinesq_elements.h:1160
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Definition: boussinesq_elements.h:1191
std::string scalar_name_paraview(const unsigned &i) const
Definition: boussinesq_elements.h:892
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Definition: boussinesq_elements.h:1039
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
Definition: refineable_elements.h:539
int local_hang_eqn(Node *const &node_pt, const unsigned &i)
Definition: refineable_elements.h:278
unsigned num_Z2_flux_terms()
Number of 'flux' terms for Z2 error estimation.
Definition: refineable_navier_stokes_elements.h:395
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Definition: refineable_navier_stokes_elements.h:403
Definition: refineable_advection_diffusion_elements.h:355
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
Definition: refineable_advection_diffusion_elements.h:395
Refineable version of Crouzeix Raviart elements. Generic class definitions.
Definition: refineable_navier_stokes_elements.h:1109
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: refineable_navier_stokes_elements.h:1178
unsigned nrecovery_order()
Definition: refineable_navier_stokes_elements.h:1157
void further_setup_hanging_nodes()
Definition: refineable_navier_stokes_elements.h:1234
void rebuild_from_sons(Mesh *&mesh_pt)
RealScalar s
Definition: level1_cplx_impl.h:130
int * m
Definition: level2_cplx_impl.h:294
#define DIM
Definition: linearised_navier_stokes_elements.h:44
double Pe
Peclet number.
Definition: axisym_heat_sphere.cc:59
double Ra
Rayleigh number.
Definition: axisym_heat_sphere.cc:62
void gravity(const double &t, const Vector< double > &xi, Vector< double > &b)
Definition: ConstraintElementsUnitTest.cpp:20
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
int error
Definition: calibrate.py:297
Definition: MortaringCantileverCompareToNonMortaring.cpp:176
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286
double Default_Physical_Constant_Value
Default value for physical constants.
Definition: multi_domain_boussinesq_elements.h:48
@ 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
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2