refineable_b_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 //////////////////////////////////////////////////////////////////////
50 //////////////////////////////////////////////////////////////////////
51 //////////////////////////////////////////////////////////////////////
52 
53 
54 //======start_of_namespace============================================
55 /// Namespace for the physical parameters in the problem
56 //====================================================================
58 {
59 
60  /// Peclet number (identically one from our non-dimensionalisation)
61  double Peclet=1.0;
62 
63  /// 1/Prandtl number
64  double Inverse_Prandtl=1.0;
65 
66  /// \short Rayleigh number, set to be greater than
67  /// the threshold for linear instability
68  double Rayleigh = 1800.0;
69 
70  /// Gravity vector
71  Vector<double> Direction_of_gravity(2);
72 
73 } // end_of_namespace
74 
75 
76 //////////////////////////////////////////////////////////////////////
77 //////////////////////////////////////////////////////////////////////
78 //////////////////////////////////////////////////////////////////////
79 
80 
81 
82 //=======start_of_problem_class=======================================
83 /// 2D Convection problem on rectangular domain, discretised
84 /// with refineable elements. The specific type
85 /// of element is specified via the template parameter.
86 //====================================================================
87 template<class ELEMENT>
88 class RefineableConvectionProblem : public Problem
89 {
90 
91 public:
92 
93  ///Constructor
95 
96  /// Destructor. Empty
98 
99  /// \short Update the problem specs before solve:
100  void actions_before_newton_solve();
101 
102  /// Update the problem after solve (empty)
104 
105  /// \short Overloaded version of the problem's access function to
106  /// the mesh. Recasts the pointer to the base Mesh object to
107  /// the actual mesh type.
108  RectangularQuadMesh<ELEMENT>* mesh_pt()
109  {
110  return dynamic_cast<RectangularQuadMesh<ELEMENT>*>(
111  Problem::mesh_pt());
112  } //end of access function to specic mesh
113 
114  /// Actions before adapt:(empty)
116 
117  /// \short Actions after adaptation,
118  /// Re-pin a single pressure degree of freedom
120  {
121  //Unpin all the pressures to avoid pinning two pressures
122  RefineableNavierStokesEquations<2>::
123  unpin_all_pressure_dofs(mesh_pt()->element_pt());
124 
125  //Pin the zero-th pressure dof in the zero-th element and set
126  // its value to zero
127  fix_pressure(0,0,0.0);
128  }
129 
130  ///Fix pressure in element e at pressure dof pdof and set to pvalue
131  void fix_pressure(const unsigned &e, const unsigned &pdof,
132  const double &pvalue)
133  {
134  //Cast to specific element and fix pressure
135  dynamic_cast<ELEMENT*>(mesh_pt()->element_pt(e))->
136  fix_pressure(pdof,pvalue);
137  } // end_of_fix_pressure
138 
139  /// \short Set the
140  /// boundary condition on the upper wall to be perturbed slightly
141  /// to force the solution into the symmetry broken state.
142  void enable_imperfection() {Imperfect = true;}
143 
144  /// \short Set the
145  /// boundary condition on the upper wall to be unperturbed.
146  void disable_imperfection() {Imperfect = false;}
147 
148  /// \short Doc the solution.
149  void doc_solution();
150 
151 private:
152 
153  /// DocInfo object
154  DocInfo Doc_info;
155 
156  /// Is the boundary condition imperfect or not
157  bool Imperfect;
158 
159 }; // end of problem class
160 
161 
162 //=======start_of_constructor=============================================
163 /// Constructor for adaptive thermal convection problem
164 //========================================================================
165 template<class ELEMENT>
167 RefineableConvectionProblem() : Imperfect(false)
168 {
169  // Set output directory
170  Doc_info.set_directory("RESLT");
171 
172  // # of elements in x-direction
173  unsigned n_x=9;
174 
175  // # of elements in y-direction
176  unsigned n_y=8;
177 
178  // Domain length in x-direction
179  double l_x=3.0;
180 
181  // Domain length in y-direction
182  double l_y=1.0;
183 
184  // Build the mesh
185  RefineableRectangularQuadMesh<ELEMENT>* cast_mesh_pt =
186  new RefineableRectangularQuadMesh<ELEMENT>(n_x,n_y,l_x,l_y);
187 
188  //Set the problem's mesh pointer
189  Problem::mesh_pt() = cast_mesh_pt;
190 
191 
192  // Create/set error estimator
193  cast_mesh_pt->spatial_error_estimator_pt()=new Z2ErrorEstimator;
194 
195  // Set error targets for adaptive refinement
196  cast_mesh_pt->max_permitted_error()=0.5e-3;
197  cast_mesh_pt->min_permitted_error()=0.5e-4;
198 
199 
200  // Set the boundary conditions for this problem: All nodes are
201  // free by default -- only need to pin the ones that have Dirichlet
202  // conditions here
203 
204  //Loop over the boundaries
205  unsigned num_bound = mesh_pt()->nboundary();
206  for(unsigned ibound=0;ibound<num_bound;ibound++)
207  {
208  //Set the maximum index to be pinned (all values by default)
209  unsigned val_max=3;
210  //If we are on the side-walls, the v-velocity and temperature
211  //satisfy natural boundary conditions, so we only pin the
212  //first value
213  if((ibound==1) || (ibound==3)) {val_max=1;}
214 
215  //Loop over the number of nodes on the boundry
216  unsigned num_nod= mesh_pt()->nboundary_node(ibound);
217  for (unsigned inod=0;inod<num_nod;inod++)
218  {
219  //Loop over the desired values stored at the nodes and pin
220  for(unsigned j=0;j<val_max;j++)
221  {
222  mesh_pt()->boundary_node_pt(ibound,inod)->pin(j);
223  }
224  }
225  }
226 
227  // Pin the zero-th pressure value in the zero-th element and
228  // set its value to zero.
229  fix_pressure(0,0,0.0);
230 
231  // Complete the build of all elements so they are fully functional
232 
233  // Loop over the elements to set up element-specific
234  // things that cannot be handled by the (argument-free!) ELEMENT
235  // constructor.
236  unsigned n_element = mesh_pt()->nelement();
237  for(unsigned i=0;i<n_element;i++)
238  {
239  // Upcast from GeneralsedElement to the present element
240  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(mesh_pt()->element_pt(i));
241 
242  // Set the Peclet number
243  el_pt->pe_pt() = &Global_Physical_Variables::Peclet;
244 
245  // Set the Peclet Strouhal number
246  el_pt->pe_st_pt() = &Global_Physical_Variables::Peclet;
247 
248  // Set the Reynolds number (1/Pr in our non-dimensionalisation)
250 
251  // Set ReSt (also 1/Pr in our non-dimensionalisation)
252  el_pt->re_st_pt() = &Global_Physical_Variables::Inverse_Prandtl;
253 
254  // Set the Rayleigh number
255  el_pt->ra_pt() = &Global_Physical_Variables::Rayleigh;
256 
257  //Set Gravity vector
259  }
260 
261  // Setup equation numbering scheme
262  cout <<"Number of equations: " << assign_eqn_numbers() << endl;
263 
264 } // end of constructor
265 
266 
267 //===================start_actions_before_newton_solve===========================
268 /// Update the problem specs before solve: (Re-)set boundary conditions
269 /// to include an imperfection (or not) depending on the control flag.
270 //========================================================================
271 template<class ELEMENT>
273 {
274  // Loop over the boundaries
275  unsigned num_bound = mesh_pt()->nboundary();
276  for(unsigned ibound=0;ibound<num_bound;ibound++)
277  {
278  // Loop over the nodes on boundary
279  unsigned num_nod=mesh_pt()->nboundary_node(ibound);
280  for(unsigned inod=0;inod<num_nod;inod++)
281  {
282  // Get pointer to node
283  Node* nod_pt=mesh_pt()->boundary_node_pt(ibound,inod);
284 
285  //Set the number of velocity components
286  unsigned vel_max=2;
287  //If we are on the side walls we only pin the x-velocity.
288  if((ibound==1) || (ibound==3)) {vel_max = 1;}
289  //Set the pinned velocities to zero
290  for(unsigned j=0;j<vel_max;j++) {nod_pt->set_value(j,0.0);}
291 
292  //If we are on the top boundary
293  if(ibound==2)
294  {
295  //Set the temperature to -0.5 (cooled)
296  nod_pt->set_value(2,-0.5);
297  //Add small velocity imperfection if desired
298  if(Imperfect)
299  {
300  //Read out the x position
301  double x = nod_pt->x(0);
302  //Set a sinusoidal perturbation in the vertical velocity
303  //This perturbation is mass conserving
304  double value = sin(2.0*3.141592654*x/3.0);
305  nod_pt->set_value(1,value);
306  }
307  }
308 
309  //If we are on the bottom boundary, set the temperature
310  //to 0.5 (heated)
311  if(ibound==0) {nod_pt->set_value(2,0.5);}
312  }
313  }
314 
315 } // end of actions before solve
316 
317 
318 
319 //====================start_of_doc_solution===============================
320 /// Doc the solution
321 //========================================================================
322 template<class ELEMENT>
324 {
325  //Declare an output stream and filename
326  ofstream some_file;
327  char filename[100];
328 
329  // Number of plot points: npts x npts
330  unsigned npts=5;
331 
332  // Output solution
333  //-----------------
334  sprintf(filename,"%s/soln%i.dat",Doc_info.directory().c_str(),
335  Doc_info.number());
336  some_file.open(filename);
337  mesh_pt()->output(some_file,npts);
338  some_file.close();
339 
340  Doc_info.number()++;
341 } // end of doc
342 
343 
344 //===============start_of_main========================================
345 /// Driver code for 2D Boussinesq convection problem with
346 /// adaptivity.
347 //====================================================================
348 int main()
349 {
350 
351  // Set the direction of gravity
354 
355  // Create the problem with 2D nine-node refineable elements.
358 
359  // Apply a perturbation to the upper boundary condition to
360  // force the solution into the symmetry-broken state.
361  problem.enable_imperfection();
362 
363  //Solve the problem with (up to) two levels of adaptation
364  problem.newton_solve(2);
365 
366  //Document the solution
367  problem.doc_solution();
368 
369  // Make the boundary conditions perfect and solve again.
370  // Now the slightly perturbed symmetry broken state computed
371  // above is used as the initial condition and the Newton solver
372  // converges to the symmetry broken solution, even without
373  // the perturbation
374  problem.disable_imperfection();
375  problem.newton_solve(2);
376  problem.doc_solution();
377 
378 } // end of main
379 
380 
381 
382 
383 
384 
385 
386 
387 
~RefineableConvectionProblem()
Destructor. Empty.
void actions_after_adapt()
Actions after adaptation, Re-pin a single pressure degree of freedom.
double Peclet
Peclet number (identically one from our non-dimensionalisation)
void actions_after_newton_solve()
Update the problem after solve (empty)
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 actions_before_adapt()
Actions before adapt:(empty)
void actions_before_newton_solve()
Update the problem specs before solve:
RectangularQuadMesh< ELEMENT > * mesh_pt()
Overloaded version of the problem&#39;s access function to the mesh. Recasts the pointer to the base Mesh...
Namespace for the physical parameters in the problem.
Vector< double > Direction_of_gravity(2)
Gravity vector.
void enable_imperfection()
Set the boundary condition on the upper wall to be perturbed slightly to force the solution into the ...
double Inverse_Prandtl
1/Prandtl number
void disable_imperfection()
Set the boundary condition on the upper wall to be unperturbed.
double Rayleigh
Rayleigh number, set to be greater than the threshold for linear instability.