two_d_adv_diff_SUPG.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-2023 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 simple 2D adv diff problem with SUPG stabilisation
27 
28 //Generic routines
29 #include "generic.h"
30 
31 // The Poisson equations
32 #include "advection_diffusion.h"
33 
34 // The mesh
35 #include "meshes/rectangular_quadmesh.h"
36 
37 using namespace std;
38 
39 using namespace oomph;
40 
41 
42 
43 //======start_of_namespace============================================
44 /// Namespace for global parameters: Unforced problem with
45 /// boundary values corresponding to a steep tanh step profile
46 /// oriented at 45 degrees across the domain.
47 //====================================================================
49 {
50 
51  /// Peclet number
52  double Peclet=200.0;
53 
54  /// Parameter for steepness of step in boundary values
55  double Alpha=50.0;
56 
57  /// Parameter for angle of step in boundary values: 45 degrees
58  double TanPhi=1.0;
59 
60  /// Some "solution" for assignment of boundary values
61  void get_boundary_values(const Vector<double>& x, Vector<double>& u)
62  {
63  u[0]=tanh(1.0-Alpha*(TanPhi*x[0]-x[1]));
64  }
65 
66  /// Zero source function
67  void source_function(const Vector<double>& x_vect, double& source)
68  {
69  source=0.0;
70  }
71 
72  /// Wind
73  void wind_function(const Vector<double>& x, Vector<double>& wind)
74  {
75  wind[0]=sin(6.0*x[1]);
76  wind[1]=cos(6.0*x[0]);
77  }
78 
79 
80 } // end of namespace
81 
82 /// ///////////////////////////////////////////////////////////////////
83 /// ///////////////////////////////////////////////////////////////////
84 /// ///////////////////////////////////////////////////////////////////
85 
86 
87 //====== start_of_problem_class=======================================
88 /// 2D AdvectionDiffusion problem on rectangular domain, discretised
89 /// with refineable 2D QAdvectionDiffusion elements. The specific type
90 /// of element is specified via the template parameter.
91 //====================================================================
92 template<class ELEMENT>
93 class SUPGAdvectionDiffusionProblem : public Problem
94 {
95 
96 public:
97 
98  /// Constructor: Pass pointer to source and wind functions, and
99  /// flag to indicate if stabilisation is to be used.
101  AdvectionDiffusionEquations<2>::AdvectionDiffusionSourceFctPt source_fct_pt,
102  AdvectionDiffusionEquations<2>::AdvectionDiffusionWindFctPt wind_fct_pt,
103  const bool& use_stabilisation);
104 
105  /// Destructor. Empty
107 
108  /// Update the problem specs before solve: Reset boundary conditions
109  /// to the values from the tanh solution and compute stabilisation
110  /// parameter.
111  void actions_before_newton_solve();
112 
113  /// Update the problem after solve (empty)
115 
116  /// Doc the solution.
117  void doc_solution();
118 
119  /// Overloaded version of the problem's access function to
120  /// the mesh. Recasts the pointer to the base Mesh object to
121  /// the actual mesh type.
122  RectangularQuadMesh<ELEMENT>* mesh_pt()
123  {
124  return dynamic_cast<RectangularQuadMesh<ELEMENT>*>(
125  Problem::mesh_pt());
126  }
127 
128 private:
129 
130  /// DocInfo object
131  DocInfo Doc_info;
132 
133  /// Pointer to source function
134  AdvectionDiffusionEquations<2>::AdvectionDiffusionSourceFctPt Source_fct_pt;
135 
136  /// Pointer to wind function
137  AdvectionDiffusionEquations<2>::AdvectionDiffusionWindFctPt Wind_fct_pt;
138 
139  /// Flag to indicate if stabilisation is to be used
141 
142 }; // end of problem class
143 
144 
145 
146 //=====start_of_constructor===============================================
147 /// Constructor for AdvectionDiffusion problem: Pass pointer to
148 /// source function and wind functions and flag to indicate
149 /// if stabilisation is to be used.
150 //========================================================================
151 template<class ELEMENT>
153  AdvectionDiffusionEquations<2>::AdvectionDiffusionSourceFctPt source_fct_pt,
154  AdvectionDiffusionEquations<2>::AdvectionDiffusionWindFctPt wind_fct_pt,
155  const bool& use_stabilisation)
156  : Source_fct_pt(source_fct_pt), Wind_fct_pt(wind_fct_pt),
157  Use_stabilisation(use_stabilisation)
158 {
159 
160  // Set output directory
161  if (use_stabilisation)
162  {
163  Doc_info.set_directory("RESLT_stabilised");
164  }
165  else
166  {
167  Doc_info.set_directory("RESLT_unstabilised");
168  }
169 
170  // Setup mesh
171 
172  // # of elements in x-direction
173  unsigned n_x=40;
174 
175  // # of elements in y-direction
176  unsigned n_y=40;
177 
178  // Domain length in x-direction
179  double l_x=1.0;
180 
181  // Domain length in y-direction
182  double l_y=2.0;
183 
184  // Build and assign mesh
185  Problem::mesh_pt() =
186  new RectangularQuadMesh<ELEMENT>(n_x,n_y,l_x,l_y);
187 
188 
189  // Set the boundary conditions for this problem: All nodes are
190  // free by default -- only need to pin the ones that have Dirichlet
191  // conditions here
192  unsigned num_bound = mesh_pt()->nboundary();
193  for(unsigned ibound=0;ibound<num_bound;ibound++)
194  {
195  unsigned num_nod= mesh_pt()->nboundary_node(ibound);
196  for (unsigned inod=0;inod<num_nod;inod++)
197  {
198  mesh_pt()->boundary_node_pt(ibound,inod)->pin(0);
199  }
200  } // end loop over boundaries
201 
202  // Complete the build of all elements so they are fully functional
203 
204  // Loop over the elements to set up element-specific
205  // things that cannot be handled by the (argument-free!) ELEMENT
206  // constructor: Pass pointer to source function
207  unsigned n_element = mesh_pt()->nelement();
208  for(unsigned i=0;i<n_element;i++)
209  {
210  // Upcast from GeneralsedElement to the present element
211  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(mesh_pt()->element_pt(i));
212 
213  //Set the source function pointer
214  el_pt->source_fct_pt() = Source_fct_pt;
215 
216  //Set the wind function pointer
217  el_pt->wind_fct_pt() = Wind_fct_pt;
218 
219  // Set the Peclet number
220  el_pt->pe_pt() = &GlobalPhysicalParameters::Peclet;
221  }
222 
223  // Setup equation numbering scheme
224  cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
225 
226 } // end of constructor
227 
228 
229 //========================================start_of_actions_before_newton_solve===
230 /// Update the problem specs before solve: (Re-)set boundary conditions
231 //========================================================================
232 template<class ELEMENT>
234 {
235  // How many boundaries are there?
236  unsigned num_bound = mesh_pt()->nboundary();
237 
238  //Loop over the boundaries
239  for(unsigned ibound=0;ibound<num_bound;ibound++)
240  {
241  // How many nodes are there on this boundary?
242  unsigned num_nod=mesh_pt()->nboundary_node(ibound);
243 
244  // Loop over the nodes on boundary
245  for (unsigned inod=0;inod<num_nod;inod++)
246  {
247  // Get pointer to node
248  Node* nod_pt=mesh_pt()->boundary_node_pt(ibound,inod);
249 
250  // Extract nodal coordinates from node:
251  Vector<double> x(2);
252  x[0]=nod_pt->x(0);
253  x[1]=nod_pt->x(1);
254 
255  // Get boundary value
256  Vector<double> u(1);
258 
259  // Assign the value to the one (and only) nodal value at this node
260  nod_pt->set_value(0,u[0]);
261  }
262  }
263 
264  // Now loop over all elements and set the stabilisation parameter
265  unsigned n_element = mesh_pt()->nelement();
266  for(unsigned i=0;i<n_element;i++)
267  {
268  // Upcast from GeneralsedElement to the present element
269  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(mesh_pt()->element_pt(i));
270 
271 
272  // Use stabilisation?
273  if (Use_stabilisation)
274  {
275  //Compute stabilisation parameter
276  el_pt->compute_stabilisation_parameter();
277  }
278  else
279  {
280  //Compute stabilisation parameter
281  el_pt->switch_off_stabilisation();
282  }
283 
284  }
285 
286 } // end of actions before solve
287 
288 
289 
290 //===============start_of_doc=============================================
291 /// Doc the solution
292 //========================================================================
293 template<class ELEMENT>
295 {
296 
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 } // end of doc
312 
313 
314 //===== start_of_main=====================================================
315 /// Driver code for 2D AdvectionDiffusion problem
316 //========================================================================
317 int main()
318 {
319 
320  //Set up the problem with stabilisation
321  {
322  bool use_stabilisation=true;
323 
324  // Create the problem with 2D nine-node elements from the
325  // QAdvectionDiffusionElement family. Pass pointer to
326  // source and wind function.
330  use_stabilisation);
331 
332  // Solve the problem
333  problem.newton_solve();
334 
335  //Output the solution
336  problem.doc_solution();
337 
338  }
339 
340 
341 
342  //Set up the problem without stabilisation
343  {
344 
345  bool use_stabilisation=false;
346 
347  // Create the problem with 2D nine-node elements from the
348  // QAdvectionDiffusionElement family. Pass pointer to
349  // source and wind function.
353  use_stabilisation);
354 
355  // Solve the problem
356  problem.newton_solve();
357 
358  //Output the solution
359  problem.doc_solution();
360 
361  }
362 
363 
364 
365 } // end of main
366 
367 
368 
369 
370 
371 
372 
373 
374 
/////////////////////////////////////////////////////////////////// /////////////////////////////////...
~SUPGAdvectionDiffusionProblem()
Destructor. Empty.
void doc_solution()
Doc the solution.
void actions_after_newton_solve()
Update the problem after solve (empty)
AdvectionDiffusionEquations< 2 >::AdvectionDiffusionWindFctPt Wind_fct_pt
Pointer to wind function.
SUPGAdvectionDiffusionProblem(AdvectionDiffusionEquations< 2 >::AdvectionDiffusionSourceFctPt source_fct_pt, AdvectionDiffusionEquations< 2 >::AdvectionDiffusionWindFctPt wind_fct_pt, const bool &use_stabilisation)
Constructor: Pass pointer to source and wind functions, and flag to indicate if stabilisation is to b...
AdvectionDiffusionEquations< 2 >::AdvectionDiffusionSourceFctPt Source_fct_pt
Pointer to source function.
bool Use_stabilisation
Flag to indicate if stabilisation is to be used.
RectangularQuadMesh< ELEMENT > * mesh_pt()
Overloaded version of the problem's access function to the mesh. Recasts the pointer to the base Mesh...
void actions_before_newton_solve()
Update the problem specs before solve: Reset boundary conditions to the values from the tanh solution...
Namespace for global parameters: Unforced problem with boundary values corresponding to a steep tanh ...
void wind_function(const Vector< double > &x, Vector< double > &wind)
Wind.
void get_boundary_values(const Vector< double > &x, Vector< double > &u)
Some "solution" for assignment of boundary values.
void source_function(const Vector< double > &x_vect, double &source)
Zero source function.
double Peclet
Peclet number.
double Alpha
Parameter for steepness of step in boundary values.
double TanPhi
Parameter for angle of step in boundary values: 45 degrees.
int main()
Driver code for 2D AdvectionDiffusion problem.