two_d_poisson_compare_solvers.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 -- compare various linear solvers
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 //===== start_of_namespace=============================================
46 /// Namespace for exact solution for Poisson equation with "sharp step"
47 //=====================================================================
48 namespace TanhSolnForPoisson
49 {
50 
51  /// Parameter for steepness of "step"
52  double Alpha=1.0;
53 
54  /// Parameter for angle Phi of "step"
55  double TanPhi=0.0;
56 
57  /// Exact solution as a Vector
58  void get_exact_u(const Vector<double>& x, Vector<double>& u)
59  {
60  u[0]=tanh(1.0-Alpha*(TanPhi*x[0]-x[1]));
61  }
62 
63  /// Source function required to make the solution above an exact solution
64  void get_source(const Vector<double>& x, double& source)
65  {
66  source = 2.0*tanh(-1.0+Alpha*(TanPhi*x[0]-x[1]))*
67  (1.0-pow(tanh(-1.0+Alpha*(TanPhi*x[0]-x[1])),2.0))*
68  Alpha*Alpha*TanPhi*TanPhi+2.0*tanh(-1.0+Alpha*(TanPhi*x[0]-x[1]))*
69  (1.0-pow(tanh(-1.0+Alpha*(TanPhi*x[0]-x[1])),2.0))*Alpha*Alpha;
70  }
71 
72 } // end of namespace
73 
74 
75 
76 
77 
78 
79 
80 
81 //====== start_of_problem_class=======================================
82 /// 2D Poisson problem on rectangular domain, discretised with
83 /// 2D QPoisson elements. The specific type of element is
84 /// specified via the template parameter.
85 //====================================================================
86 template<class ELEMENT>
87 class PoissonProblem : public Problem
88 {
89 
90 public:
91 
92  /// \short Constructor: Pass pointer to source function and flag to indicate
93  /// if order of elements should be messed up. nel_1d is the number of
94  /// elements along the 1d mesh edges -- total number of elements is
95  /// nel_1d x nel_1d.
96  PoissonProblem(PoissonEquations<2>::PoissonSourceFctPt source_fct_pt,
97  const unsigned& nel_1d,
98  const bool& mess_up_order);
99 
100  /// Destructor (empty)
102 
103  /// \short Update the problem specs before solve: Reset boundary conditions
104  /// to the values from the exact solution.
105  void actions_before_newton_solve();
106 
107  /// Update the problem after solve (empty)
109 
110  /// \short Doc the solution. DocInfo object stores flags/labels for where the
111  /// output gets written to
112  void doc_solution(DocInfo& doc_info);
113 
114 private:
115 
116  /// Pointer to source function
117  PoissonEquations<2>::PoissonSourceFctPt Source_fct_pt;
118 
119 }; // end of problem class
120 
121 
122 
123 
124 //=====start_of_constructor===============================================
125 /// Constructor for Poisson problem: Pass pointer to source function
126 /// and a flag that specifies if the order of the elements should
127 /// be messed up. nel_1d is the number of elements along
128 /// the 1d mesh edges -- total number of elements is nel_1d x nel_1d.
129 //========================================================================
130 template<class ELEMENT>
132 PoissonProblem(PoissonEquations<2>::PoissonSourceFctPt source_fct_pt,
133  const unsigned& nel_1d, const bool& mess_up_order)
134  : Source_fct_pt(source_fct_pt)
135 {
136 
137  // Setup mesh
138 
139  // # of elements in x-direction
140  unsigned n_x=nel_1d;
141 
142  // # of elements in y-direction
143  unsigned n_y=nel_1d;
144 
145  // Domain length in x-direction
146  double l_x=1.0;
147 
148  // Domain length in y-direction
149  double l_y=2.0;
150 
151  // Build and assign mesh
152  Problem::mesh_pt() = new SimpleRectangularQuadMesh<ELEMENT>(n_x,n_y,l_x,l_y);
153 
154 
155 
156  // Mess up order of elements in the Problem's mesh:
157  //-------------------------------------------------
158  if (mess_up_order)
159  {
160  unsigned n_element=mesh_pt()->nelement();
161 
162  // Make copy
163  Vector<ELEMENT*> tmp_element_pt(n_element);
164  for (unsigned e=0;e<n_element;e++)
165  {
166  tmp_element_pt[e]=dynamic_cast<ELEMENT*>(mesh_pt()->element_pt(e));
167  }
168 
169  // Reorder
170  unsigned e_half=unsigned(0.5*double(n_element));
171  unsigned e_lo=e_half-3;
172  unsigned e_up=n_element-e_lo;
173  unsigned count=0;
174  for (unsigned e=0;e<e_lo;e++)
175  {
176  mesh_pt()->element_pt(count)=tmp_element_pt[e];
177  count++;
178  mesh_pt()->element_pt(count)=tmp_element_pt[n_element-e-1];
179  count++;
180  }
181  for (unsigned e=e_lo;e<e_up;e++)
182  {
183  mesh_pt()->element_pt(count)=tmp_element_pt[e];
184  count++;
185  }
186  }
187 
188  // Set the boundary conditions for this problem: All nodes are
189  // free by default -- only need to pin the ones that have Dirichlet conditions
190  // here.
191  unsigned num_bound = mesh_pt()->nboundary();
192  for(unsigned ibound=0;ibound<num_bound;ibound++)
193  {
194  unsigned num_nod= mesh_pt()->nboundary_node(ibound);
195  for (unsigned inod=0;inod<num_nod;inod++)
196  {
197  mesh_pt()->boundary_node_pt(ibound,inod)->pin(0);
198  }
199  }
200 
201  // Complete the build of all elements so they are fully functional
202 
203  // Loop over the elements to set up element-specific
204  // things that cannot be handled by the (argument-free!) ELEMENT
205  // constructor: Pass pointer to source function
206  unsigned n_element = mesh_pt()->nelement();
207  for(unsigned i=0;i<n_element;i++)
208  {
209  // Upcast from GeneralsedElement to the present element
210  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(mesh_pt()->element_pt(i));
211 
212  //Set the source function pointer
213  el_pt->source_fct_pt() = Source_fct_pt;
214  }
215 
216 
217  // Setup equation numbering scheme
218  cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
219 
220 } // end of constructor
221 
222 
223 
224 
225 //========================================start_of_actions_before_newton_solve===
226 /// Update the problem specs before solve: (Re-)set boundary conditions
227 /// to the values from the exact solution.
228 //========================================================================
229 template<class ELEMENT>
231 {
232  // How many boundaries are there?
233  unsigned num_bound = mesh_pt()->nboundary();
234 
235  //Loop over the boundaries
236  for(unsigned ibound=0;ibound<num_bound;ibound++)
237  {
238  // How many nodes are there on this boundary?
239  unsigned num_nod=mesh_pt()->nboundary_node(ibound);
240 
241  // Loop over the nodes on boundary
242  for (unsigned inod=0;inod<num_nod;inod++)
243  {
244  // Get pointer to node
245  Node* nod_pt=mesh_pt()->boundary_node_pt(ibound,inod);
246 
247  // Extract nodal coordinates from node:
248  Vector<double> x(2);
249  x[0]=nod_pt->x(0);
250  x[1]=nod_pt->x(1);
251 
252  // Compute the value of the exact solution at the nodal point
253  Vector<double> u(1);
255 
256  // Assign the value to the one (and only) nodal value at this node
257  nod_pt->set_value(0,u[0]);
258  }
259  }
260 } // end of actions before solve
261 
262 
263 
264 //===============start_of_doc=============================================
265 /// Doc the solution: doc_info contains labels/output directory etc.
266 //========================================================================
267 template<class ELEMENT>
268 void PoissonProblem<ELEMENT>::doc_solution(DocInfo& doc_info)
269 {
270 
271  ofstream some_file;
272  char filename[100];
273 
274  // Number of plot points: npts x npts
275  unsigned npts=5;
276 
277  // Output solution
278  //-----------------
279  sprintf(filename,"%s/soln%i.dat",doc_info.directory().c_str(),
280  doc_info.number());
281  some_file.open(filename);
282  mesh_pt()->output(some_file,npts);
283  some_file.close();
284 
285 
286  // Output exact solution
287  //----------------------
288  sprintf(filename,"%s/exact_soln%i.dat",doc_info.directory().c_str(),
289  doc_info.number());
290  some_file.open(filename);
291  mesh_pt()->output_fct(some_file,npts,TanhSolnForPoisson::get_exact_u);
292  some_file.close();
293 
294  // Doc error and return of the square of the L2 error
295  //---------------------------------------------------
296  double error,norm;
297  sprintf(filename,"%s/error%i.dat",doc_info.directory().c_str(),
298  doc_info.number());
299  some_file.open(filename);
300  mesh_pt()->compute_error(some_file,TanhSolnForPoisson::get_exact_u,
301  error,norm);
302  some_file.close();
303 
304  // Doc L2 error and norm of solution
305  cout << "\nNorm of error : " << sqrt(error) << std::endl;
306  cout << "Norm of solution: " << sqrt(norm) << std::endl << std::endl;
307 
308 } // end of doc
309 
310 
311 
312 
313 
314 
315 //===== start_of_run======================================================
316 /// Build and run problem with specified linear solver. Also pass
317 /// flag to specify if order of elements should be messed up.
318 /// nel_1d is the number of elements along the 1D mesh edge.
319 /// Total number of elements is nel_1d x nel_1d.
320 //========================================================================
321 void run(const string& dir_name,
322  LinearSolver* linear_solver_pt,
323  const unsigned nel_1d,
324  bool mess_up_order)
325 {
326 
327  //Set up the problem
328  //------------------
329 
330  // Create the problem with 2D nine-node elements from the
331  // QPoissonElement family. Pass pointer to source function
332  // and flag to specify if order of elements should be messed up.
334  problem(&TanhSolnForPoisson::get_source,nel_1d,mess_up_order);
335 
336 
337  /// Set linear solver:
338  //--------------------
339  problem.linear_solver_pt() = linear_solver_pt;
340 
341 
342  // Create label for output
343  //------------------------
344  DocInfo doc_info;
345 
346  // Set output directory
347  doc_info.set_directory(dir_name);
348 
349  // Step number
350  doc_info.number()=0;
351 
352  // Check if we're ready to go:
353  //----------------------------
354  cout << "\n\n\nProblem self-test:";
355  if (problem.self_test()==0)
356  {
357  cout << "passed: Problem can be solved." << std::endl;
358  }
359  else
360  {
361  throw OomphLibError("Self test failed",
362  OOMPH_CURRENT_FUNCTION,
363  OOMPH_EXCEPTION_LOCATION);
364  }
365 
366 
367  // Set the orientation of the "step" to 45 degrees
369 
370  // Initial value for the steepness of the "step"
372 
373  // Solve the problem
374  problem.newton_solve();
375 
376  //Output the solution
377  problem.doc_solution(doc_info);
378 
379 
380 } //end of run
381 
382 
383 
384 
385 
386 //===== start_of_main=====================================================
387 /// Driver code for 2D Poisson problem -- compare different solvers
388 //========================================================================
389 int main()
390 {
391 
392  // Pointer to linear solver
393  LinearSolver* linear_solver_pt;
394 
395  // Result directory
396  string dir_name;
397 
398  // Cpu start/end times
399  clock_t cpu_start,cpu_end;
400 
401  // Storage for cpu times with and without messed up order
402  Vector<double> superlu_cr_cpu(2);
403  Vector<double> superlu_cc_cpu(2);
404  Vector<double> hsl_ma42_cpu(2);
405  Vector<double> hsl_ma42_reordered_cpu(2);
406  Vector<double> dense_lu_cpu(2);
407  Vector<double> fd_lu_cpu(2);
408 
409  // Flag to indicate if order should be messed up:
410  bool mess_up_order;
411 
412  // Number of elements along 1D edge of mesh; total number of elements
413  // is nel_1d x nel_1d
414  unsigned nel_1d;
415 
416  // Do run with and without messed up elements
417  for (unsigned mess=0;mess<2;mess++)
418  {
419 
420  // Mess up order?
421  if (mess==0)
422  {
423  mess_up_order=true;
424  }
425  else
426  {
427  mess_up_order=false;
428  }
429 
430  /// Use SuperLU with compressed row storage
431  //-----------------------------------------
432 
433  cout << std::endl;
434  cout << " Use SuperLU with compressed row storage: " << std::endl;
435  cout << "========================================= " << std::endl;
436  cout << std::endl;
437 
438  // Start cpu clock
439  cpu_start=clock();
440 
441  // Build linear solver
442  linear_solver_pt = new SuperLUSolver;
443 
444  /// Use compressed row storage
445  static_cast<SuperLUSolver*>(linear_solver_pt)
446  ->use_compressed_row_for_superlu_serial();
447 
448  /// Switch on full doc
449  static_cast<SuperLUSolver*>(linear_solver_pt)->enable_doc_stats();
450 
451  // Choose result directory
452  dir_name="RESLT_cr";
453 
454  // Number of elements along 1D edge of mesh; total number of elements
455  // is nel_1d x nel_1d
456  nel_1d=20;
457 
458  // Run it
459  run(dir_name,linear_solver_pt,nel_1d,mess_up_order);
460 
461  // Note: solver does not have to be deleted -- it's killed automatically
462  // in the problem destructor.
463 
464  // End cpu clock
465  cpu_end=clock();
466 
467  // Timing
468  superlu_cr_cpu[mess]=double(cpu_end-cpu_start)/CLOCKS_PER_SEC;
469 
470 
471 
472 
473  /// Use SuperLU with compressed column storage
474  //--------------------------------------------
475 
476 
477  cout << std::endl;
478  cout << " Use SuperLU with compressed column storage: " << std::endl;
479  cout << "============================================ " << std::endl;
480  cout << std::endl;
481 
482 
483  // Start cpu clock
484  cpu_start=clock();
485 
486  // Build linear solver
487  linear_solver_pt = new SuperLUSolver;
488 
489  /// Use compressed row storage
490  static_cast<SuperLUSolver*>(linear_solver_pt)
491  ->use_compressed_column_for_superlu_serial();
492 
493  /// Switch on full doc
494  static_cast<SuperLUSolver*>(linear_solver_pt)->enable_doc_stats();
495 
496  // Choose result directory
497  dir_name="RESLT_cc";
498 
499  // Number of elements along 1D edge of mesh; total number of elements
500  // is nel_1d x nel_1d
501  nel_1d=20;
502 
503  // Run it
504  run(dir_name,linear_solver_pt,nel_1d,mess_up_order);
505 
506  // Note: solver does not have to be deleted -- it's killed automatically
507  // in the problem destructor.
508 
509  // End cpu clock
510  cpu_end=clock();
511 
512  // Timing
513  superlu_cc_cpu[mess]=double(cpu_end-cpu_start)/CLOCKS_PER_SEC;
514 
515 
516 #ifdef HAVE_HSL_SOURCES
517 
518  /// Use HSL frontal solver MA42 without element re-ordering
519  //---------------------------------------------------------
520 
521  cout << std::endl;
522  cout << " Use HSL frontal solver MA42 without element re-ordering: " << std::endl;
523  cout << "========================================================= " << std::endl;
524  cout << std::endl;
525 
526  // Start cpu clock
527  cpu_start=clock();
528 
529  // Build linear solver
530  linear_solver_pt = new HSL_MA42;
531 
532  /// Switch on full doc
533  static_cast<HSL_MA42*>(linear_solver_pt)->enable_doc_stats();
534 
535  // Choose result directory
536  dir_name="RESLT_frontal";
537 
538 
539  // Number of elements along 1D edge of mesh; total number of elements
540  // is nel_1d x nel_1d
541  nel_1d=20;
542 
543  // Run it
544  run(dir_name,linear_solver_pt,nel_1d,mess_up_order);
545 
546  // Note: solver does not have to be deleted -- it's killed automatically
547  // in the problem destructor.
548 
549 
550  // End cpu clock
551  cpu_end=clock();
552 
553  // Timing
554  hsl_ma42_cpu[mess]=double(cpu_end-cpu_start)/CLOCKS_PER_SEC;
555 
556 
557 
558  /// Use HSL frontal solver MA42 with element reordering by Sloan's algorithm
559  //--------------------------------------------------------------------------
560 
561  cout << std::endl;
562  cout << " Use HSL frontal solver MA42 with element re-ordering: " << std::endl;
563  cout << "====================================================== " << std::endl;
564  cout << std::endl;
565 
566  // Start cpu clock
567  cpu_start=clock();
568 
569  // Build linear solver
570  linear_solver_pt = new HSL_MA42;
571 
572  /// Switch on full doc
573  static_cast<HSL_MA42*>(linear_solver_pt)->enable_doc_stats();
574 
575 
576  /// Switch on re-ordering
577  static_cast<HSL_MA42*>(linear_solver_pt)->enable_reordering();
578 
579  // Choose result directory
580  dir_name="RESLT_frontal_reordered";
581 
582 
583  // Number of elements along 1D edge of mesh; total number of elements
584  // is nel_1d x nel_1d
585  nel_1d=20;
586 
587  // Run it
588  run(dir_name,linear_solver_pt,nel_1d,mess_up_order);
589 
590  // Note: solver does not have to be deleted -- it's killed automatically
591  // in the problem destructor.
592 
593 
594  // End cpu clock
595  cpu_end=clock();
596 
597  // Timing
598  hsl_ma42_reordered_cpu[mess]=double(cpu_end-cpu_start)/CLOCKS_PER_SEC;
599 
600 
601 #endif
602 
603  /// Use dense matrix LU decomposition
604  //-----------------------------------
605 
606  cout << std::endl;
607  cout << " Use dense matrix LU decomposition: " << std::endl;
608  cout << "=================================== " << std::endl;
609  cout << std::endl;
610 
611  // Start cpu clock
612  cpu_start=clock();
613 
614  // Build linear solver
615  linear_solver_pt = new DenseLU;
616 
617  // Choose result directory
618  dir_name="RESLT_dense_LU";
619 
620  // Number of elements along 1D edge of mesh; total number of elements
621  // is nel_1d x nel_1d
622  nel_1d=4;
623 
624  // Run it
625  run(dir_name,linear_solver_pt,nel_1d,mess_up_order);
626 
627  // Note: solver does not have to be deleted -- it's killed automatically
628  // in the problem destructor.
629 
630 
631  // End cpu clock
632  cpu_end=clock();
633 
634  // Timing
635  dense_lu_cpu[mess]=double(cpu_end-cpu_start)/CLOCKS_PER_SEC;
636 
637 
638 
639 
640  /// Use dense matrix LU decomposition
641  //-----------------------------------
642 
643  cout << std::endl;
644  cout << " Use dense FD-ed matrix LU decomposition: " << std::endl;
645  cout << "========================================= " << std::endl;
646  cout << std::endl;
647 
648  // Start cpu clock
649  cpu_start=clock();
650 
651  // Build linear solver
652  linear_solver_pt = new FD_LU;
653 
654  // Choose result directory
655  dir_name="RESLT_FD_LU";
656 
657  // Number of elements along 1D edge of mesh; total number of elements
658  // is nel_1d x nel_1d
659  nel_1d=4;
660 
661  // Run it
662  run(dir_name,linear_solver_pt,nel_1d,mess_up_order);
663 
664  // Note: solver does not have to be deleted -- it's killed automatically
665  // in the problem destructor.
666 
667  // End cpu clock
668  cpu_end=clock();
669 
670  // Timing
671  fd_lu_cpu[mess]=double(cpu_end-cpu_start)/CLOCKS_PER_SEC;
672 
673  }
674 
675 
676  // Doc timings with and without messed up elements
677  for (unsigned mess=0;mess<2;mess++)
678  {
679  cout << std::endl << std::endl << std::endl ;
680  if (mess==0)
681  {
682  cout << "TIMINGS WITH MESSED UP ORDERING OF ELEMENTS: " << std::endl;
683  cout << "============================================ " << std::endl;
684 
685  }
686  else
687  {
688  cout << "TIMINGS WITHOUT MESSED UP ORDERING OF ELEMENTS: " << std::endl;
689  cout << "=============================================== " << std::endl;
690  }
691 
692  cout << "CPU time with SuperLU compressed row : "
693  << superlu_cr_cpu[mess] << std::endl;
694  cout << "CPU time with SuperLU compressed col : "
695  << superlu_cc_cpu[mess] << std::endl;
696 #ifdef HAVE_HSL_SOURCES
697  cout << "CPU time with MA42 frontal solver : "
698  << hsl_ma42_cpu[mess] << std::endl;
699  cout << "CPU time with MA42 frontal solver (incl. reordering) : "
700  << hsl_ma42_reordered_cpu[mess] << std::endl;
701 #endif
702  cout << "CPU time with dense LU solver (fewer elements!) : "
703  << dense_lu_cpu[mess] << std::endl;
704  cout << "CPU time with dense LU solver & FD (fewer elements!) : "
705  << fd_lu_cpu[mess] << std::endl;
706  cout << std::endl << std::endl << std::endl ;
707  }
708 
709 
710 
711 
712 } //end of main
713 
714 
715 
716 
717 
int main()
Driver code for 2D Poisson problem – compare different solvers.
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.
void run(const string &dir_name, LinearSolver *linear_solver_pt, const unsigned nel_1d, bool mess_up_order)
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.