37 #include "advection_diffusion.h" 38 #include "navier_stokes.h" 39 #include "multi_physics.h" 42 #include "meshes/rectangular_quadmesh.h" 45 using namespace oomph;
78 template<
class NST_ELEMENT,
class AD_ELEMENT>
103 set_boundary_conditions(time_pt()->time());
108 const double &pvalue)
111 dynamic_cast<NST_ELEMENT*
>(nst_mesh_pt()->element_pt(e))->
112 fix_pressure(pdof,pvalue);
119 void set_boundary_conditions(
const double &time);
124 return dynamic_cast<RectangularQuadMesh<NST_ELEMENT>*
>(Nst_mesh_pt);
130 return dynamic_cast<RectangularQuadMesh<AD_ELEMENT>*
>(Adv_diff_mesh_pt);
151 template<
class NST_ELEMENT,
class AD_ELEMENT>
156 add_time_stepper_pt(
new BDF<2>);
159 Doc_info.set_directory(
"RESLT");
175 new RectangularQuadMesh<NST_ELEMENT>(n_x,n_y,l_x,l_y,time_stepper_pt());
177 new RectangularQuadMesh<AD_ELEMENT>(n_x,n_y,l_x,l_y,time_stepper_pt());
184 unsigned num_bound = nst_mesh_pt()->nboundary();
185 for(
unsigned ibound=0;ibound<num_bound;ibound++)
191 unsigned num_nod= nst_mesh_pt()->nboundary_node(ibound);
192 for (
unsigned inod=0;inod<num_nod;inod++)
196 if ((ibound==1) || (ibound==3))
202 for(
unsigned j=0;j<val_max;j++)
204 nst_mesh_pt()->boundary_node_pt(ibound,inod)->pin(j);
211 fix_pressure(0,0,0.0);
214 num_bound = adv_diff_mesh_pt()->nboundary();
215 for(
unsigned ibound=0;ibound<num_bound;ibound++)
221 unsigned num_nod= adv_diff_mesh_pt()->nboundary_node(ibound);
222 for (
unsigned inod=0;inod<num_nod;inod++)
227 if ((ibound==1) || (ibound==3))
232 for(
unsigned j=0;j<val_max;j++)
234 adv_diff_mesh_pt()->boundary_node_pt(ibound,inod)->pin(j);
245 unsigned n_nst_element = nst_mesh_pt()->nelement();
246 for(
unsigned i=0;i<n_nst_element;i++)
249 NST_ELEMENT *el_pt =
dynamic_cast<NST_ELEMENT*
> 250 (nst_mesh_pt()->element_pt(i));
268 el_pt->ignore_external_geometric_data();
271 el_pt->disable_ALE();
275 unsigned n_ad_element = adv_diff_mesh_pt()->nelement();
276 for(
unsigned i=0;i<n_ad_element;i++)
279 AD_ELEMENT *el_pt =
dynamic_cast<AD_ELEMENT*
> 280 (adv_diff_mesh_pt()->element_pt(i));
289 el_pt->disable_ALE();
295 el_pt->ignore_external_geometric_data();
299 add_sub_mesh(Nst_mesh_pt);
300 add_sub_mesh(Adv_diff_mesh_pt);
304 Multi_domain_functions::
305 setup_multi_domain_interactions<NST_ELEMENT,AD_ELEMENT>
306 (
this,nst_mesh_pt(),adv_diff_mesh_pt());
309 cout <<
"Number of equations: " << assign_eqn_numbers() << endl;
319 template<
class NST_ELEMENT,
class AD_ELEMENT>
324 unsigned num_bound=nst_mesh_pt()->nboundary();
325 for(
unsigned ibound=0;ibound<num_bound;ibound++)
328 unsigned num_nod=nst_mesh_pt()->nboundary_node(ibound);
329 for(
unsigned inod=0;inod<num_nod;inod++)
332 Node* nod_pt=nst_mesh_pt()->boundary_node_pt(ibound,inod);
339 if((ibound==1) || (ibound==3)) {vel_max = 1;}
342 for(
unsigned j=0;j<vel_max;j++) {nod_pt->set_value(j,0.0);}
348 double epsilon = 0.01;
351 double x = nod_pt->x(0);
355 double value = sin(2.0*MathematicalConstants::Pi*x/3.0)*
356 epsilon*time*exp(-time);
357 nod_pt->set_value(1,value);
364 num_bound=adv_diff_mesh_pt()->nboundary();
365 for(
unsigned ibound=0;ibound<num_bound;ibound++)
368 unsigned num_nod=adv_diff_mesh_pt()->nboundary_node(ibound);
369 for(
unsigned inod=0;inod<num_nod;inod++)
372 Node* nod_pt=adv_diff_mesh_pt()->boundary_node_pt(ibound,inod);
376 if(ibound==2) {nod_pt->set_value(0,-0.5);}
380 if(ibound==0) {nod_pt->set_value(0,0.5);}
390 template<
class NST_ELEMENT,
class AD_ELEMENT>
401 sprintf(filename,
"%s/fluid_soln%i.dat",Doc_info.directory().c_str(),
403 some_file.open(filename);
404 nst_mesh_pt()->output(some_file,npts);
408 sprintf(filename,
"%s/temperature_soln%i.dat",Doc_info.directory().c_str(),
410 some_file.open(filename);
411 adv_diff_mesh_pt()->output(some_file,npts);
422 int main(
int argc,
char **argv)
436 QAdvectionDiffusionElement<2,3> > ,
438 QCrouzeixRaviartElement<2> >
446 QAdvectionDiffusionBoussinesqElement<2> >
455 problem.steady_newton_solve();
464 problem.assign_initial_values_impulsive(dt);
467 unsigned n_steps = 200;
469 problem.refine_uniformly();
473 if(argc > 1) {n_steps = 5;}
476 for(
unsigned i=0;i<n_steps;++i)
478 problem.unsteady_newton_solve(dt);
void doc_solution()
Doc the solution.
double Peclet
Peclet number (identically one from our non-dimensionalisation)
RectangularQuadMesh< NST_ELEMENT > * nst_mesh_pt()
Access function to the Navier-Stokes mesh.
void actions_before_newton_solve()
Update the problem specs before solve (empty)
RectangularQuadMesh< AD_ELEMENT > * Adv_diff_mesh_pt
Mesh of advection diffusion elements.
void actions_before_adapt()
Actions before adapt:(empty)
RectangularQuadMesh< AD_ELEMENT > * adv_diff_mesh_pt()
Access function to the Advection-Diffusion mesh.
ConvectionProblem()
Constructor.
int main(int argc, char **argv)
Driver code for 2D Boussinesq convection problem.
void fix_pressure(const unsigned &e, const unsigned &pdof, const double &pvalue)
Fix pressure in element e at pressure dof pdof and set to pvalue.
~ConvectionProblem()
Destructor. Empty.
Namespace for the physical parameters in the problem.
Vector< double > Direction_of_gravity(2)
Gravity vector.
void actions_after_newton_solve()
Update the problem after solve (empty)
double Inverse_Prandtl
1/Prandtl number
void actions_before_implicit_timestep()
Actions before the timestep (update the the time-dependent boundary conditions)
double Rayleigh
Rayleigh number, set to be greater than the threshold for linear instability.
RectangularQuadMesh< NST_ELEMENT > * Nst_mesh_pt
Mesh of Navier Stokes elements.
void set_boundary_conditions(const double &time)
Set the boundary conditions.