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