two_d_poisson2.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 simple 2D poisson problem
31 
32 //Generic routines
33 #include "generic.h"
34 
35 // The Poisson equations
36 #include "poisson.h"
37 
38 // Include the templated mesh header only -- the mesh is
39 // precompiled and instantiated with the required element type
40 // in the separate file, two_d_poisson2_mesh.cc to avoid
41 // recompilation.
42 #include "meshes/simple_rectangular_quadmesh.template.h"
43 
44 using namespace std;
45 
46 using namespace oomph;
47 
48 //===== start_of_namespace=============================================
49 /// Namespace for exact solution for Poisson equation with "sharp step"
50 //=====================================================================
51 namespace TanhSolnForPoisson
52 {
53 
54  /// Parameter for steepness of "step"
55  double Alpha=1.0;
56 
57  /// Parameter for angle Phi of "step"
58  double TanPhi=0.0;
59 
60  /// Exact solution as a Vector
61  void get_exact_u(const Vector<double>& x, Vector<double>& u)
62  {
63  u[0]=tanh(1.0-Alpha*(TanPhi*x[0]-x[1]));
64  }
65 
66  /// Source function required to make the solution above an exact solution
67  void get_source(const Vector<double>& x, double& source)
68  {
69  source = 2.0*tanh(-1.0+Alpha*(TanPhi*x[0]-x[1]))*
70  (1.0-pow(tanh(-1.0+Alpha*(TanPhi*x[0]-x[1])),2.0))*
71  Alpha*Alpha*TanPhi*TanPhi+2.0*tanh(-1.0+Alpha*(TanPhi*x[0]-x[1]))*
72  (1.0-pow(tanh(-1.0+Alpha*(TanPhi*x[0]-x[1])),2.0))*Alpha*Alpha;
73  }
74 
75 } // end of namespace
76 
77 
78 
79 
80 
81 
82 
83 
84 //====== start_of_problem_class=======================================
85 /// 2D Poisson problem on rectangular domain, discretised with
86 /// 2D QPoisson elements. The specific type of element is
87 /// specified via the template parameter.
88 //====================================================================
89 template<class ELEMENT>
90 class PoissonProblem : public Problem
91 {
92 
93 public:
94 
95  /// Constructor: Pass pointer to source function
96  PoissonProblem(PoissonEquations<2>::PoissonSourceFctPt source_fct_pt);
97 
98  /// Destructor (empty)
100 
101  /// \short Update the problem specs before solve: Reset boundary conditions
102  /// to the values from the exact solution.
103  void actions_before_newton_solve();
104 
105  /// Update the problem after solve (empty)
107 
108  /// \short Doc the solution. DocInfo object stores flags/labels for where the
109  /// output gets written to
110  void doc_solution(DocInfo& doc_info);
111 
112 private:
113 
114  /// Pointer to source function
115  PoissonEquations<2>::PoissonSourceFctPt Source_fct_pt;
116 
117 }; // end of problem class
118 
119 
120 
121 
122 //=====start_of_constructor===============================================
123 /// Constructor for Poisson problem: Pass pointer to source function.
124 //========================================================================
125 template<class ELEMENT>
127  PoissonProblem(PoissonEquations<2>::PoissonSourceFctPt source_fct_pt)
128  : Source_fct_pt(source_fct_pt)
129 {
130 
131  // Setup mesh
132 
133  // # of elements in x-direction
134  unsigned n_x=4;
135 
136  // # of elements in y-direction
137  unsigned n_y=4;
138 
139  // Domain length in x-direction
140  double l_x=1.0;
141 
142  // Domain length in y-direction
143  double l_y=2.0;
144 
145  // Build and assign mesh
146  Problem::mesh_pt() = new SimpleRectangularQuadMesh<ELEMENT>(n_x,n_y,l_x,l_y);
147 
148  // Set the boundary conditions for this problem: All nodes are
149  // free by default -- only need to pin the ones that have Dirichlet conditions
150  // here.
151  unsigned num_bound = mesh_pt()->nboundary();
152  for(unsigned ibound=0;ibound<num_bound;ibound++)
153  {
154  unsigned num_nod= mesh_pt()->nboundary_node(ibound);
155  for (unsigned inod=0;inod<num_nod;inod++)
156  {
157  mesh_pt()->boundary_node_pt(ibound,inod)->pin(0);
158  }
159  }
160 
161  // Complete the build of all elements so they are fully functional
162 
163  // Loop over the elements to set up element-specific
164  // things that cannot be handled by the (argument-free!) ELEMENT
165  // constructor: Pass pointer to source function
166  unsigned n_element = mesh_pt()->nelement();
167  for(unsigned i=0;i<n_element;i++)
168  {
169  // Upcast from GeneralsedElement to the present element
170  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(mesh_pt()->element_pt(i));
171 
172  //Set the source function pointer
173  el_pt->source_fct_pt() = Source_fct_pt;
174  }
175 
176 
177  // Setup equation numbering scheme
178  cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
179 
180 } // end of constructor
181 
182 
183 
184 
185 //========================================start_of_actions_before_newton_solve===
186 /// Update the problem specs before solve: (Re-)set boundary conditions
187 /// to the values from the exact solution.
188 //========================================================================
189 template<class ELEMENT>
191 {
192  // How many boundaries are there?
193  unsigned num_bound = mesh_pt()->nboundary();
194 
195  //Loop over the boundaries
196  for(unsigned ibound=0;ibound<num_bound;ibound++)
197  {
198  // How many nodes are there on this boundary?
199  unsigned num_nod=mesh_pt()->nboundary_node(ibound);
200 
201  // Loop over the nodes on boundary
202  for (unsigned inod=0;inod<num_nod;inod++)
203  {
204  // Get pointer to node
205  Node* nod_pt=mesh_pt()->boundary_node_pt(ibound,inod);
206 
207  // Extract nodal coordinates from node:
208  Vector<double> x(2);
209  x[0]=nod_pt->x(0);
210  x[1]=nod_pt->x(1);
211 
212  // Compute the value of the exact solution at the nodal point
213  Vector<double> u(1);
215 
216  // Assign the value to the one (and only) nodal value at this node
217  nod_pt->set_value(0,u[0]);
218  }
219  }
220 } // end of actions before solve
221 
222 
223 
224 //===============start_of_doc=============================================
225 /// Doc the solution: doc_info contains labels/output directory etc.
226 //========================================================================
227 template<class ELEMENT>
228 void PoissonProblem<ELEMENT>::doc_solution(DocInfo& doc_info)
229 {
230 
231  ofstream some_file;
232  char filename[100];
233 
234  // Number of plot points: npts x npts
235  unsigned npts=5;
236 
237  // Output solution
238  //-----------------
239  sprintf(filename,"%s/soln%i.dat",doc_info.directory().c_str(),
240  doc_info.number());
241  some_file.open(filename);
242  mesh_pt()->output(some_file,npts);
243  some_file.close();
244 
245 
246  // Output exact solution
247  //----------------------
248  sprintf(filename,"%s/exact_soln%i.dat",doc_info.directory().c_str(),
249  doc_info.number());
250  some_file.open(filename);
251  mesh_pt()->output_fct(some_file,npts,TanhSolnForPoisson::get_exact_u);
252  some_file.close();
253 
254  // Doc error and return of the square of the L2 error
255  //---------------------------------------------------
256  double error,norm;
257  sprintf(filename,"%s/error%i.dat",doc_info.directory().c_str(),
258  doc_info.number());
259  some_file.open(filename);
260  mesh_pt()->compute_error(some_file,TanhSolnForPoisson::get_exact_u,
261  error,norm);
262  some_file.close();
263 
264  // Doc L2 error and norm of solution
265  cout << "\nNorm of error : " << sqrt(error) << std::endl;
266  cout << "Norm of solution: " << sqrt(norm) << std::endl << std::endl;
267 
268 } // end of doc
269 
270 
271 
272 
273 
274 
275 //===== start_of_main=====================================================
276 /// Driver code for 2D Poisson problem
277 //========================================================================
278 int main()
279 {
280 
281  //Set up the problem
282  //------------------
283 
284  // Create the problem with 2D nine-node elements from the
285  // QPoissonElement family. Pass pointer to source function.
288 
289  // Create label for output
290  //------------------------
291  DocInfo doc_info;
292 
293  // Set output directory
294  doc_info.set_directory("RESLT");
295 
296  // Step number
297  doc_info.number()=0;
298 
299  // Check if we're ready to go:
300  //----------------------------
301  cout << "\n\n\nProblem self-test ";
302  if (problem.self_test()==0)
303  {
304  cout << "passed: Problem can be solved." << std::endl;
305  }
306  else
307  {
308  throw OomphLibError("Self test failed",
309  OOMPH_CURRENT_FUNCTION,
310  OOMPH_EXCEPTION_LOCATION);
311  }
312 
313 
314  // Set the orientation of the "step" to 45 degrees
316 
317  // Initial value for the steepness of the "step"
319 
320 
321  // Do a couple of solutions for different forcing functions
322  //---------------------------------------------------------
323  unsigned nstep=4;
324  for (unsigned istep=0;istep<nstep;istep++)
325  {
326 
327  // Increase the steepness of the step:
329 
330  cout << "\n\nSolving for TanhSolnForPoisson::Alpha="
331  << TanhSolnForPoisson::Alpha << std::endl << std::endl;
332 
333  // Solve the problem
334  problem.newton_solve();
335 
336  //Output the solution
337  problem.doc_solution(doc_info);
338 
339  //Increment counter for solutions
340  doc_info.number()++;
341 
342  }
343 
344 } //end of main
345 
346 
347 
348 
349 
350 
351 
352 
353 
~PoissonProblem()
Destructor (empty)
void actions_after_newton_solve()
Update the problem after solve (empty)
void doc_solution(DocInfo &doc_info)
Doc the solution. DocInfo object stores flags/labels for where the output gets written to...
void get_exact_u(const Vector< double > &x, Vector< double > &u)
Exact solution as a Vector.
double Alpha
Parameter for steepness of "step".
void actions_before_newton_solve()
Update the problem specs before solve: Reset boundary conditions to the values from the exact solutio...
double TanPhi
Parameter for angle Phi of "step".
Namespace for exact solution for Poisson equation with "sharp step".
PoissonProblem(PoissonEquations< 2 >::PoissonSourceFctPt source_fct_pt)
Constructor: Pass pointer to source function.
void get_source(const Vector< double > &x, double &source)
Source function required to make the solution above an exact solution.
int main()
Driver code for 2D Poisson problem.