33 #include "advection_diffusion.h"
34 #include "navier_stokes.h"
35 #include "multi_physics.h"
38 #include "meshes/rectangular_quadmesh.h"
41 using namespace oomph;
74 template<
class NST_ELEMENT,
class AD_ELEMENT>
99 set_boundary_conditions(time_pt()->time());
104 const double &pvalue)
107 dynamic_cast<NST_ELEMENT*
>(nst_mesh_pt()->element_pt(e))->
108 fix_pressure(pdof,pvalue);
120 return dynamic_cast<RectangularQuadMesh<NST_ELEMENT>*
>(Nst_mesh_pt);
126 return dynamic_cast<RectangularQuadMesh<AD_ELEMENT>*
>(Adv_diff_mesh_pt);
147 template<
class NST_ELEMENT,
class AD_ELEMENT>
152 add_time_stepper_pt(
new BDF<2>);
155 Doc_info.set_directory(
"RESLT");
171 new RectangularQuadMesh<NST_ELEMENT>(n_x,n_y,l_x,l_y,time_stepper_pt());
173 new RectangularQuadMesh<AD_ELEMENT>(n_x,n_y,l_x,l_y,time_stepper_pt());
180 unsigned num_bound = nst_mesh_pt()->nboundary();
181 for(
unsigned ibound=0;ibound<num_bound;ibound++)
187 unsigned num_nod= nst_mesh_pt()->nboundary_node(ibound);
188 for (
unsigned inod=0;inod<num_nod;inod++)
192 if ((ibound==1) || (ibound==3))
198 for(
unsigned j=0;j<val_max;j++)
200 nst_mesh_pt()->boundary_node_pt(ibound,inod)->pin(j);
207 fix_pressure(0,0,0.0);
210 num_bound = adv_diff_mesh_pt()->nboundary();
211 for(
unsigned ibound=0;ibound<num_bound;ibound++)
217 unsigned num_nod= adv_diff_mesh_pt()->nboundary_node(ibound);
218 for (
unsigned inod=0;inod<num_nod;inod++)
223 if ((ibound==1) || (ibound==3))
228 for(
unsigned j=0;j<val_max;j++)
230 adv_diff_mesh_pt()->boundary_node_pt(ibound,inod)->pin(j);
241 unsigned n_nst_element = nst_mesh_pt()->nelement();
242 for(
unsigned i=0;i<n_nst_element;i++)
245 NST_ELEMENT *el_pt =
dynamic_cast<NST_ELEMENT*
>
246 (nst_mesh_pt()->element_pt(i));
264 el_pt->ignore_external_geometric_data();
267 el_pt->disable_ALE();
271 unsigned n_ad_element = adv_diff_mesh_pt()->nelement();
272 for(
unsigned i=0;i<n_ad_element;i++)
275 AD_ELEMENT *el_pt =
dynamic_cast<AD_ELEMENT*
>
276 (adv_diff_mesh_pt()->element_pt(i));
285 el_pt->disable_ALE();
291 el_pt->ignore_external_geometric_data();
295 add_sub_mesh(Nst_mesh_pt);
296 add_sub_mesh(Adv_diff_mesh_pt);
300 Multi_domain_functions::
301 setup_multi_domain_interactions<NST_ELEMENT,AD_ELEMENT>
302 (
this,nst_mesh_pt(),adv_diff_mesh_pt());
305 cout <<
"Number of equations: " << assign_eqn_numbers() << endl;
315 template<
class NST_ELEMENT,
class AD_ELEMENT>
320 unsigned num_bound=nst_mesh_pt()->nboundary();
321 for(
unsigned ibound=0;ibound<num_bound;ibound++)
324 unsigned num_nod=nst_mesh_pt()->nboundary_node(ibound);
325 for(
unsigned inod=0;inod<num_nod;inod++)
328 Node* nod_pt=nst_mesh_pt()->boundary_node_pt(ibound,inod);
335 if((ibound==1) || (ibound==3)) {vel_max = 1;}
338 for(
unsigned j=0;j<vel_max;j++) {nod_pt->set_value(j,0.0);}
344 double epsilon = 0.01;
347 double x = nod_pt->x(0);
351 double value = sin(2.0*MathematicalConstants::Pi*x/3.0)*
352 epsilon*time*exp(-time);
353 nod_pt->set_value(1,value);
360 num_bound=adv_diff_mesh_pt()->nboundary();
361 for(
unsigned ibound=0;ibound<num_bound;ibound++)
364 unsigned num_nod=adv_diff_mesh_pt()->nboundary_node(ibound);
365 for(
unsigned inod=0;inod<num_nod;inod++)
368 Node* nod_pt=adv_diff_mesh_pt()->boundary_node_pt(ibound,inod);
372 if(ibound==2) {nod_pt->set_value(0,-0.5);}
376 if(ibound==0) {nod_pt->set_value(0,0.5);}
386 template<
class NST_ELEMENT,
class AD_ELEMENT>
397 sprintf(filename,
"%s/fluid_soln%i.dat",Doc_info.directory().c_str(),
399 some_file.open(filename);
400 nst_mesh_pt()->output(some_file,npts);
404 sprintf(filename,
"%s/temperature_soln%i.dat",Doc_info.directory().c_str(),
406 some_file.open(filename);
407 adv_diff_mesh_pt()->output(some_file,npts);
418 int main(
int argc,
char **argv)
432 QAdvectionDiffusionElement<2,3> > ,
434 QCrouzeixRaviartElement<2> >
442 QAdvectionDiffusionBoussinesqElement<2> >
451 problem.steady_newton_solve();
454 problem.doc_solution();
460 problem.assign_initial_values_impulsive(dt);
463 unsigned n_steps = 200;
465 problem.refine_uniformly();
469 if(argc > 1) {n_steps = 5;}
472 for(
unsigned i=0;i<n_steps;++i)
474 problem.unsteady_newton_solve(dt);
475 problem.doc_solution();
/////////////////////////////////////////////////////////////////// /////////////////////////////////...
ConvectionProblem()
Constructor.
ConvectionProblem()
Constructor.
RectangularQuadMesh< AD_ELEMENT > * adv_diff_mesh_pt()
Access function to the Advection-Diffusion mesh.
void actions_before_newton_solve()
Update the problem specs before solve (empty)
void actions_before_implicit_timestep()
Actions before the timestep (update the the time-dependent boundary conditions)
void set_boundary_conditions(const double &time)
Set the boundary conditions.
RectangularQuadMesh< NST_ELEMENT > * Nst_mesh_pt
Mesh of Navier Stokes elements.
RectangularQuadMesh< NST_ELEMENT > * nst_mesh_pt()
Access function to the Navier-Stokes mesh.
void doc_solution()
Doc the solution.
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.
void set_boundary_conditions(const double &time)
Set the boundary conditions.
void actions_after_newton_solve()
Update the problem after solve (empty)
~ConvectionProblem()
Destructor. Empty.
void doc_solution()
Doc the solution.
RectangularQuadMesh< AD_ELEMENT > * Adv_diff_mesh_pt
Mesh of advection diffusion elements.
void actions_before_adapt()
Actions before adapt:(empty)
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
////////////////////////////////////////////////////////////////////// //////////////////////////////...
int main(int argc, char **argv)
Driver code for 2D Boussinesq convection problem.
Namespace for the physical parameters in the problem.
Vector< double > Direction_of_gravity(2)
Gravity vector.
double Rayleigh
Rayleigh number, set to be greater than the threshold for linear instability.
double Inverse_Prandtl
1/Prandtl number
double Peclet
Peclet number (identically one from our non-dimensionalisation)