mpi/solvers/TrilinosSolver_test.cc File Reference
#include "generic.h"
#include "poisson.h"
#include "meshes/rectangular_quadmesh.h"

Classes

class  PoissonProblem< ELEMENT >
 Micky mouse Poisson problem. More...
 

Namespaces

 TanhSolnForPoisson
 Namespace for exact solution for Poisson equation with "sharp step".
 

Functions

void TanhSolnForPoisson::get_exact_u (const Vector< double > &x, Vector< double > &u)
 Exact solution as a Vector. More...
 
void TanhSolnForPoisson::source_function (const Vector< double > &x, double &source)
 Source function required to make the solution above an exact solution. More...
 
int main (int argc, char **argv)
 Driver code for 2D Poisson problem. More...
 

Function Documentation

◆ main()

int main ( int argc  ,
char **  argv 
)

Driver code for 2D Poisson problem.

316 {
317 #ifdef OOMPH_HAS_MPI
318  MPI_Helpers::init(argc,argv);
319 #endif
320 
321 
322  // Note: there is a memory leak in trilinos itself as
323  // demonstrated by this bit of code!
324  if (false)
325  {
326  // set up a parameter list
327  Teuchos::ParameterList ml_parameters;
328  ML_Epetra::SetDefaults("SA", ml_parameters);
329 
330  std::cout << "hello" << std::endl;
331  exit(0);
332  }
333 
334  // Set the orientation of the "step" to 45 degrees
336 
337  // Initial value for the steepness of the "step"
339 
340  // Create the problem with 2D nine-node elements from the
341  // QPoissonElement family. Pass pointer to source function.
344 
345 
346  DocInfo doc_info;
347  doc_info.set_directory("RESLT");
348 
349  // File to report the number of Newton iterations
350  ofstream conv_file;
351  char filename[100];
352  sprintf(filename,"%s/conv.dat",doc_info.directory().c_str());
353  conv_file.open(filename);
354 
355 
356  // Trilinos solvers with the ML preconditioner
357  //--------------------------------------------
358  {
359 
360  cout << "==================================="
361  << "=================================" << endl;
362  cout << "Testing Trilinos solver preconditioned with "
363  << "TrilinosMLPreconditioner" << endl;
364  cout << "======================================="
365  << "=============================" << endl;
366 
367  // Create a Trilinos Solver
368  TrilinosAztecOOSolver* linear_solver_pt = new TrilinosAztecOOSolver;
369 
370  // Create the Trilinos ML preconditioner
371  TrilinosMLPreconditioner* preconditioner_pt = new TrilinosMLPreconditioner;
372 
373  // Set the preconditioner pointer
374  linear_solver_pt->preconditioner_pt() = preconditioner_pt;
375 
376  // Set linear solver
377  problem.linear_solver_pt() = linear_solver_pt;
378 
379  linear_solver_pt->enable_doc_time();
380 
381  // Loop over specific built-in Trilinos solvers
382  for (unsigned s=0; s<=2; s++)
383  {
384  // Select the linear solver
385  switch(s)
386  {
387  case 0:
388  linear_solver_pt->solver_type() = TrilinosAztecOOSolver::CG;
389  cout << "Using CG as solver" << endl;
390  cout << "==================" << endl;
391  break;
392  case 1:
393  linear_solver_pt->solver_type() = TrilinosAztecOOSolver::GMRES;
394  cout << "Using GMRES as solver" << endl;
395  cout << "=====================" << endl;
396  break;
397  case 2:
398  linear_solver_pt->solver_type() = TrilinosAztecOOSolver::BiCGStab;
399  cout << "Using BiCGStab as solver" << endl;
400  cout << "========================" << endl;
401  break;
402  default: cout << "Error selecting solver" << endl;
403  }
404 
405  // Solve the problem
406  problem.newton_solve();
407 
408  // Doc
409  problem.doc_solution(doc_info);
410  doc_info.number()++;
411  oomph_info << "Number of linear solver iterations: "
412  << dynamic_cast<TrilinosAztecOOSolver*>(
413  linear_solver_pt)->iterations() << std::endl << std::endl;
414  conv_file << "Number of linear solver iterations: "
415  << dynamic_cast<TrilinosAztecOOSolver*>(
416  linear_solver_pt)->iterations() << std::endl;
417 
418  } // end of loop over preconditioner
419 
420  // delete the Trilinos solver
421  delete linear_solver_pt;
422  delete preconditioner_pt;
423  }
424 
425 
426 
427  // Trilinos solvers with the IFPACK preconditioner
428  //------------------------------------------------
429  {
430 
431  cout << "================================="
432  << "=======================================" << endl;
433  cout << "Testing Trilinos solver preconditioned"
434  << " with TrilinosIFPACKPreconditioner" << endl;
435  cout << "==========================================="
436  << "=============================" << endl;
437 
438  // Create a Trilinos Solver
439  TrilinosAztecOOSolver* linear_solver_pt = new TrilinosAztecOOSolver;
440 
441  // create the Trilinos IFPACK preconditioner
442  TrilinosIFPACKPreconditioner* preconditioner_pt =
444 
445  // set the preconditioner pointer
446  linear_solver_pt->preconditioner_pt() = preconditioner_pt;
447 
448  // Set linear solver
449  problem.linear_solver_pt() = linear_solver_pt;
450 
451  linear_solver_pt->enable_doc_time();
452 
453  // Loop over specific built-in Trilinos solvers
454  for (unsigned s=0; s<=2; s++)
455  {
456  // Select the linear solver
457  switch(s)
458  {
459  case 0:
460  linear_solver_pt->solver_type() = TrilinosAztecOOSolver::CG;
461  cout << "Using CG as solver" << endl;
462  cout << "==================" << endl;
463  break;
464  case 1:
465  linear_solver_pt->solver_type() = TrilinosAztecOOSolver::GMRES;
466  cout << "Using GMRES as solver" << endl;
467  cout << "=====================" << endl;
468  break;
469  case 2:
470  linear_solver_pt->solver_type() = TrilinosAztecOOSolver::BiCGStab;
471  cout << "Using BiCGStab as solver" << endl;
472  cout << "========================" << endl;
473  break;
474  default: cout << "Error selecting solver" << endl;
475  }
476 
477  // Solve the problem
478  problem.newton_solve();
479 
480  // Doc
481  problem.doc_solution(doc_info);
482  doc_info.number()++;
483  oomph_info << "Number of linear solver iterations: "
484  << dynamic_cast<TrilinosAztecOOSolver*>(
485  linear_solver_pt)->iterations() << std::endl << std::endl;
486  conv_file << "Number of linear solver iterations: "
487  << dynamic_cast<TrilinosAztecOOSolver*>(
488  linear_solver_pt)->iterations() << std::endl;
489 
490  } // end of loop over preconditioner
491 
492  // delete the Trilinos solver
493  delete linear_solver_pt;
494  delete preconditioner_pt;
495  }
496 
497 
498 
499 
500 
501  // Test the oomph-lib iterative linear solvers with Trilinos preconditioners
502  //--------------------------------------------------------------------------
503 
504  cout
505  << "======================================"
506  << "=================================="
507  << endl;
508  cout
509  << "Testing oomph-lib iterative linear solvers "
510  << "with Trilinos preconditioners"
511  << endl;
512  cout
513  << "======================================="
514  << "================================="
515  << endl << endl << endl;
516 
517 
518 
519  // Pointer to linear solver
520  LinearSolver* linear_solver_pt=0;
521 
522  // Loop over two Trilonos preconditioners
523  for (unsigned i_prec=0;i_prec<2;i_prec++)
524  {
525 
526  // Pointer to preconditioner
527  Preconditioner* preconditioner_pt=0;
528 
529  if (i_prec==0)
530  {
531  // Create a IFPACK preconditioner
532  preconditioner_pt = new TrilinosIFPACKPreconditioner;
533  cout << "Using TrilinosIFPACKPreconditioner as preconditioner" << endl;
534  cout << "=====================================================" << endl;
535  }
536  else
537  {
538  // Create a ML preconditioner
539  preconditioner_pt = new TrilinosMLPreconditioner;
540  cout << "Using TrilinosMLPreconditioner as preconditioner" << endl;
541  cout << "=================================================" << endl;
542  }
543 
544 
545  // Test CG
546  //--------
547  cout << "CG iterative solver" << endl;
548  cout << "-------------------" << endl;
549  linear_solver_pt = new CG<CRDoubleMatrix>;
550  static_cast<IterativeLinearSolver*>(linear_solver_pt)->tolerance()=1.0e-10;
551  problem.linear_solver_pt() = linear_solver_pt;
552 
553  // set the preconditioner in the iterative solver
554  static_cast<IterativeLinearSolver*>(linear_solver_pt)->
555  preconditioner_pt()= preconditioner_pt;
556 
557  // Solve the problem
558  problem.newton_solve();
559 
560 
561  // Doc
562  problem.doc_solution(doc_info);
563  doc_info.number()++;
564  oomph_info << "Number of linear solver iterations: "
565  << dynamic_cast<CG<CRDoubleMatrix>*>(
566  linear_solver_pt)->iterations() << std::endl << std::endl;
567  conv_file << "Number of linear solver iterations: "
568  << dynamic_cast<CG<CRDoubleMatrix>*>(
569  linear_solver_pt)->iterations() << std::endl;
570 
571  delete linear_solver_pt;
572 
573 
574  // Test BiCGStab
575  //--------------
576  cout << "BiCGStab iterative solver" << endl;
577  cout << "-------------------------" << endl;
578  linear_solver_pt = new BiCGStab<CRDoubleMatrix>;
579  static_cast<IterativeLinearSolver*>(linear_solver_pt)->tolerance()=1.0e-10;
580  problem.linear_solver_pt() = linear_solver_pt;
581 
582  // set the preconditioner in the iterative solver
583  static_cast<IterativeLinearSolver*>(linear_solver_pt)->preconditioner_pt()
584  = preconditioner_pt;
585 
586  // Solve the problem
587  problem.newton_solve();
588 
589  // Doc
590  problem.doc_solution(doc_info);
591  doc_info.number()++;
592  oomph_info << "Number of linear solver iterations: "
593  << dynamic_cast<BiCGStab<CRDoubleMatrix>*>(
594  linear_solver_pt)->iterations() << std::endl << std::endl;
595  conv_file << "Number of linear solver iterations: "
596  << dynamic_cast<BiCGStab<CRDoubleMatrix>*>(
597  linear_solver_pt)->iterations() << std::endl;
598 
599  delete linear_solver_pt;
600 
601 
602  // Test GMRES
603  //-----------
604  cout << "GMRES iterative solver" << endl;
605  cout << "----------------------" << endl;
606  linear_solver_pt = new GMRES<CRDoubleMatrix>;
607  static_cast<IterativeLinearSolver*>(linear_solver_pt)->tolerance()=1.0e-10;
608  problem.linear_solver_pt() = linear_solver_pt;
609 
610  // set the preconditioner in the iterative solver
611  static_cast<IterativeLinearSolver*>(linear_solver_pt)->preconditioner_pt()
612  = preconditioner_pt;
613 
614  // Solve the problem
615  problem.newton_solve();
616 
617  // Doc
618  problem.doc_solution(doc_info);
619  doc_info.number()++;
620  oomph_info << "Number of linear solver iterations: "
621  << dynamic_cast<GMRES<CRDoubleMatrix>*>(
622  linear_solver_pt)->iterations() << std::endl << std::endl;
623  conv_file << "Number of linear solver iterations: "
624  << dynamic_cast<GMRES<CRDoubleMatrix>*>(
625  linear_solver_pt)->iterations() << std::endl;
626 
627  delete linear_solver_pt;
628 
629  // Kill preconditioner
630  delete preconditioner_pt;
631  preconditioner_pt=0;
632 
633  } // end of loop over preconditioners
634 
635 
636  // Trilinos solvers with the OOMPH-lib
637  //------------------------------------------------
638  {
639 
640  cout << "=================================="
641  << "=======================================" << endl;
642  cout << "Testing Trilinos solver preconditioned with"
643  << " MatrixBasedDiagPreconditioner" << endl;
644  cout << "======================================="
645  << "==================================" << endl;
646 
647  // Create a Trilinos Solver
648  TrilinosAztecOOSolver* linear_solver_pt = new TrilinosAztecOOSolver;
649 
650  // create the oomphlib diagonal preconditioner
651  Preconditioner* preconditioner_pt =
653 
654  // set the preconditioner pointer
655  linear_solver_pt->preconditioner_pt() = preconditioner_pt;
656 
657  // Set linear solver
658  problem.linear_solver_pt() = linear_solver_pt;
659 
660  linear_solver_pt->enable_doc_time();
661 
662  // Loop over specific built-in Trilinos solvers
663  for (unsigned s=0; s<=2; s++)
664  {
665  // Select the linear solver
666  switch(s)
667  {
668  case 0:
669  linear_solver_pt->solver_type() = TrilinosAztecOOSolver::CG;
670  cout << "Using CG as solver" << endl;
671  cout << "==================" << endl;
672  break;
673  case 1:
674  linear_solver_pt->solver_type() = TrilinosAztecOOSolver::GMRES;
675  cout << "Using GMRES as solver" << endl;
676  cout << "=====================" << endl;
677  break;
678  case 2:
679  linear_solver_pt->solver_type() = TrilinosAztecOOSolver::BiCGStab;
680  cout << "Using BiCGStab as solver" << endl;
681  cout << "========================" << endl;
682  break;
683  default: cout << "Error selecting solver" << endl;
684  }
685 
686  // Solve the problem
687  problem.newton_solve();
688 
689  // Doc
690  problem.doc_solution(doc_info);
691  doc_info.number()++;
692  oomph_info << "Number of linear solver iterations: "
693  << dynamic_cast<TrilinosAztecOOSolver*>(
694  linear_solver_pt)->iterations() << std::endl << std::endl;
695  conv_file << "Number of linear solver iterations: "
696  << dynamic_cast<TrilinosAztecOOSolver*>(
697  linear_solver_pt)->iterations() << std::endl;
698 
699  } // end of loop over preconditioner
700 
701  // delete the Trilinos solver
702  delete linear_solver_pt;
703 
704  // Delete oomph-lib preconditioner
705  delete preconditioner_pt;
706  }
707 
708  conv_file.close();
709 
710  // This code is not used, just compiled -- don't delete as it's
711  // extracted by doxygen for a tutorial
712  bool never_get_here=true;
713  if (!never_get_here)
714  {
715  // Create oomph-lib linear solver
716  IterativeLinearSolver* linear_solver_pt=new GMRES<CRDoubleMatrix>;
717 
718  // Create Trilinos IFPACK preconditioner as oomph-lib Preconditioner
719  Preconditioner* preconditioner_pt=new TrilinosIFPACKPreconditioner;
720 
721  // Pass pointer to preconditioner to oomph-lib IterativeLinearSolver
722  linear_solver_pt->preconditioner_pt()=preconditioner_pt;
723  }
724 
725 #ifdef OOMPH_HAS_MPI
726 MPI_Helpers::finalize();
727 #endif
728 
729 
730 } //end of main
Micky mouse Poisson problem.
Definition: HypreSolver_test.cc:81
The conjugate gradient method.
Definition: iterative_linear_solver.h:410
The conjugate gradient method.
Definition: iterative_linear_solver.h:284
Definition: oomph_utilities.h:499
std::string directory() const
Output directory.
Definition: oomph_utilities.h:524
void set_directory(const std::string &directory)
Definition: oomph_utilities.cc:298
unsigned & number()
Number used (e.g.) for labeling output files.
Definition: oomph_utilities.h:554
The GMRES method.
Definition: iterative_linear_solver.h:1227
Definition: iterative_linear_solver.h:54
Preconditioner *& preconditioner_pt()
Access function to preconditioner.
Definition: iterative_linear_solver.h:95
Definition: linear_solver.h:68
void enable_doc_time()
Enable documentation of solve times.
Definition: linear_solver.h:110
Matrix-based diagonal preconditioner.
Definition: general_purpose_preconditioners.h:49
Definition: preconditioner.h:54
Definition: trilinos_solver.h:267
unsigned & solver_type()
Access function to Solver_type.
Definition: trilinos_solver.h:442
Definition: trilinos_preconditioners.h:268
Definition: trilinos_preconditioners.h:152
RealScalar s
Definition: level1_cplx_impl.h:130
string filename
Definition: MergeRestartFiles.py:39
double TanPhi
Parameter for angle Phi of "step".
Definition: HypreSolver_test.cc:51
void source_function(const Vector< double > &x, double &source)
Source function required to make the solution above an exact solution.
Definition: HypreSolver_test.cc:60
double Alpha
Parameter for steepness of step.
Definition: extrude_with_macro_element_representation.cc:185
OomphInfo oomph_info
Definition: oomph_definitions.cc:319
Constructor for SteadyAxisymAdvectionDiffusion problem
Definition: steady_axisym_advection_diffusion.cc:213

References TanhSolnForPoisson::Alpha, oomph::DocInfo::directory(), oomph::LinearSolver::enable_doc_time(), MergeRestartFiles::filename, oomph::DocInfo::number(), oomph::oomph_info, oomph::IterativeLinearSolver::preconditioner_pt(), problem, s, oomph::DocInfo::set_directory(), oomph::TrilinosAztecOOSolver::solver_type(), TanhSolnForPoisson::source_function(), and TanhSolnForPoisson::TanPhi.