multi_domain_boussinesq_convection.cc
Go to the documentation of this file.
1//LIC// ====================================================================
2//LIC// This file forms part of oomph-lib, the object-oriented,
3//LIC// multi-physics finite-element library, available
4//LIC// at http://www.oomph-lib.org.
5//LIC//
6//LIC// Copyright (C) 2006-2025 Matthias Heil and Andrew Hazel
7//LIC//
8//LIC// This library is free software; you can redistribute it and/or
9//LIC// modify it under the terms of the GNU Lesser General Public
10//LIC// License as published by the Free Software Foundation; either
11//LIC// version 2.1 of the License, or (at your option) any later version.
12//LIC//
13//LIC// This library is distributed in the hope that it will be useful,
14//LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15//LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16//LIC// Lesser General Public License for more details.
17//LIC//
18//LIC// You should have received a copy of the GNU Lesser General Public
19//LIC// License along with this library; if not, write to the Free Software
20//LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21//LIC// 02110-1301 USA.
22//LIC//
23//LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24//LIC//
25//LIC//====================================================================
26//Driver for a multi-physics problem that couples a Navier--Stokes
27//mesh to an advection diffusion mesh, giving Boussinesq convection
28
29
30//Oomph-lib headers, we require the generic, advection-diffusion,
31//and navier-stokes elements.
32#include "generic.h"
33#include "advection_diffusion.h"
34#include "navier_stokes.h"
35#include "multi_physics.h"
36
37// Both meshes are the standard rectangular quadmesh
38#include "meshes/rectangular_quadmesh.h"
39
40// Use the oomph and std namespaces
41using namespace oomph;
42using namespace std;
43
44
45//======start_of_namespace============================================
46/// Namespace for the physical parameters in the problem
47//====================================================================
49{
50 /// Peclet number (identically one from our non-dimensionalisation)
51 double Peclet=1.0;
52
53 /// 1/Prandtl number
54 double Inverse_Prandtl=1.0;
55
56 /// Rayleigh number, set to be greater than
57 /// the threshold for linear instability
58 double Rayleigh = 1800.0;
59
60 /// Gravity vector
61 Vector<double> Direction_of_gravity(2);
62
63} // end_of_namespace
64
65//////////////////////////////////////////////////////////////////////
66//////////////////////////////////////////////////////////////////////
67//////////////////////////////////////////////////////////////////////
68
69//====== start_of_problem_class=======================================
70/// 2D Convection problem on two rectangular domains, discretised
71/// with Navier-Stokes and Advection-Diffusion elements. The specific type
72/// of elements is specified via the template parameters.
73//====================================================================
74template<class NST_ELEMENT,class AD_ELEMENT>
75class ConvectionProblem : public Problem
76{
77
78public:
79
80 /// Constructor
82
83 /// Destructor. Empty
85
86 /// Update the problem specs before solve (empty)
88
89 /// Update the problem after solve (empty)
91
92 /// Actions before adapt:(empty)
94
95 /// Actions before the timestep (update the the time-dependent
96 /// boundary conditions)
98 {
99 set_boundary_conditions(time_pt()->time());
100 }
101
102 /// Fix pressure in element e at pressure dof pdof and set to pvalue
103 void fix_pressure(const unsigned &e, const unsigned &pdof,
104 const double &pvalue)
105 {
106 //Cast to specific element and fix pressure
107 dynamic_cast<NST_ELEMENT*>(nst_mesh_pt()->element_pt(e))->
108 fix_pressure(pdof,pvalue);
109 } // end_of_fix_pressure
110
111 /// Doc the solution.
113
114 /// Set the boundary conditions
115 void set_boundary_conditions(const double &time);
116
117 /// Access function to the Navier-Stokes mesh
118 RectangularQuadMesh<NST_ELEMENT>* nst_mesh_pt()
119 {
120 return dynamic_cast<RectangularQuadMesh<NST_ELEMENT>*>(Nst_mesh_pt);
121 }
122
123 /// Access function to the Advection-Diffusion mesh
124 RectangularQuadMesh<AD_ELEMENT>* adv_diff_mesh_pt()
125 {
126 return dynamic_cast<RectangularQuadMesh<AD_ELEMENT>*>(Adv_diff_mesh_pt);
127 }
128
129private:
130
131 /// DocInfo object
132 DocInfo Doc_info;
133
134protected:
135
136 /// Mesh of Navier Stokes elements
137 RectangularQuadMesh<NST_ELEMENT>* Nst_mesh_pt;
138
139 /// Mesh of advection diffusion elements
140 RectangularQuadMesh<AD_ELEMENT>* Adv_diff_mesh_pt;
141
142}; // end of problem class
143
144//===========start_of_constructor=========================================
145/// Constructor for convection problem
146//========================================================================
147template<class NST_ELEMENT,class AD_ELEMENT>
149{
150
151 //Allocate a timestepper
152 add_time_stepper_pt(new BDF<2>);
153
154 // Set output directory
155 Doc_info.set_directory("RESLT");
156
157 // # of elements in x-direction
158 unsigned n_x=8;
159
160 // # of elements in y-direction
161 unsigned n_y=8;
162
163 // Domain length in x-direction
164 double l_x=3.0;
165
166 // Domain length in y-direction
167 double l_y=1.0;
168
169 // Build two standard rectangular quadmesh
170 Nst_mesh_pt =
171 new RectangularQuadMesh<NST_ELEMENT>(n_x,n_y,l_x,l_y,time_stepper_pt());
172 Adv_diff_mesh_pt =
173 new RectangularQuadMesh<AD_ELEMENT>(n_x,n_y,l_x,l_y,time_stepper_pt());
174
175 // Set the boundary conditions for this problem: All nodes are
176 // free by default -- only need to pin the ones that have Dirichlet
177 // conditions here
178
179 //Loop over the boundaries
180 unsigned num_bound = nst_mesh_pt()->nboundary();
181 for(unsigned ibound=0;ibound<num_bound;ibound++)
182 {
183 //Set the maximum index to be pinned (both velocity values by default)
184 unsigned val_max=2;
185
186 //Loop over the number of nodes on the boundry
187 unsigned num_nod= nst_mesh_pt()->nboundary_node(ibound);
188 for (unsigned inod=0;inod<num_nod;inod++)
189 {
190 //If we are on the side-walls, the v-velocity satisfies natural
191 //boundary conditions, so we only pin the u-velocity
192 if ((ibound==1) || (ibound==3))
193 {
194 val_max=1;
195 }
196
197 //Loop over the desired values stored at the nodes and pin
198 for(unsigned j=0;j<val_max;j++)
199 {
200 nst_mesh_pt()->boundary_node_pt(ibound,inod)->pin(j);
201 }
202 }
203 }
204
205 //Pin the zero-th pressure dof in element 0 and set its value to
206 //zero:
207 fix_pressure(0,0,0.0);
208
209 //Loop over the boundaries of the AD mesh
210 num_bound = adv_diff_mesh_pt()->nboundary();
211 for(unsigned ibound=0;ibound<num_bound;ibound++)
212 {
213 //Set the maximum index to be pinned (the temperature value by default)
214 unsigned val_max=1;
215
216 //Loop over the number of nodes on the boundry
217 unsigned num_nod= adv_diff_mesh_pt()->nboundary_node(ibound);
218 for (unsigned inod=0;inod<num_nod;inod++)
219 {
220 //If we are on the side-walls, the temperature
221 //satisfies natural boundary conditions, so we don't pin anything
222 // in this mesh
223 if ((ibound==1) || (ibound==3))
224 {
225 val_max=0;
226 }
227 //Loop over the desired values stored at the nodes and pin
228 for(unsigned j=0;j<val_max;j++)
229 {
230 adv_diff_mesh_pt()->boundary_node_pt(ibound,inod)->pin(j);
231 }
232 }
233 }
234
235
236 // Complete the build of all elements so they are fully functional
237
238 // Loop over the elements to set up element-specific
239 // things that cannot be handled by the (argument-free!) ELEMENT
240 // constructors.
241 unsigned n_nst_element = nst_mesh_pt()->nelement();
242 for(unsigned i=0;i<n_nst_element;i++)
243 {
244 // Upcast from GeneralsedElement to the present element
245 NST_ELEMENT *el_pt = dynamic_cast<NST_ELEMENT*>
246 (nst_mesh_pt()->element_pt(i));
247
248 // Set the Reynolds number (1/Pr in our non-dimensionalisation)
250
251 // Set ReSt (also 1/Pr in our non-dimensionalisation)
253
254 // Set the Rayleigh number
255 el_pt->ra_pt() = &Global_Physical_Variables::Rayleigh;
256
257 //Set Gravity vector
259
260 // We can ignore the external geometric data in the "external"
261 // advection diffusion element when computing the Jacobian matrix
262 // because the interaction does not involve spatial gradients of
263 // the temperature (and also because the mesh isn't moving!)
264 el_pt->ignore_external_geometric_data();
265
266 //The mesh is fixed, so we can disable ALE
267 el_pt->disable_ALE();
268
269 }
270
271 unsigned n_ad_element = adv_diff_mesh_pt()->nelement();
272 for(unsigned i=0;i<n_ad_element;i++)
273 {
274 // Upcast from GeneralsedElement to the present element
275 AD_ELEMENT *el_pt = dynamic_cast<AD_ELEMENT*>
276 (adv_diff_mesh_pt()->element_pt(i));
277
278 // Set the Peclet number
279 el_pt->pe_pt() = &Global_Physical_Variables::Peclet;
280
281 // Set the Peclet number multiplied by the Strouhal number
282 el_pt->pe_st_pt() =&Global_Physical_Variables::Peclet;
283
284 //The mesh is fixed, so we can disable ALE
285 el_pt->disable_ALE();
286
287 // We can ignore the external geometric data in the "external"
288 // advection diffusion element when computing the Jacobian matrix
289 // because the interaction does not involve spatial gradients of
290 // the temperature (and also because the mesh isn't moving!)
291 el_pt->ignore_external_geometric_data();
292 }
293
294 // combine the submeshes
295 add_sub_mesh(Nst_mesh_pt);
296 add_sub_mesh(Adv_diff_mesh_pt);
297 build_global_mesh();
298
299 // Set sources
300 Multi_domain_functions::
301 setup_multi_domain_interactions<NST_ELEMENT,AD_ELEMENT>
302 (this,nst_mesh_pt(),adv_diff_mesh_pt());
303
304 // Setup equation numbering scheme
305 cout <<"Number of equations: " << assign_eqn_numbers() << endl;
306
307} // end of constructor
308
309
310
311//===========start_of_set_boundary_conditions================
312/// Set the boundary conditions as a function of continuous
313/// time
314//===========================================================
315template<class NST_ELEMENT,class AD_ELEMENT>
317 const double &time)
318{
319 // Loop over all the boundaries on the NST mesh
320 unsigned num_bound=nst_mesh_pt()->nboundary();
321 for(unsigned ibound=0;ibound<num_bound;ibound++)
322 {
323 // Loop over the nodes on boundary
324 unsigned num_nod=nst_mesh_pt()->nboundary_node(ibound);
325 for(unsigned inod=0;inod<num_nod;inod++)
326 {
327 // Get pointer to node
328 Node* nod_pt=nst_mesh_pt()->boundary_node_pt(ibound,inod);
329
330 //Set the number of velocity components to be pinned
331 //(both by default)
332 unsigned vel_max=2;
333
334 //If we are on the side walls we only set the x-velocity.
335 if((ibound==1) || (ibound==3)) {vel_max = 1;}
336
337 //Set the pinned velocities to zero on NST mesh
338 for(unsigned j=0;j<vel_max;j++) {nod_pt->set_value(j,0.0);}
339
340 //If we are on the top boundary
341 if(ibound==2)
342 {
343 //Add small velocity imperfection if desired
344 double epsilon = 0.01;
345
346 //Read out the x position
347 double x = nod_pt->x(0);
348
349 //Set a sinusoidal perturbation in the vertical velocity
350 //This perturbation is mass conserving
351 double value = sin(2.0*MathematicalConstants::Pi*x/3.0)*
352 epsilon*time*exp(-time);
353 nod_pt->set_value(1,value);
354 }
355
356 }
357 }
358
359 // Loop over all the boundaries on the AD mesh
360 num_bound=adv_diff_mesh_pt()->nboundary();
361 for(unsigned ibound=0;ibound<num_bound;ibound++)
362 {
363 // Loop over the nodes on boundary
364 unsigned num_nod=adv_diff_mesh_pt()->nboundary_node(ibound);
365 for(unsigned inod=0;inod<num_nod;inod++)
366 {
367 // Get pointer to node
368 Node* nod_pt=adv_diff_mesh_pt()->boundary_node_pt(ibound,inod);
369
370 //If we are on the top boundary, set the temperature
371 //to -0.5 (cooled)
372 if(ibound==2) {nod_pt->set_value(0,-0.5);}
373
374 //If we are on the bottom boundary, set the temperature
375 //to 0.5 (heated)
376 if(ibound==0) {nod_pt->set_value(0,0.5);}
377 }
378 }
379
380
381} // end_of_set_boundary_conditions
382
383//===============start_doc_solution=======================================
384/// Doc the solution
385//========================================================================
386template<class NST_ELEMENT,class AD_ELEMENT>
388{
389 //Declare an output stream and filename
390 ofstream some_file;
391 char filename[100];
392
393 // Number of plot points: npts x npts
394 unsigned npts=5;
395
396 // Output Navier-Stokes solution
397 snprintf(filename, sizeof(filename), "%s/fluid_soln%i.dat",Doc_info.directory().c_str(),
398 Doc_info.number());
399 some_file.open(filename);
400 nst_mesh_pt()->output(some_file,npts);
401 some_file.close();
402
403 // Output advection diffusion solution
404 snprintf(filename, sizeof(filename), "%s/temperature_soln%i.dat",Doc_info.directory().c_str(),
405 Doc_info.number());
406 some_file.open(filename);
407 adv_diff_mesh_pt()->output(some_file,npts);
408 some_file.close();
409
410 Doc_info.number()++;
411
412} // end of doc
413
414
415//=======start_of_main================================================
416/// Driver code for 2D Boussinesq convection problem
417//====================================================================
418int main(int argc, char **argv)
419{
420
421 // Set the direction of gravity
424
425#define NEW
426//#undef NEW
427#ifdef NEW
428
429 //Construct our problem
432 QAdvectionDiffusionElement<2,3> > ,
434 QCrouzeixRaviartElement<2> >
435 >
436 problem;
437
438#else
439
440//Construct our problem
442 QAdvectionDiffusionBoussinesqElement<2> >
443 problem;
444
445#endif
446
447 // Apply the boundary condition at time zero
448 problem.set_boundary_conditions(0.0);
449
450 //Perform a single steady Newton solve
451 problem.steady_newton_solve();
452
453 //Document the solution
454 problem.doc_solution();
455
456 //Set the timestep
457 double dt = 0.1;
458
459 //Initialise the value of the timestep and set an impulsive start
460 problem.assign_initial_values_impulsive(dt);
461
462 //Set the number of timesteps to our default value
463 unsigned n_steps = 200;
464
465 problem.refine_uniformly();
466
467 //If we have a command line argument, perform fewer steps
468 //(used for self-test runs)
469 if(argc > 1) {n_steps = 5;}
470
471 //Perform n_steps timesteps
472 for(unsigned i=0;i<n_steps;++i)
473 {
474 problem.unsteady_newton_solve(dt);
475 problem.doc_solution();
476 }
477
478} // end of main
479
480
481
482
483
484
485
486
487
2D Convection problem on rectangular domain, discretised with refineable elements....
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.
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)
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)
DocInfo Doc_info
DocInfo object.
RectangularQuadMesh< NST_ELEMENT > * nst_mesh_pt()
Access function to the Navier-Stokes mesh.
Build AdvectionDiffusionBoussinesqElement that inherits from ElementWithExternalElement so that it ca...
Build NavierStokesBoussinesqElement that inherits from ElementWithExternalElement so that it can "com...
int main()
Driver code for 2D Boussinesq convection problem with adaptivity.
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)