two_d_adv_diff_adapt2.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 adaptive mesh refinement
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 /// Namespace for physical parameters and boundary conditions
43 /// "sharp" step
44 //====================================================================
46 {
47 
48  /// Peclet number
49  double Peclet=200.0;
50 
51  /// Parameter for steepness of step in tanh profile
52  double Alpha;
53 
54  /// Parameter for angle of step in tanh profile
55  double TanPhi;
56 
57  /// Tanh profile for assignment of boundary conditons as a Vector
58  void tanh_profile(const Vector<double>& x, Vector<double>& u)
59  {
60  u[0]=tanh(1.0-Alpha*(TanPhi*x[0]-x[1]));
61  }
62 
63  /// Tanh profile for assignment of boundary conditons as a Vector
64  void tanh_profile(const Vector<double>& x, double& u)
65  {
66  u=tanh(1.0-Alpha*(TanPhi*x[0]-x[1]));
67  }
68 
69  /// Zero source function
70  void source_function(const Vector<double>& x_vect, double& source)
71  {
72  //double x=x_vect[0];
73  //double y=x_vect[1];
74  source = 0.0;
75  }
76 
77  /// Wind: Recirculating cells
78  void wind_function(const Vector<double>& x, Vector<double>& wind)
79  {
80  wind[0]=sin(6.0*x[1]);
81  wind[1]=cos(6.0*x[0]);
82  }
83 
84 } // end of namespace
85 
86 /// ///////////////////////////////////////////////////////////////////
87 /// ///////////////////////////////////////////////////////////////////
88 /// ///////////////////////////////////////////////////////////////////
89 
90 //====== start_of_problem_class=======================================
91 /// 2D AdvectionDiffusion problem on rectangular domain, discretised
92 /// with refineable 2D QAdvectionDiffusion elements. The specific
93 /// type of element is specified via the template parameter.
94 //====================================================================
95 template<class ELEMENT>
96 class RefineableAdvectionDiffusionProblem : public Problem
97 {
98 
99 public:
100 
101  /// Constructor: Pass pointer to source function
103  AdvectionDiffusionEquations<2>::AdvectionDiffusionSourceFctPt source_fct_pt,
104  AdvectionDiffusionEquations<2>::AdvectionDiffusionWindFctPt wind_fct_pt);
105 
106  /// Destructor (empty)
108 
109  /// Update the problem specs before solve: Reset boundary conditions
110  /// to the values from the exact solution.
112 
113  /// Update the problem after solve (empty)
115 
116  /// Doc the solution. DocInfo object stores flags/labels for where the
117  /// output gets written to
118  void doc_solution(DocInfo& doc_info);
119 
120  /// Overloaded version of the problem's access function to
121  /// the mesh. Recasts the pointer to the base Mesh object to
122  /// the actual mesh type.
123  RefineableRectangularQuadMesh<ELEMENT>* mesh_pt()
124  {
125  return dynamic_cast<RefineableRectangularQuadMesh<ELEMENT>*>(
126  Problem::mesh_pt());
127  }
128 
129 private:
130 
131  /// Pointer to source function
132  AdvectionDiffusionEquations<2>::AdvectionDiffusionSourceFctPt Source_fct_pt;
133 
134  /// Pointer to wind function
135  AdvectionDiffusionEquations<2>::AdvectionDiffusionWindFctPt Wind_fct_pt;
136 
137 }; // end of problem class
138 
139 
140 
141 
142 //=====start_of_constructor===============================================
143 /// Constructor for AdvectionDiffusion problem: Pass pointer to source
144 /// function.
145 //========================================================================
146 template<class ELEMENT>
148  AdvectionDiffusionEquations<2>::AdvectionDiffusionSourceFctPt source_fct_pt,
149  AdvectionDiffusionEquations<2>::AdvectionDiffusionWindFctPt wind_fct_pt)
150  : Source_fct_pt(source_fct_pt), Wind_fct_pt(wind_fct_pt)
151 {
152 
153  // Setup mesh
154 
155  // # of elements in x-direction
156  unsigned n_x=4;
157 
158  // # of elements in y-direction
159  unsigned n_y=4;
160 
161  // Domain length in x-direction
162  double l_x=1.0;
163 
164  // Domain length in y-direction
165  double l_y=2.0;
166 
167  // Build and assign mesh
168  Problem::mesh_pt() =
169  new RefineableRectangularQuadMesh<ELEMENT>(n_x,n_y,l_x,l_y);
170 
171  // Create/set error estimator
172  mesh_pt()->spatial_error_estimator_pt()=new Z2ErrorEstimator;
173 
174  // Set the boundary conditions for this problem: All nodes are
175  // free by default -- only need to pin the ones that have Dirichlet conditions
176  // here.
177  unsigned num_bound = mesh_pt()->nboundary();
178  for(unsigned ibound=0;ibound<num_bound;ibound++)
179  {
180  unsigned num_nod= mesh_pt()->nboundary_node(ibound);
181  for (unsigned inod=0;inod<num_nod;inod++)
182  {
183  mesh_pt()->boundary_node_pt(ibound,inod)->pin(0);
184  }
185  }
186 
187  // Complete the build of all elements so they are fully functional
188 
189  // Loop over the elements to set up element-specific
190  // things that cannot be handled by the (argument-free!) ELEMENT
191  // constructor: Pass pointer to source function
192  unsigned n_element = mesh_pt()->nelement();
193  for(unsigned i=0;i<n_element;i++)
194  {
195  // Upcast from GeneralsedElement to the present element
196  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(mesh_pt()->element_pt(i));
197 
198  //Set the source function pointer
199  el_pt->source_fct_pt() = Source_fct_pt;
200 
201  //Set the wind function pointer
202  el_pt->wind_fct_pt() = Wind_fct_pt;
203 
204  // Set the Peclet number
205  el_pt->pe_pt() = &TanhSolnForAdvectionDiffusion::Peclet;
206  }
207 
208  // Setup equation numbering scheme
209  cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
210 
211 } // end of constructor
212 
213 //========================================start_of_actions_before_newton_solve===
214 /// Update the problem specs before solve: (Re-)set boundary conditions
215 /// to the values from the exact solution.
216 //========================================================================
217 template<class ELEMENT>
219 {
220  // How many boundaries are there?
221  unsigned num_bound = mesh_pt()->nboundary();
222 
223  //Loop over the boundaries
224  for(unsigned ibound=0;ibound<num_bound;ibound++)
225  {
226  // How many nodes are there on this boundary?
227  unsigned num_nod=mesh_pt()->nboundary_node(ibound);
228 
229  // Loop over the nodes on boundary
230  for (unsigned inod=0;inod<num_nod;inod++)
231  {
232  // Get pointer to node
233  Node* nod_pt=mesh_pt()->boundary_node_pt(ibound,inod);
234 
235  // Extract nodal coordinates from node:
236  Vector<double> x(2);
237  x[0]=nod_pt->x(0);
238  x[1]=nod_pt->x(1);
239 
240  // Compute the value of the exact solution at the nodal point
241  Vector<double> u(1);
243 
244  // Assign the value to the one (and only) nodal value at this node
245  nod_pt->set_value(0,u[0]);
246  }
247  }
248 } // end of actions before solve
249 
250 //===============start_of_doc=============================================
251 /// Doc the solution: doc_info contains labels/output directory etc.
252 //========================================================================
253 template<class ELEMENT>
255  DocInfo& doc_info)
256 {
257 
258  ofstream some_file;
259  char filename[100];
260 
261  // Number of plot points: npts x npts
262  unsigned npts=5;
263 
264  // Output solution
265  //-----------------
266  sprintf(filename,"%s/soln%i.dat",doc_info.directory().c_str(),
267  doc_info.number());
268  some_file.open(filename);
269  mesh_pt()->output(some_file,npts);
270  some_file.close();
271 
272 } // end of doc
273 
274 //===== start_of_main=====================================================
275 /// Driver code for 2D AdvectionDiffusion problem
276 //========================================================================
277 int main()
278 {
279 
280  //Set up the problem
281  //------------------
282 
283  // Create the problem with 2D nine-node refineable elements from the
284  // RefineableQuadAdvectionDiffusionElement family. Pass pointer to
285  // source and wind function.
286  RefineableAdvectionDiffusionProblem<RefineableQAdvectionDiffusionElement<2,
289 
290  // Create label for output
291  //------------------------
292  DocInfo doc_info;
293 
294  // Set output directory
295  doc_info.set_directory("RESLT");
296 
297  // Step number
298  doc_info.number()=0;
299 
300  // Check if we're ready to go:
301  //----------------------------
302  cout << "\n\n\nProblem self-test ";
303  if (problem.self_test()==0)
304  {
305  cout << "passed: Problem can be solved." << std::endl;
306  }
307  else
308  {
309  throw OomphLibError("Self test failed",
310  OOMPH_CURRENT_FUNCTION,
311  OOMPH_EXCEPTION_LOCATION);
312  }
313 
314  // Set the orientation of the "step" to 45 degrees
316 
317  // Choose a large value for the steepness of the "step"
319 
320  // Doc (default) refinement targets
321  //----------------------------------
322  problem.mesh_pt()->doc_adaptivity_targets(cout);
323 
324 
325 
326  // Solve/doc the problem on the initial, very coarse mesh
327  //-------------------------------------------------------
328 
329  // Solve the problem
330  problem.newton_solve();
331 
332  //Output solution
333  problem.doc_solution(doc_info);
334 
335  //Increment counter for solutions
336  doc_info.number()++; // end of Solve/doc very coarse mesh
337 
338  // Now do (up to) four rounds of fully automatic adapation in response to
339  //-----------------------------------------------------------------------
340  // error estimate
341  //---------------
342  unsigned max_solve=4;
343  for (unsigned isolve=0;isolve<max_solve;isolve++)
344  {
345  // Adapt problem/mesh
346  problem.adapt();
347 
348  // Re-solve the problem if the adaptation has changed anything
349  if ((problem.mesh_pt()->nrefined() !=0)||
350  (problem.mesh_pt()->nunrefined()!=0))
351  {
352  problem.newton_solve();
353  }
354  else
355  {
356  cout << "Mesh wasn't adapted --> we'll stop here" << std::endl;
357  break;
358  }
359 
360  //Output solution
361  problem.doc_solution(doc_info);
362 
363  //Increment counter for solutions
364  doc_info.number()++;
365  } // end (up to) four rounds of fully automatic adapation
366 
367 
368 } //end of main
369 
370 
371 
372 
373 
374 
375 
376 
377 
/////////////////////////////////////////////////////////////////// /////////////////////////////////...
RefineableAdvectionDiffusionProblem(AdvectionDiffusionEquations< 2 >::AdvectionDiffusionSourceFctPt source_fct_pt, AdvectionDiffusionEquations< 2 >::AdvectionDiffusionWindFctPt wind_fct_pt)
Constructor: Pass pointer to source function.
void actions_before_newton_solve()
Update the problem specs before solve: Reset boundary conditions to the values from the exact solutio...
AdvectionDiffusionEquations< 2 >::AdvectionDiffusionSourceFctPt Source_fct_pt
Pointer to source function.
void actions_after_newton_solve()
Update the problem after solve (empty)
AdvectionDiffusionEquations< 2 >::AdvectionDiffusionWindFctPt Wind_fct_pt
Pointer to wind function.
RefineableRectangularQuadMesh< ELEMENT > * mesh_pt()
Overloaded version of the problem's access function to the mesh. Recasts the pointer to the base Mesh...
~RefineableAdvectionDiffusionProblem()
Destructor (empty)
Namespace for exact solution for AdvectionDiffusion equation with "sharp" step.
double TanPhi
Parameter for angle of step.
void tanh_profile(const Vector< double > &x, Vector< double > &u)
Tanh profile for assignment of boundary conditons as a Vector.
double Alpha
Parameter for steepness of step.
void source_function(const Vector< double > &x_vect, double &source)
Source function required to make the solution above an exact solution.
void wind_function(const Vector< double > &x, Vector< double > &wind)
Wind.
void tanh_profile(const Vector< double > &x, double &u)
Tanh profile for assignment of boundary conditons as a Vector.
int main()
Driver code for 2D AdvectionDiffusion problem.