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