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 the Navier--Stokes
31 //equations to the advection diffusion equations,
32 //giving Boussinesq convection
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 // The mesh is our 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 rectangular domain, discretised
75 /// with refineable elements. The specific type
76 /// of element is specified via the template parameter.
77 //====================================================================
78 template<class 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<ELEMENT*>(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 Overloaded version of the problem's access function to
122  /// the mesh. Recasts the pointer to the base Mesh object to
123  /// the actual mesh type.
124  RectangularQuadMesh<ELEMENT>* mesh_pt()
125  {
126  return dynamic_cast<RectangularQuadMesh<ELEMENT>*>(
127  Problem::mesh_pt());
128  }
129 
130 private:
131 
132  /// DocInfo object
133  DocInfo Doc_info;
134 
135 }; // end of problem class
136 
137 //===========start_of_constructor=========================================
138 /// \short Constructor for convection problem
139 //========================================================================
140 template<class ELEMENT>
142 {
143  //Allocate a timestepper
144  add_time_stepper_pt(new BDF<2>);
145 
146  // Set output directory
147  Doc_info.set_directory("RESLT");
148 
149  // # of elements in x-direction
150  unsigned n_x=8;
151 
152  // # of elements in y-direction
153  unsigned n_y=8;
154 
155  // Domain length in x-direction
156  double l_x=3.0;
157 
158  // Domain length in y-direction
159  double l_y=1.0;
160 
161  // Build a standard rectangular quadmesh
162  Problem::mesh_pt() =
163  new RectangularQuadMesh<ELEMENT>(n_x,n_y,l_x,l_y,time_stepper_pt());
164 
165  // Set the boundary conditions for this problem: All nodes are
166  // free by default -- only need to pin the ones that have Dirichlet
167  // conditions here
168 
169  //Loop over the boundaries
170  unsigned num_bound = mesh_pt()->nboundary();
171  for(unsigned ibound=0;ibound<num_bound;ibound++)
172  {
173  //Set the maximum index to be pinned (all values by default)
174  unsigned val_max=3;
175  //If we are on the side-walls, the v-velocity and temperature
176  //satisfy natural boundary conditions, so we only pin the
177  //first value
178  if((ibound==1) || (ibound==3)) {val_max=1;}
179 
180  //Loop over the number of nodes on the boundry
181  unsigned num_nod= mesh_pt()->nboundary_node(ibound);
182  for (unsigned inod=0;inod<num_nod;inod++)
183  {
184  //Loop over the desired values stored at the nodes and pin
185  for(unsigned j=0;j<val_max;j++)
186  {
187  mesh_pt()->boundary_node_pt(ibound,inod)->pin(j);
188  }
189  }
190  }
191 
192  //Pin the zero-th pressure dof in element 0 and set its value to
193  //zero:
194  fix_pressure(0,0,0.0);
195 
196  // Complete the build of all elements so they are fully functional
197 
198  // Loop over the elements to set up element-specific
199  // things that cannot be handled by the (argument-free!) ELEMENT
200  // constructor.
201  unsigned n_element = mesh_pt()->nelement();
202  for(unsigned i=0;i<n_element;i++)
203  {
204  // Upcast from GeneralsedElement to the present element
205  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(mesh_pt()->element_pt(i));
206 
207  // Set the Peclet number
208  el_pt->pe_pt() = &Global_Physical_Variables::Peclet;
209 
210  // Set the Peclet number multiplied by the Strouhal number
211  el_pt->pe_st_pt() =&Global_Physical_Variables::Peclet;
212 
213  // Set the Reynolds number (1/Pr in our non-dimensionalisation)
215 
216  // Set ReSt (also 1/Pr in our non-dimensionalisation)
217  el_pt->re_st_pt() = &Global_Physical_Variables::Inverse_Prandtl;
218 
219  // Set the Rayleigh number
220  el_pt->ra_pt() = &Global_Physical_Variables::Rayleigh;
221 
222  //Set Gravity vector
224 
225  //The mesh is fixed, so we can disable ALE
226  el_pt->disable_ALE();
227  }
228 
229  // Setup equation numbering scheme
230  cout <<"Number of equations: " << assign_eqn_numbers() << endl;
231 
232 } // end of constructor
233 
234 
235 
236 //===========start_of_set_boundary_conditions================
237 /// Set the boundary conditions as a function of continuous
238 /// time
239 //===========================================================
240 template<class ELEMENT>
242  const double &time)
243 {
244  // Loop over the boundaries
245  unsigned num_bound = mesh_pt()->nboundary();
246  for(unsigned ibound=0;ibound<num_bound;ibound++)
247  {
248  // Loop over the nodes on boundary
249  unsigned num_nod=mesh_pt()->nboundary_node(ibound);
250  for(unsigned inod=0;inod<num_nod;inod++)
251  {
252  // Get pointer to node
253  Node* nod_pt=mesh_pt()->boundary_node_pt(ibound,inod);
254 
255  //Set the number of velocity components
256  unsigned vel_max=2;
257 
258  //If we are on the side walls we only set the x-velocity.
259  if((ibound==1) || (ibound==3)) {vel_max = 1;}
260 
261  //Set the pinned velocities to zero
262  for(unsigned j=0;j<vel_max;j++) {nod_pt->set_value(j,0.0);}
263 
264  //If we are on the top boundary
265  if(ibound==2)
266  {
267  //Set the temperature to -0.5 (cooled)
268  nod_pt->set_value(2,-0.5);
269 
270  //Add small velocity imperfection if desired
271  double epsilon = 0.01;
272 
273  //Read out the x position
274  double x = nod_pt->x(0);
275 
276  //Set a sinusoidal perturbation in the vertical velocity
277  //This perturbation is mass conserving
278  double value = sin(2.0*MathematicalConstants::Pi*x/3.0)*
279  epsilon*time*exp(-time);
280  nod_pt->set_value(1,value);
281  }
282 
283  //If we are on the bottom boundary, set the temperature
284  //to 0.5 (heated)
285  if(ibound==0) {nod_pt->set_value(2,0.5);}
286  }
287  }
288 } // end_of_set_boundary_conditions
289 
290 //===============start_doc_solution=======================================
291 /// Doc the solution
292 //========================================================================
293 template<class ELEMENT>
295 {
296  //Declare an output stream and filename
297  ofstream some_file;
298  char filename[100];
299 
300  // Number of plot points: npts x npts
301  unsigned npts=5;
302 
303  // Output solution
304  //-----------------
305  sprintf(filename,"%s/soln%i.dat",Doc_info.directory().c_str(),
306  Doc_info.number());
307  some_file.open(filename);
308  mesh_pt()->output(some_file,npts);
309  some_file.close();
310 
311  Doc_info.number()++;
312 } // end of doc
313 
314 
315 //=======start_of_main================================================
316 /// Driver code for 2D Boussinesq convection problem
317 //====================================================================
318 int main(int argc, char **argv)
319 {
320 
321  // Set the direction of gravity
324 
325  //Construct our problem
327 
328  // Apply the boundary condition at time zero
329  problem.set_boundary_conditions(0.0);
330 
331  //Perform a single steady Newton solve
332  problem.steady_newton_solve();
333 
334  //Document the solution
335  problem.doc_solution();
336 
337  //Set the timestep
338  double dt = 0.1;
339 
340  //Initialise the value of the timestep and set an impulsive start
341  problem.assign_initial_values_impulsive(dt);
342 
343  //Set the number of timesteps to our default value
344  unsigned n_steps = 200;
345 
346  //If we have a command line argument, perform fewer steps
347  //(used for self-test runs)
348  if(argc > 1) {n_steps = 5;}
349 
350  //Perform n_steps timesteps
351  for(unsigned i=0;i<n_steps;++i)
352  {
353  problem.unsteady_newton_solve(dt);
354  problem.doc_solution();
355  }
356 
357 } // end of main
358 
359 
360 
361 
362 
363 
364 
365 
366 
void doc_solution()
Doc the solution.
double Peclet
Peclet number (identically one from our non-dimensionalisation)
void actions_before_newton_solve()
Update the problem specs before solve (empty)
void actions_before_adapt()
Actions before adapt:(empty)
RectangularQuadMesh< ELEMENT > * mesh_pt()
Overloaded version of the problem&#39;s access function to the mesh. Recasts the pointer to the base Mesh...
ConvectionProblem()
Constructor.
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.
DocInfo Doc_info
DocInfo object.
int main(int argc, char **argv)
Driver code for 2D Boussinesq convection 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.
void set_boundary_conditions(const double &time)
Set the boundary conditions.