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