mesh_from_triangle_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 code for a simple test poisson problem using a mesh
31 // generated from an input file generated by the triangle mesh generator
32 // Triangle.
33 
34 //Generic routines
35 #include "generic.h"
36 
37 // Poisson equations
38 #include "poisson.h"
39 
40 // The mesh
41 #include "meshes/triangle_mesh.h"
42 
43 using namespace std;
44 
45 using namespace oomph;
46 
47 //====================================================================
48 /// Namespace for exact solution for Poisson equation with sharp step
49 //====================================================================
51 {
52 
53  /// Parameter for steepness of step
54  double Alpha;
55 
56  /// Parameter for angle of step
57  double Beta;
58 
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*(Beta*x[0]-x[1]));
64  }
65 
66 
67  /// Exact solution as a scalar
68  void get_exact_u(const Vector<double>& x, double& u)
69  {
70  u=tanh(1.0-Alpha*(Beta*x[0]-x[1]));
71  }
72 
73 
74  /// Source function to make it an exact solution
75  void get_source(const Vector<double>& x, double& source)
76  {
77  source = 2.0*tanh(-1.0+Alpha*(Beta*x[0]-x[1]))*
78  (1.0-pow(tanh(-1.0+Alpha*(Beta*x[0]-x[1])),2.0))*
79  Alpha*Alpha*Beta*Beta+2.0*tanh(-1.0+Alpha*(Beta*x[0]-x[1]))*
80  (1.0-pow(tanh(-1.0+Alpha*(Beta*x[0]-x[1])),2.0))*Alpha*Alpha;
81  }
82 
83 }
84 
85 
86 
87 
88 
89 
90 
91 //====================================================================
92 /// Micky mouse Poisson problem.
93 //====================================================================
94 
95 // Poisson problem
96 template<class ELEMENT>
97 class PoissonProblem : public Problem
98 {
99 
100 public:
101 
102 
103  /// \short Constructor: Pass pointer to source function and names of
104  /// two triangle input files
105  PoissonProblem(PoissonEquations<2>::PoissonSourceFctPt source_fct_pt,
106  const string& node_file_name,
107  const string& element_file_name,
108  const string& poly_file_name);
109 
110  /// Destructor (empty)
112 
113  /// Update the problem specs before solve: (Re)set boundary conditions
115  {
116  //Loop over the boundaries
117  unsigned num_bound = mesh_pt()->nboundary();
118  for(unsigned ibound=0;ibound<num_bound;ibound++)
119  {
120  // Loop over the nodes on boundary
121  unsigned num_nod=mesh_pt()->nboundary_node(ibound);
122  for (unsigned inod=0;inod<num_nod;inod++)
123  {
124  Node* nod_pt=mesh_pt()->boundary_node_pt(ibound,inod);
125  double u;
126  Vector<double> x(2);
127  x[0]=nod_pt->x(0);
128  x[1]=nod_pt->x(1);
130  nod_pt->set_value(0,u);
131  }
132  }
133  }
134 
135  /// Update the problem specs before solve (empty)
137  {}
138 
139 #ifdef ADAPT
140  /// Actions performed after the adaptation steps
141  void actions_after_adapt();
142 #endif
143 
144 #ifdef ADAPT
145  /// Access function for the specific mesh
147  {
148  return dynamic_cast<RefineableTriangleMesh<ELEMENT>*>(Problem::mesh_pt());
149  }
150 #else
151  /// Access function for the specific mesh
153  {
154  return dynamic_cast<TriangleMesh<ELEMENT>*>(Problem::mesh_pt());
155  }
156 #endif
157 
158  /// Doc the solution
159  void doc_solution(DocInfo& doc_info);
160 
161 private:
162 
163  /// Pointer to source function
164  PoissonEquations<2>::PoissonSourceFctPt Source_fct_pt;
165 
166 #ifdef ADAPT
167  /// Pointer to the bulk mesh
169 
170  /// Error estimator
171  Z2ErrorEstimator* error_estimator_pt;
172 #endif
173 
174 };
175 
176 
177 
178 //========================================================================
179 /// Constructor for Poisson problem
180 //========================================================================
181 template<class ELEMENT>
183  PoissonProblem(PoissonEquations<2>::PoissonSourceFctPt source_fct_pt,
184  const string& node_file_name,
185  const string& element_file_name,
186  const string& poly_file_name)
187  : Source_fct_pt(source_fct_pt)
188 {
189 
190  // Setup parameters for exact tanh solution
191 
192  // Steepness of step
194 
195  // Orientation of step
197 
198 #ifdef ADAPT
199  //Create mesh
200  Bulk_mesh_pt = new RefineableTriangleMesh<ELEMENT>(node_file_name,
201  element_file_name,
202  poly_file_name);
203 
204  // Set error estimator for bulk mesh
205  Z2ErrorEstimator* error_estimator_pt=new Z2ErrorEstimator;
206  Bulk_mesh_pt->spatial_error_estimator_pt() = error_estimator_pt;
207 
208  // Set the problem mesh pointer to the bulk mesh
209  Problem::mesh_pt() = Bulk_mesh_pt;
210 
211 #else
212  //Create mesh
213  Problem::mesh_pt() = new TriangleMesh<ELEMENT>(node_file_name,
214  element_file_name,
215  poly_file_name);
216 #endif
217 
218  // Set the boundary conditions for this problem: All nodes are
219  // free by default -- just pin the ones that have Dirichlet conditions
220  // here.
221  unsigned num_bound = mesh_pt()->nboundary();
222  for(unsigned ibound=0;ibound<num_bound;ibound++)
223  {
224  unsigned num_nod= mesh_pt()->nboundary_node(ibound);
225  for (unsigned inod=0;inod<num_nod;inod++)
226  {
227  mesh_pt()->boundary_node_pt(ibound,inod)->pin(0);
228  }
229  }
230 
231  // Complete the build of all elements so they are fully functional
232 
233  //Find number of elements in mesh
234  unsigned n_element = mesh_pt()->nelement();
235 
236  // Loop over the elements to set up element-specific
237  // things that cannot be handled by constructor
238  for(unsigned i=0;i<n_element;i++)
239  {
240  // Upcast from GeneralElement to the present element
241  ELEMENT* el_pt = dynamic_cast<ELEMENT*>(mesh_pt()->element_pt(i));
242 
243  //Set the source function pointer
244  el_pt->source_fct_pt() = Source_fct_pt;
245  }
246 
247  // Setup equation numbering scheme
248  cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
249 
250 }
251 
252 #ifdef ADAPT
253 //========================================================================
254 /// Actions performed after the adaptation steps
255 //========================================================================
256 template<class ELEMENT>
258 {
259  // Since the mesh adaptation strategy is based on mesh re-generation
260  // we need to pin the nodes in the new regenerated mesh, and pass the
261  // source function to the elements
262 
263  // Set the boundary conditions for this problem: All nodes are free
264  // by default -- just pin the ones that have Dirichlet conditions
265  // here.
266  unsigned num_bound = mesh_pt()->nboundary();
267  for(unsigned ibound=0;ibound<num_bound;ibound++)
268  {
269  unsigned num_nod= mesh_pt()->nboundary_node(ibound);
270  for (unsigned inod=0;inod<num_nod;inod++)
271  {
272  mesh_pt()->boundary_node_pt(ibound,inod)->pin(0);
273  }
274  }
275 
276  // Complete the build of all elements so they are fully functional
277 
278  //Find number of elements in mesh
279  unsigned n_element = mesh_pt()->nelement();
280 
281  // Loop over the elements to set up element-specific
282  // things that cannot be handled by constructor
283  for(unsigned i=0;i<n_element;i++)
284  {
285  // Upcast from GeneralElement to the present element
286  ELEMENT* el_pt = dynamic_cast<ELEMENT*>(mesh_pt()->element_pt(i));
287 
288  //Set the source function pointer
289  el_pt->source_fct_pt() = Source_fct_pt;
290  }
291 
292  // Setup equation numbering scheme
293  cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
294 }
295 #endif
296 
297 
298 //========================================================================
299 /// Doc the solution
300 //========================================================================
301 template<class ELEMENT>
303 {
304 
305  ofstream some_file;
306  char filename[100];
307 
308  // Number of plot points
309  unsigned npts;
310  npts=5;
311 
312 
313 
314  // Output boundaries
315  //------------------
316  sprintf(filename,"%s/boundaries.dat",doc_info.directory().c_str());
317  some_file.open(filename);
318  mesh_pt()->output_boundaries(some_file);
319  some_file.close();
320 
321 
322  // Output solution
323  //----------------
324  sprintf(filename,"%s/soln%i.dat",doc_info.directory().c_str(),
325  doc_info.number());
326  some_file.open(filename);
327  mesh_pt()->output(some_file,npts);
328  some_file.close();
329 
330 
331  // Output exact solution
332  //----------------------
333  sprintf(filename,"%s/exact_soln%i.dat",doc_info.directory().c_str(),
334  doc_info.number());
335  some_file.open(filename);
336  mesh_pt()->output_fct(some_file,npts,TanhSolnForPoisson::get_exact_u);
337  some_file.close();
338 
339 
340  // Doc error
341  //----------
342  double error,norm;
343  sprintf(filename,"%s/error%i.dat",doc_info.directory().c_str(),
344  doc_info.number());
345  some_file.open(filename);
346  mesh_pt()->compute_error(some_file,TanhSolnForPoisson::get_exact_u,
347  error,norm);
348  some_file.close();
349  cout << "error: " << sqrt(error) << std::endl;
350  cout << "norm : " << sqrt(norm) << std::endl << std::endl;
351 
352 
353  // Get norm of solution
354  sprintf(filename,"%s/norm%i.dat",doc_info.directory().c_str(),
355  doc_info.number());
356  some_file.open(filename);
357  double norm_soln=0.0;
358  mesh_pt()->compute_norm(norm_soln);
359  some_file << sqrt(norm_soln) << std::endl;
360  cout << "Norm of computed solution: " << sqrt(norm_soln) << endl;
361  some_file.close();
362 
363 }
364 
365 
366 
367 
368 
369 //========================================================================
370 /// Demonstrate how to solve Poisson problem
371 //========================================================================
372 int main(int argc, char* argv[])
373 {
374  // Store command line arguments
375  CommandLineArgs::setup(argc,argv);
376 
377  // Check number of command line arguments: Need exactly two.
378  if (argc!=4)
379  {
380  std::string error_message =
381  "Wrong number of command line arguments.\n";
382  error_message +=
383  "Must specify the following file names \n";
384  error_message +=
385  "filename.node then filename.ele then filename.poly\n";
386 
387  throw OomphLibError(error_message,
388  OOMPH_CURRENT_FUNCTION,
389  OOMPH_EXCEPTION_LOCATION);
390  }
391 
392 #ifdef ADAPT
393  const unsigned max_adapt = 3;
394 #endif
395 
396  // Convert arguments to strings that specify the input file names
397  string node_file_name(argv[1]);
398  string element_file_name(argv[2]);
399  string poly_file_name(argv[3]);
400 
401  // Label for output
402  DocInfo doc_info;
403 
404  // Output directory
405  doc_info.set_directory("RESLT");
406 
407 #ifndef ADAPT
408  // Do the problem with cubic elements
409  //-----------------------------------
410  {
411  cout << std::endl << "Cubic elements" << std::endl;
412  cout << "==============" << std::endl << std::endl;
413 
414  //Set up the problem
416  problem(&TanhSolnForPoisson::get_source,node_file_name,
417  element_file_name,
418  poly_file_name);
419 
420  // Solve the problem
421  problem.newton_solve();
422 
423  //Output solution
424  problem.doc_solution(doc_info);
425 
426  //Increment counter for solutions
427  doc_info.number()++;
428  }
429 #endif
430 
431 
432  // Do the problem with quadratic elements
433  //---------------------------------------
434  {
435  cout << std::endl << "Quadratic elements" << std::endl;
436  cout << "===================" << std::endl << std::endl;
437 
438 #ifdef ADAPT
439  //Set up the problem
441  problem(&TanhSolnForPoisson::get_source,node_file_name,
442  element_file_name,
443  poly_file_name);
444 #else
445  //Set up the problem
447  problem(&TanhSolnForPoisson::get_source,node_file_name,
448  element_file_name,
449  poly_file_name);
450 #endif
451 
452 #ifdef ADAPT
453  // Solve the problem with adaptation
454  problem.newton_solve(max_adapt);
455 #else
456  // Solve the problem
457  problem.newton_solve();
458 #endif
459 
460  //Output solution
461  problem.doc_solution(doc_info);
462 
463  //Increment counter for solutions
464  doc_info.number()++;
465  }
466 
467 
468 
469  // Do the problem with linear elements
470  //------------------------------------
471  {
472  cout << std::endl << "Linear elements" << std::endl;
473  cout << "===============" << std::endl << std::endl;
474 
475 #ifdef ADAPT
476  //Set up the problem
478  problem(&TanhSolnForPoisson::get_source,node_file_name,
479  element_file_name,
480  poly_file_name);
481 #else
482  //Set up the problem
484  problem(&TanhSolnForPoisson::get_source,node_file_name,
485  element_file_name,
486  poly_file_name);
487 #endif
488 
489 #ifdef ADAPT
490  // Solve the problem with adaptation
491  problem.newton_solve(max_adapt);
492 #else
493  // Solve the problem
494  problem.newton_solve();
495 #endif
496 
497  //Output solution
498  problem.doc_solution(doc_info);
499 
500  //Increment counter for solutions
501  doc_info.number()++;
502  }
503 
504 }
505 
506 
507 
int main(int argc, char *argv[])
Demonstrate how to solve Poisson problem.
~PoissonProblem()
Destructor (empty)
void actions_after_newton_solve()
Update the problem specs before solve (empty)
void doc_solution(DocInfo &doc_info)
Doc the solution.
void get_exact_u(const Vector< double > &x, Vector< double > &u)
Exact solution as a Vector.
RefineableTriangleMesh< ELEMENT > * Bulk_mesh_pt
Pointer to the bulk mesh.
Z2ErrorEstimator * error_estimator_pt
Error estimator.
RefineableTriangleMesh< ELEMENT > * mesh_pt()
Access function for the specific mesh.
double Alpha
Parameter for steepness of step.
void actions_before_newton_solve()
Update the problem specs before solve: (Re)set boundary conditions.
PoissonEquations< 2 >::PoissonSourceFctPt Source_fct_pt
Pointer to source function.
PoissonProblem(PoissonEquations< 2 >::PoissonSourceFctPt source_fct_pt, const string &node_file_name, const string &element_file_name, const string &poly_file_name)
Constructor: Pass pointer to source function and names of two triangle input files.
void get_exact_u(const Vector< double > &x, double &u)
Exact solution as a scalar.
Namespace for exact solution for Poisson equation with sharp step.
TriangleMesh< ELEMENT > * mesh_pt()
Access function for the specific mesh.
double Beta
Parameter for angle of step.
Micky mouse Poisson problem.
Unstructured refineable Triangle Mesh.
void get_source(const Vector< double > &x, double &source)
Source function to make it an exact solution.
void actions_after_adapt()
Actions performed after the adaptation steps.