30 #ifndef OOMPH_MULTI_DOMAIN_BOUSSINESQ_ELEMENTS_HEADER
31 #define OOMPH_MULTI_DOMAIN_BOUSSINESQ_ELEMENTS_HEADER
36 #include "advection_diffusion.h"
37 #include "navier_stokes.h"
45 namespace MultiDomainBoussinesqHelper
62 template<
class NST_ELEMENT,
class AD_ELEMENT>
64 :
public virtual NST_ELEMENT,
65 public virtual ElementWithExternalElement
72 : NST_ELEMENT(), ElementWithExternalElement()
78 this->set_ninteraction(1);
82 const double&
ra()
const
99 NST_ELEMENT::further_build();
104 cast_father_element_pt =
dynamic_cast<
106 this->father_element_pt());
110 this->
Ra_pt = cast_father_element_pt->
ra_pt();
113 if (!cast_father_element_pt->external_geometric_data_is_included())
115 this->ignore_external_geometric_data();
125 const Vector<double>& s,
126 const Vector<double>& x,
127 Vector<double>& body_force)
130 const unsigned interaction = 0;
134 const AD_ELEMENT* adv_diff_el_pt =
135 dynamic_cast<AD_ELEMENT*
>(external_element_pt(interaction, ipt));
138 const double interpolated_t = adv_diff_el_pt->interpolated_u_adv_diff(
139 external_element_local_coord(interaction, ipt));
143 Vector<double> gravity(NST_ELEMENT::g());
146 const unsigned n_dim = this->dim();
147 for (
unsigned i = 0; i < n_dim; i++)
149 body_force[i] = -gravity[i] * interpolated_t *
ra();
157 DenseMatrix<double>& jacobian)
160 NST_ELEMENT::fill_in_contribution_to_jacobian(residuals, jacobian);
162 #ifdef USE_FD_FOR_DERIVATIVES_WRT_EXTERNAL_DATA_IN_MULTI_DOMAIN_BOUSSINESQ
165 this->fill_in_jacobian_from_external_interaction_by_fd(residuals,
180 Vector<double>& residuals,
181 DenseMatrix<double>& jacobian,
182 DenseMatrix<double>& mass_matrix)
187 FiniteElement::fill_in_contribution_to_jacobian_and_mass_matrix(
188 residuals, jacobian, mass_matrix);
196 DenseMatrix<double>& result,
197 Vector<unsigned>& global_eqn_number);
203 DenseMatrix<double>& jacobian)
208 const unsigned n_dim = this->dim();
209 unsigned u_nodal_nst[n_dim];
210 for (
unsigned i = 0; i < n_dim; i++)
212 u_nodal_nst[i] = this->u_index_nst(i);
216 const unsigned n_node = this->nnode();
219 Shape psif(n_node), testf(n_node);
220 DShape dpsifdx(n_node, n_dim), dtestfdx(n_node, n_dim);
223 const unsigned n_intpt = this->integral_pt()->nweight();
226 int local_eqn = 0, local_unknown = 0;
229 HangInfo* hang_info_pt = 0;
232 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
235 double w = this->integral_pt()->weight(ipt);
238 double J = this->dshape_and_dtest_eulerian_at_knot_nst(
239 ipt, psif, dpsifdx, testf, dtestfdx);
248 DenseMatrix<double> dbody_dexternal_element_data;
252 Vector<unsigned> global_eqn_number_of_external_element_data;
257 dbody_dexternal_element_data,
258 global_eqn_number_of_external_element_data);
261 const unsigned n_external_element_data =
262 global_eqn_number_of_external_element_data.size();
265 for (
unsigned l = 0; l < n_node; l++)
270 unsigned n_master = 1;
271 double hang_weight = 1.0;
274 bool is_node_hanging = this->node_pt(l)->is_hanging();
279 hang_info_pt = this->node_pt(l)->hanging_pt();
280 n_master = hang_info_pt->nmaster();
289 for (
unsigned m = 0; m < n_master; m++)
295 hang_weight = hang_info_pt->master_weight(m);
305 for (
unsigned i = 0; i < n_dim; i++)
311 local_eqn = this->local_hang_eqn(
312 hang_info_pt->master_node_pt(m), u_nodal_nst[i]);
317 local_eqn = this->nodal_local_eqn(l, u_nodal_nst[i]);
323 for (
unsigned l2 = 0; l2 < n_external_element_data; l2++)
327 local_unknown = this->local_eqn_number(
328 global_eqn_number_of_external_element_data[l2]);
329 if (local_unknown >= 0)
332 jacobian(local_eqn, local_unknown) +=
333 dbody_dexternal_element_data(i, l2) * testf(l) *
346 std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list)
const
349 NST_ELEMENT::get_dof_numbers_for_unknowns(dof_lookup_list);
355 return NST_ELEMENT::ndof_types();
374 template<
class AD_ELEMENT,
class NST_ELEMENT>
376 :
public virtual AD_ELEMENT,
377 public virtual ElementWithExternalElement
382 : AD_ELEMENT(), ElementWithExternalElement()
385 this->set_ninteraction(1);
400 void output(std::ostream& outfile,
const unsigned& nplot)
403 unsigned n_dim = this->dim();
404 Vector<double> s(n_dim);
407 outfile << this->tecplot_zone_string(nplot);
410 unsigned num_plot_points = this->nplot_points(nplot);
411 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
414 this->get_s_plot(iplot, nplot, s);
417 for (
unsigned i = 0; i < n_dim; i++)
419 outfile << this->interpolated_x(s, i) <<
" ";
423 outfile << AD_ELEMENT::interpolated_u_adv_diff(s) << std::endl;
425 outfile << std::endl;
428 this->write_tecplot_zone_footer(outfile, nplot);
436 FiniteElement::output(outfile);
442 FiniteElement::output(file_pt);
446 void output(FILE* file_pt,
const unsigned& n_plot)
448 FiniteElement::output(file_pt, n_plot);
457 const Vector<double>& s,
458 const Vector<double>& x,
459 Vector<double>& wind)
const
462 unsigned interaction = 0;
465 NST_ELEMENT* nst_el_pt =
466 dynamic_cast<NST_ELEMENT*
>(external_element_pt(interaction, ipt));
469 nst_el_pt->interpolated_u_nst(
470 external_element_local_coord(interaction, ipt), wind);
477 DenseMatrix<double>& jacobian)
480 AD_ELEMENT::fill_in_contribution_to_jacobian(residuals, jacobian);
482 #ifdef USE_FD_FOR_DERIVATIVES_WRT_EXTERNAL_DATA_IN_MULTI_DOMAIN_BOUSSINESQ
485 this->fill_in_jacobian_from_external_interaction_by_fd(residuals,
502 Vector<std::set<FiniteElement*>>
const& external_elements_pt,
503 std::set<std::pair<Data*, unsigned>>& paired_interaction_data);
508 Vector<double>& residuals,
509 DenseMatrix<double>& jacobian,
510 DenseMatrix<double>& mass_matrix)
515 FiniteElement::fill_in_contribution_to_jacobian_and_mass_matrix(
516 residuals, jacobian, mass_matrix);
525 Vector<double>& result,
526 Vector<unsigned>& global_eqn_number)
529 unsigned interaction = 0;
532 NST_ELEMENT* source_el_pt =
533 dynamic_cast<NST_ELEMENT*
>(external_element_pt(interaction, ipt));
538 source_el_pt->dinterpolated_u_nst_ddata(
539 external_element_local_coord(interaction, ipt),
549 DenseMatrix<double>& jacobian)
553 const unsigned u_nodal_adv_diff = this->u_index_adv_diff();
556 const unsigned n_node = this->nnode();
559 const unsigned n_dim = this->dim();
562 Shape psi(n_node), test(n_node);
563 DShape dpsidx(n_node, n_dim), dtestdx(n_node, n_dim);
566 const unsigned n_intpt = this->integral_pt()->nweight();
569 int local_eqn = 0, local_unknown = 0;
572 HangInfo* hang_info_pt = 0;
575 const double peclet = this->pe();
578 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
581 double w = this->integral_pt()->weight(ipt);
584 double J = this->dshape_and_dtest_eulerian_at_knot_adv_diff(
585 ipt, psi, dpsidx, test, dtestdx);
591 Vector<double> interpolated_dudx(n_dim, 0.0);
593 for (
unsigned l = 0; l < n_node; l++)
596 for (
unsigned j = 0; j < n_dim; j++)
598 interpolated_dudx[j] +=
599 this->nodal_value(l, u_nodal_adv_diff) * dpsidx(l, j);
605 Vector<double> dwind_dexternal_element_data;
609 Vector<unsigned> global_eqn_number_of_external_element_data;
612 for (
unsigned i2 = 0; i2 < n_dim; i2++)
618 dwind_dexternal_element_data,
619 global_eqn_number_of_external_element_data);
623 const unsigned n_external_element_data =
624 global_eqn_number_of_external_element_data.size();
627 for (
unsigned l = 0; l < n_node; l++)
631 unsigned n_master = 1;
632 double hang_weight = 1.0;
635 bool is_node_hanging = this->node_pt(l)->is_hanging();
640 hang_info_pt = this->node_pt(l)->hanging_pt();
641 n_master = hang_info_pt->nmaster();
650 for (
unsigned m = 0; m < n_master; m++)
656 hang_weight = hang_info_pt->master_weight(m);
668 local_eqn = this->local_hang_eqn(
669 hang_info_pt->master_node_pt(m), u_nodal_adv_diff);
674 local_eqn = this->nodal_local_eqn(l, u_nodal_adv_diff);
680 for (
unsigned l2 = 0; l2 < n_external_element_data; l2++)
684 local_unknown = this->local_eqn_number(
685 global_eqn_number_of_external_element_data[l2]);
686 if (local_unknown >= 0)
689 jacobian(local_eqn, local_unknown) -=
690 peclet * dwind_dexternal_element_data[l2] *
691 interpolated_dudx[i2] * test(l) * hang_weight * W;
703 std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list)
const
706 unsigned n_node = this->nnode();
709 std::pair<unsigned, unsigned> dof_lookup;
712 for (
unsigned n = 0; n < n_node; n++)
715 unsigned nv = this->node_pt(n)->nvalue();
718 for (
unsigned v = 0; v < nv; v++)
721 int local_eqn_number = this->nodal_local_eqn(n, v);
726 if (local_eqn_number >= 0)
730 dof_lookup.first = this->eqn_number(local_eqn_number);
733 dof_lookup.second = 0;
736 dof_lookup_list.push_front(dof_lookup);
762 template<
class AD_ELEMENT,
class NST_ELEMENT>
765 Vector<std::set<FiniteElement*>>
const& external_elements_pt,
766 std::set<std::pair<Data*, unsigned>>& paired_interaction_data)
769 const unsigned interaction = 0;
773 for (std::set<FiniteElement*>::iterator it =
774 external_elements_pt[interaction].begin();
775 it != external_elements_pt[interaction].end();
779 NST_ELEMENT* external_fluid_el_pt =
dynamic_cast<NST_ELEMENT*
>(*it);
782 unsigned nnod = external_fluid_el_pt->nnode();
783 for (
unsigned j = 0; j < nnod; j++)
786 Data* veloc_data_pt = external_fluid_el_pt->node_pt(j);
789 const unsigned n_dim = this->dim();
790 for (
unsigned i = 0; i < n_dim; i++)
793 unsigned val = external_fluid_el_pt->u_index_nst(i);
797 paired_interaction_data.insert(std::make_pair(veloc_data_pt, val));
814 template<
class NST_ELEMENT,
class AD_ELEMENT>
816 :
public virtual NST_ELEMENT,
817 public virtual ElementWithExternalElement
828 : NST_ELEMENT(), ElementWithExternalElement()
833 this->set_ninteraction(1);
837 const double&
ra()
const
853 const Vector<double>& s,
854 const Vector<double>& x,
855 Vector<double>& result);
861 DenseMatrix<double>& result,
862 Vector<unsigned>& global_eqn_number);
868 DenseMatrix<double>& jacobian)
870 #ifdef USE_FD_JACOBIAN_NST_IN_MULTI_DOMAIN_BOUSSINESQ
873 ElementWithExternalElement::fill_in_contribution_to_jacobian(residuals,
879 NST_ELEMENT::fill_in_contribution_to_jacobian(residuals, jacobian);
890 Vector<double>& residuals,
891 DenseMatrix<double>& jacobian,
892 DenseMatrix<double>& mass_matrix)
897 FiniteElement::fill_in_contribution_to_jacobian_and_mass_matrix(
898 residuals, jacobian, mass_matrix);
904 DenseMatrix<double>& jacobian)
908 const unsigned n_dim = this->dim();
909 Vector<unsigned> u_nodal_nst(n_dim);
910 for (
unsigned i = 0; i < n_dim; i++)
912 u_nodal_nst[i] = this->u_index_nst(i);
916 const unsigned n_node = this->nnode();
919 Shape psif(n_node), testf(n_node);
920 DShape dpsifdx(n_node, n_dim), dtestfdx(n_node, n_dim);
923 const unsigned n_intpt = this->integral_pt()->nweight();
926 int local_eqn = 0, local_unknown = 0;
929 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
932 double w = this->integral_pt()->weight(ipt);
935 double J = this->dshape_and_dtest_eulerian_at_knot_nst(
936 ipt, psif, dpsifdx, testf, dtestfdx);
945 DenseMatrix<double> dbody_dexternal_element_data;
949 Vector<unsigned> global_eqn_number_of_external_element_data;
954 dbody_dexternal_element_data,
955 global_eqn_number_of_external_element_data);
958 const unsigned n_external_element_data =
959 global_eqn_number_of_external_element_data.size();
962 for (
unsigned l = 0; l < n_node; l++)
969 for (
unsigned i = 0; i < n_dim; i++)
972 local_eqn = this->nodal_local_eqn(l, u_nodal_nst[i]);
976 for (
unsigned l2 = 0; l2 < n_external_element_data; l2++)
980 local_unknown = this->local_eqn_number(
981 global_eqn_number_of_external_element_data[l2]);
982 if (local_unknown >= 0)
985 jacobian(local_eqn, local_unknown) +=
986 dbody_dexternal_element_data(i, l2) * testf(l) * W;
1001 template<
class NST_ELEMENT,
class AD_ELEMENT>
1004 const unsigned& ipt,
1005 const Vector<double>& s,
1006 const Vector<double>& x,
1007 Vector<double>& result)
1010 unsigned interaction = 0;
1013 const double interpolated_t =
1014 dynamic_cast<AD_ELEMENT*
>(external_element_pt(interaction, ipt))
1015 ->interpolated_u_adv_diff(
1016 external_element_local_coord(interaction, ipt));
1020 Vector<double> gravity(NST_ELEMENT::g());
1023 const unsigned n_dim = this->dim();
1024 for (
unsigned i = 0; i < n_dim; i++)
1026 result[i] = -gravity[i] * interpolated_t * ra();
1032 template<
class NST_ELEMENT,
class AD_ELEMENT>
1034 :
public virtual FaceGeometry<NST_ELEMENT>
1045 template<
class NST_ELEMENT,
class AD_ELEMENT>
1048 :
public virtual FaceGeometry<FaceGeometry<NST_ELEMENT>>
1069 template<
class AD_ELEMENT,
class NST_ELEMENT>
1071 :
public virtual AD_ELEMENT,
1072 public virtual ElementWithExternalElement
1077 : AD_ELEMENT(), ElementWithExternalElement()
1080 this->set_ninteraction(1);
1088 const Vector<double>& s,
1089 const Vector<double>& x,
1090 Vector<double>& wind)
const;
1104 void output(std::ostream& outfile,
const unsigned& nplot)
1107 const unsigned n_dim = this->dim();
1110 Vector<double> s(n_dim);
1113 outfile << this->tecplot_zone_string(nplot);
1116 unsigned num_plot_points = this->nplot_points(nplot);
1117 for (
unsigned iplot = 0; iplot < num_plot_points; iplot++)
1120 this->get_s_plot(iplot, nplot, s);
1123 for (
unsigned i = 0; i < n_dim; i++)
1125 outfile << this->interpolated_x(s, i) <<
" ";
1129 outfile << AD_ELEMENT::interpolated_u_adv_diff(s) << std::endl;
1131 outfile << std::endl;
1134 this->write_tecplot_zone_footer(outfile, nplot);
1142 FiniteElement::output(outfile);
1148 FiniteElement::output(file_pt);
1152 void output(FILE* file_pt,
const unsigned& n_plot)
1154 FiniteElement::output(file_pt, n_plot);
1161 const unsigned& ipt,
1163 Vector<double>& result,
1164 Vector<unsigned>& global_eqn_number);
1170 DenseMatrix<double>& jacobian)
1172 #ifdef USE_FD_JACOBIAN_IN_MULTI_DOMAIN_BOUSSINESQ
1175 ElementWithExternalElement::fill_in_contribution_to_jacobian(residuals,
1181 AD_ELEMENT::fill_in_contribution_to_jacobian(residuals, jacobian);
1192 Vector<double>& residuals,
1193 DenseMatrix<double>& jacobian,
1194 DenseMatrix<double>& mass_matrix)
1199 FiniteElement::fill_in_contribution_to_jacobian_and_mass_matrix(
1200 residuals, jacobian, mass_matrix);
1207 DenseMatrix<double>& jacobian)
1210 const unsigned u_nodal_adv_diff = this->u_index_adv_diff();
1213 const unsigned n_node = this->nnode();
1216 const unsigned n_dim = this->dim();
1219 Shape psi(n_node), test(n_node);
1220 DShape dpsidx(n_node, n_dim), dtestdx(n_node, n_dim);
1223 const unsigned n_intpt = this->integral_pt()->nweight();
1226 int local_eqn = 0, local_unknown = 0;
1229 const double peclet = this->pe();
1232 for (
unsigned ipt = 0; ipt < n_intpt; ipt++)
1235 double w = this->integral_pt()->weight(ipt);
1238 double J = this->dshape_and_dtest_eulerian_at_knot_adv_diff(
1239 ipt, psi, dpsidx, test, dtestdx);
1245 Vector<double> interpolated_dudx(n_dim, 0.0);
1247 for (
unsigned l = 0; l < n_node; l++)
1250 for (
unsigned j = 0; j < n_dim; j++)
1252 interpolated_dudx[j] +=
1253 this->raw_nodal_value(l, u_nodal_adv_diff) * dpsidx(l, j);
1259 Vector<double> dwind_dexternal_element_data;
1263 Vector<unsigned> global_eqn_number_of_external_element_data;
1267 for (
unsigned i2 = 0; i2 < n_dim; i2++)
1273 dwind_dexternal_element_data,
1274 global_eqn_number_of_external_element_data);
1277 const unsigned n_external_element_data =
1278 global_eqn_number_of_external_element_data.size();
1281 for (
unsigned l = 0; l < n_node; l++)
1287 local_eqn = this->nodal_local_eqn(l, u_nodal_adv_diff);
1291 for (
unsigned l2 = 0; l2 < n_external_element_data; l2++)
1295 local_unknown = this->local_eqn_number(
1296 global_eqn_number_of_external_element_data[l2]);
1297 if (local_unknown >= 0)
1300 jacobian(local_eqn, local_unknown) -=
1301 peclet * dwind_dexternal_element_data[l2] *
1302 interpolated_dudx[i2] * test(l) * W;
1319 template<
class AD_ELEMENT,
class NST_ELEMENT>
1322 const Vector<double>& s,
1323 const Vector<double>& x,
1324 Vector<double>& wind)
const
1327 unsigned interaction = 0;
1330 NST_ELEMENT* source_el_pt =
1331 dynamic_cast<NST_ELEMENT*
>(external_element_pt(interaction, ipt));
1334 source_el_pt->interpolated_u_nst(
1335 external_element_local_coord(interaction, ipt), wind);
1342 template<
class AD_ELEMENT,
class NST_ELEMENT>
1345 const unsigned& ipt,
1347 Vector<double>& result,
1348 Vector<unsigned>& global_eqn_number)
1351 unsigned interaction = 0;
1354 NST_ELEMENT* source_el_pt =
1355 dynamic_cast<NST_ELEMENT*
>(external_element_pt(interaction, ipt));
1360 source_el_pt->dinterpolated_u_nst_ddata(
1361 external_element_local_coord(interaction, ipt),
1375 template<
class NST_ELEMENT,
class AD_ELEMENT>
1378 const unsigned& ipt,
1379 DenseMatrix<double>& result,
1380 Vector<unsigned>& global_eqn_number)
1384 Vector<double> gravity(NST_ELEMENT::g());
1387 unsigned interaction = 0;
1390 Vector<double> du_adv_diff_ddata;
1393 AD_ELEMENT* source_el_pt =
1394 dynamic_cast<AD_ELEMENT*
>(external_element_pt(interaction, ipt));
1396 if (source_el_pt == 0)
1398 throw OomphLibError(
"External element could not be cast to AD_ELEMENT\n",
1399 OOMPH_CURRENT_FUNCTION,
1400 OOMPH_EXCEPTION_LOCATION);
1405 source_el_pt->dinterpolated_u_adv_diff_ddata(
1406 external_element_local_coord(interaction, ipt),
1411 unsigned n_external_element_data = du_adv_diff_ddata.size();
1414 const unsigned n_dim = this->dim();
1415 result.resize(n_dim, n_external_element_data);
1418 for (
unsigned i = 0; i < n_dim; i++)
1421 for (
unsigned n = 0; n < n_external_element_data; n++)
1423 result(i, n) = -gravity[i] * du_adv_diff_ddata[n] * ra();
1433 template<
class NST_ELEMENT,
class AD_ELEMENT>
1436 const unsigned& ipt,
1437 DenseMatrix<double>& result,
1438 Vector<unsigned>& global_eqn_number)
1442 Vector<double> gravity(NST_ELEMENT::g());
1445 unsigned interaction = 0;
1448 Vector<double> du_adv_diff_ddata;
1451 AD_ELEMENT* source_el_pt =
1452 dynamic_cast<AD_ELEMENT*
>(external_element_pt(interaction, ipt));
1454 if (source_el_pt == 0)
1456 throw OomphLibError(
"External element could not be cast to AD_ELEMENT\n",
1457 OOMPH_CURRENT_FUNCTION,
1458 OOMPH_EXCEPTION_LOCATION);
1463 source_el_pt->dinterpolated_u_adv_diff_ddata(
1464 external_element_local_coord(interaction, ipt),
1469 unsigned n_external_element_data = du_adv_diff_ddata.size();
1472 const unsigned n_dim = this->dim();
1473 result.resize(n_dim, n_external_element_data);
1476 for (
unsigned i = 0; i < n_dim; i++)
1479 for (
unsigned n = 0; n < n_external_element_data; n++)
1481 result(i, n) = -gravity[i] * du_adv_diff_ddata[n] * ra();
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
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.
AdvectionDiffusionBoussinesqElement()
Constructor: call the underlying constructors.
void get_wind_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Overload the wind function in the advection-diffusion equations. This provides the coupling from the ...
void get_dwind_adv_diff_dexternal_element_data(const unsigned &ipt, const unsigned &i, Vector< double > &result, Vector< unsigned > &global_eqn_number)
Fill in the derivatives of the wind with respect to the external unknowns.
void fill_in_off_diagonal_block_analytic(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute the contribution of the external degrees of freedom (velocities) on the AdvectionDiffusion eq...
void output(std::ostream &outfile)
Overload the standard output function with the broken default.
void output(FILE *file_pt)
C-style output function: Broken default.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function: Broken default.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute the element's residual vector and the Jacobian matrix. Jacobian is computed by finite-differe...
void output(std::ostream &outfile, const unsigned &nplot)
Output function: Output x, y, theta at Nplot^DIM plot points.
FaceGeometry()
Constructor calls the constructor of the NST_ELEMENT (Only the Intel compiler seems to need this!...
FaceGeometry()
Constructor calls the constructor of the NST_ELEMENT (Only the Intel compiler seems to need this!...
////////////////////////////////////////////////////////////////////// //////////////////////////////...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute the element's residual vector and the Jacobian matrix. Jacobian is computed by finite-differe...
double * Ra_pt
Pointer to a private data member, the Rayleigh number.
void get_body_force_nst(const double &time, const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &result)
Overload get_body_force_nst to get the temperature "body force" from the "source" AdvectionDiffusion ...
const double & ra() const
Access function for the Rayleigh number (const version)
void fill_in_off_diagonal_block_analytic(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute the contribution of the external degrees of freedom (temperatures) on the Navier-Stokes equat...
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.
void get_dbody_force_nst_dexternal_element_data(const unsigned &ipt, DenseMatrix< double > &result, Vector< unsigned > &global_eqn_number)
Fill in the derivatives of the body force with respect to the external unknowns.
NavierStokesBoussinesqElement()
Constructor: call the underlying constructors and initialise the pointer to the Rayleigh number to po...
double *& ra_pt()
Access function for the pointer to the Rayleigh number.
//////////////////////////////////////////////////////////////// ////////////////////////////////////...
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Classify dofs for use in block preconditioner.
void identify_all_field_data_for_external_interaction(Vector< std::set< FiniteElement * >> const &external_elements_pt, std::set< std::pair< Data *, unsigned >> &paired_interaction_data)
Overload the function that must return all field data involved in the interaction with the external (...
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.
RefineableAdvectionDiffusionBoussinesqElement()
Constructor: call the underlying constructors.
void get_dwind_adv_diff_dexternal_element_data(const unsigned &ipt, const unsigned &i, Vector< double > &result, Vector< unsigned > &global_eqn_number)
Fill in the derivatives of the wind with respect to the external unknowns.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function: Broken default.
void get_wind_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Overload the wind function in the advection-diffusion equations. This provides the coupling from the ...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute the element's residual vector and the Jacobian matrix.
void output(std::ostream &outfile, const unsigned &nplot)
Output function: Output x, y, theta at Nplot^DIM plot points.
void fill_in_off_diagonal_block_analytic(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute the contribution of the external degrees of freedom (velocities) on the advection-diffusion e...
void output(FILE *file_pt)
C-style output function: Broken default.
unsigned ndof_types() const
Specify number of dof types for use in block preconditioner.
void output(std::ostream &outfile)
Overload the standard output function with the broken default.
////////////////////////////////////////////////////////////////////// //////////////////////////////...
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.
void get_dbody_force_nst_dexternal_element_data(const unsigned &ipt, DenseMatrix< double > &result, Vector< unsigned > &global_eqn_number)
Fill in the derivatives of the body force with respect to the external unknowns.
double * Ra_pt
Pointer to a private data member, the Rayleigh number.
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Classify dof numbers as in underlying element.
double *& ra_pt()
Access function for the pointer to the Rayleigh number.
void get_body_force_nst(const double &time, const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &body_force)
Overload get_body_force_nst() to return the temperature-dependent buoyancy force, using the temperatu...
void further_build()
Call the underlying single-physics element's further_build() functions and make sure that the pointer...
unsigned ndof_types() const
Get number of dof types from underlying element.
void fill_in_off_diagonal_block_analytic(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute the contribution of the external degrees of freedom (temperatures) on the Navier-Stokes equat...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute the element's residual vector and the Jacobian matrix.
const double & ra() const
Access function for the Rayleigh number (const version)
RefineableNavierStokesBoussinesqElement()
Constructor: call the underlying constructors and initialise the pointer to the Rayleigh number to po...
double Default_Physical_Constant_Value
Default value for physical constants.