hypre_solver.h
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-2023 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 #ifndef OOMPH_HYPRE_SOLVER_HEADER
27 #define OOMPH_HYPRE_SOLVER_HEADER
28 
29 // Headers required for hypre
30 #include "_hypre_utilities.h"
31 #include "HYPRE.h"
32 #include "HYPRE_parcsr_ls.h"
33 #include "HYPRE_krylov.h"
34 #include "HYPRE_IJ_mv.h"
35 #include "HYPRE_parcsr_mv.h"
36 
37 // OOMPH-LIB includes
39 #include "linear_solver.h"
40 #include "preconditioner.h"
41 
42 // hypre's global error flag
43 extern int hypre__global_error;
44 
45 
46 namespace oomph
47 {
48  //==================================================================
49  /// Helper functions for use with the Hypre library
50  //==================================================================
51  namespace HypreHelpers
52  {
53  /// Default for AMG strength (0.25 recommended for 2D problems;
54  /// larger (0.5-0.75, say) for 3D
55  extern double AMG_strength;
56 
57  /// Default AMG coarsening strategy. Coarsening types include:
58  /// 0 = CLJP (parallel coarsening using independent sets)
59  /// 1 = classical RS with no boundary treatment (not recommended
60  /// in parallel)
61  /// 3 = modified RS with 3rd pass to add C points on the boundaries
62  /// 6 = Falgout (uses 1 then CLJP using interior coarse points as
63  /// first independent set)
64  /// 8 = PMIS (parallel coarsening using independent sets - lower
65  /// complexities than 0, maybe also slower convergence)
66  /// 10= HMIS (one pass RS on each processor then PMIS on interior
67  /// coarse points as first independent set)
68  /// 11= One pass RS on each processor (not recommended)
69  extern unsigned AMG_coarsening;
70 
71  /// AMG interpolation truncation factor
72  extern double AMG_truncation; // 0.0;
73 
74  /// Helper function to check the Hypre error flag, return the
75  /// message associated with any error, and reset the error flag to zero.
76  int check_HYPRE_error_flag(std::ostringstream& message);
77 
78  /// Helper function to create a HYPRE_IJVector and HYPRE_ParVector.
79  /// + If no MPI then serial vectors are created
80  /// + If MPI and serial input vector then distributed hypre vectors are
81  /// created
82  /// + If MPI and distributed input vector the distributed output vectors
83  /// are created.
84  extern void create_HYPRE_Vector(const DoubleVector& oomph_vec,
85  const LinearAlgebraDistribution* dist_pt,
86  HYPRE_IJVector& hypre_ij_vector,
87  HYPRE_ParVector& hypre_par_vector);
88 
89  /// Helper function to create an empty HYPRE_IJVector and
90  /// HYPRE_ParVector.
91  /// + If no MPI then serial vectors are created
92  /// + If MPI and serial distribution then distributed hypre vectors are
93  /// created
94  /// + If MPI and distributed input distribution the distributed output
95  /// vectors are created.
96  void create_HYPRE_Vector(const LinearAlgebraDistribution* oomph_vec,
97  HYPRE_IJVector& hypre_ij_vector,
98  HYPRE_ParVector& hypre_par_vector);
99 
100  /// Helper function to create a serial HYPRE_IJMatrix and
101  /// HYPRE_ParCSRMatrix from a CRDoubleMatrix
102  void create_HYPRE_Matrix(CRDoubleMatrix* oomph_matrix,
103  HYPRE_IJMatrix& hypre_ij_matrix,
104  HYPRE_ParCSRMatrix& hypre_par_matrix,
105  LinearAlgebraDistribution* dist_pt);
106 
107  /// Helper function to set Euclid options using a command line
108  /// like array.
109  void euclid_settings_helper(const bool& use_block_jacobi,
110  const bool& use_row_scaling,
111  const bool& use_ilut,
112  const int& level,
113  const double& drop_tol,
114  const int& print_level,
115  HYPRE_Solver& euclid_object);
116 
117  } // namespace HypreHelpers
118 
119 
120  //=====================================================================
121  /// An interface class to the suite of Hypre solvers and preconditioners
122  /// to allow use of:
123  ///
124  /// BoomerAMG (AMG), CG, GMRES or BiCGStab, Euclid (ILU) or
125  /// ParaSails (Approximate inverse)
126  ///
127  /// Hypre's Krylov subspace solvers (CG, GMRES and BiCGStab) may be
128  /// preconditioned using:
129  ///
130  /// BoomerAMG, Euclid or ParaSails
131  ///
132  //====================================================================
134  {
135  public:
136  /// Constructor
138  {
139 #ifdef PARANOID
140 #ifndef HYPRE_SEQUENTIAL
141  // For the MPI version of Hypre, check MPI_Helpers::setup has been called
143  {
144  std::ostringstream error_message;
145  error_message << "When using the MPI version of Hypre please first\n"
146  << "call function MPI_Helpers::setup()\n";
147  throw OomphLibError(error_message.str(),
148  OOMPH_CURRENT_FUNCTION,
149  OOMPH_EXCEPTION_LOCATION);
150  }
151 #endif
152 #endif
153 
154  // setup the distribution
156 
157  // These keep track of which solver and preconditioner
158  // (if any) currently exist
161 
162  // Do we want to output info and results of timings?
163  Output_info = true;
164 
165  // General control paramaters
166  Tolerance = 1e-10;
167  Max_iter = 100;
170 
171  // Default AMG control parameters -- these seem OK;
172  // see hypre documenation for details.
174  if (MPI_Helpers::communicator_pt()->nproc() > 1)
175  {
176  // Jacobi in parallel
178  }
179  else
180  {
181  // Gauss Seidel in serial
183  }
186  AMG_damping = 1.0;
190  AMG_max_levels = 100;
191  AMG_max_row_sum = 1.0;
192  AMG_print_level = 0;
193 
194  // Parameters for using Euclicd as an AMG smoother (defaults copied
195  // from the normal defaults listed in the manual).
200  AMGEuclidSmoother_drop_tol = 0; // No dropping
202 
203  // Print level for CG, GMRES and BiCGStab
204  Krylov_print_level = 0;
205 
206  // Set ParaSails default values
207  ParaSails_symmetry = 2;
208  ParaSails_nlevel = 1;
209  ParaSails_thresh = 0.1;
210  ParaSails_filter = 0.1;
211 
212  // Set Euclid default values
213  Euclid_droptol = 0.0;
214  Euclid_rowScale = false;
215  Euclid_using_ILUT = false;
216  Euclid_using_BJ = false;
217  Euclid_level = 1;
218  Euclid_print_level = 0;
219 
220  // Set to true to periodically check the hypre error flag
221  // and output any messages
222  Hypre_error_messages = false;
223  }
224 
225  /// Destructor.
227  {
228  // call function to delete solver data
230 
231  // delete teh oomph-lib distribution
232  delete Hypre_distribution_pt;
233  }
234 
235  /// Broken copy constructor.
236  HypreInterface(const HypreInterface&) = delete;
237 
238  /// Broken assignment operator.
239  void operator=(const HypreInterface&) = delete;
240 
241  /// Turn on the hypre error messages
243  {
244  Hypre_error_messages = true;
245  }
246 
247  /// Turn off hypre error messages
249  {
250  Hypre_error_messages = false;
251  }
252 
253  /// Enumerated flag to define which Hypre methods are used
254  /// CAREFUL: DON'T CHANGE THE ORDER OF THESE!
256  {
257  CG,
263  None
264  };
265 
266  /// Function to return value of which solver (if any) is currently stored.
267  unsigned existing_solver()
268  {
269  return Existing_solver;
270  }
271 
272  /// Function return value of which preconditioner (if any) is stored.
274  {
276  }
277 
278  //??ds comment, write access functions and move to protected?
285 
286  protected:
287  /// Function deletes all solver data.
288  void hypre_clean_up_memory();
289 
290  /// Function which sets values of First_global_row,
291  /// Last_global_row and other partitioning data and creates the distributed
292  /// Hypre matrix (stored in Matrix_ij/Matrix_par) from the CRDoubleMatrix.
293  void hypre_matrix_setup(CRDoubleMatrix* matrix_pt);
294 
295  /// Sets up the data required for to use as an oomph-lib
296  /// LinearSolver or Preconditioner. This must be called after
297  /// the Hypre matrix has been generated using hypre_matrix_setup(...).
298  void hypre_solver_setup();
299 
300  /// Helper function performs a solve if any solver
301  /// exists.
302  void hypre_solve(const DoubleVector& rhs, DoubleVector& solution);
303 
304  /// Flag is true to output info and results of timings
306 
307  /// Maximum number of iterations used in solver.
308  unsigned Max_iter;
309 
310  /// Tolerance used to terminate solver.
311  double Tolerance;
312 
313  /// Hypre method flag. Valid values are specified in enumeration
314  unsigned Hypre_method;
315 
316  /// Preconditioner method flag used with Hypre's PCG,
317  /// GMRES and BiCGStab in solve(...) or resolve(...). Valid values
318  /// are BoomerAMG, Euclid, ParaSails or None (all enumerated above),
319  /// for any other value no preconditioner is set.
321 
322  /// Used to set the Hypre printing level for AMG
323  /// 0: no printout
324  /// 1: print setup information
325  /// 2: print solve information
326  /// 3: print setup and solve information
327  unsigned AMG_print_level;
328 
329  /// Maximum number of levels used in AMG
330  unsigned AMG_max_levels;
331 
332  /// Parameter to identify diagonally dominant parts of the matrix in AMG
334 
335  /// Flag to determine whether simple smoothers (determined by the
336  /// AMG_simple_smoother flag) or complex smoothers (determined by the
337  /// AMG_complex_smoother flag are used in AMG
339 
340  /// Simple smoothing methods used in BoomerAMG. Relaxation types
341  /// include:
342  /// 0 = Jacobi
343  /// 1 = Gauss-Seidel, sequential
344  /// (very slow in parallel!)
345  /// 2 = Gauss-Seidel, interior points in parallel, boundary sequential
346  /// (slow in parallel!)
347  /// 3 = hybrid Gauss-Seidel or SOR, forward solve
348  /// 4 = hybrid Gauss-Seidel or SOR, backward solve
349  /// 6 = hybrid symmetric Gauss-Seidel or SSOR
350  /// To use these methods set AMG_using_simple_smoothing to true
352 
353  /// Complex smoothing methods used in BoomerAMG. Relaxation types
354  /// are:
355  /// 6 = Schwarz
356  /// 7 = Pilut
357  /// 8 = ParaSails
358  /// 9 = Euclid
359  /// To use these methods set AMG_using_simple_smoothing to false
361 
362  /// The number of smoother iterations to apply
364 
365  /// Damping factor for BoomerAMG smoothed Jacobi or hybrid SOR
366  double AMG_damping;
367 
368  /// Connection strength threshold parameter for BoomerAMG
369  double AMG_strength;
370 
371  /// Interpolation truncation factor for BoomerAMG
373 
374  /// AMG coarsening strategy. Coarsening types include:
375  /// 0 = CLJP (parallel coarsening using independent sets)
376  /// 1 = classical RS with no boundary treatment (not recommended
377  /// in parallel)
378  /// 3 = modified RS with 3rd pass to add C points on the boundaries
379  /// 6 = Falgout (uses 1 then CLJP using interior coarse points as
380  /// first independent set)
381  /// 8 = PMIS (parallel coarsening using independent sets - lower
382  /// complexities than 0, maybe also slower convergence)
383  /// 10= HMIS (one pass RS on each processor then PMIS on interior
384  /// coarse points as first independent set)
385  /// 11= One pass RS on each processor (not recommended)
386  unsigned AMG_coarsening;
387 
388  /// ParaSails symmetry flag, used to inform ParaSails of
389  /// Symmetry of definitenss of problem and type of ParaSails
390  /// preconditioner:
391  /// 0 = nonsymmetric and/or indefinite problem, nonsymmetric preconditioner
392  /// 1 = SPD problem, and SPD (factored preconditioner)
393  /// 2 = nonsymmetric, definite problem and SDP (factored preconditoner)
395 
396  /// ParaSails nlevel parameter
398 
399  /// ParaSails thresh parameter
401 
402  /// ParaSails filter parameter
404 
405  /// Euclid drop tolerance for ILU(k) and ILUT factorization
407 
408  /// Flag to switch on Euclid row scaling
410 
411  /// Flag to determine if ILUT (if true) or ILU(k) is used in Euclid
413 
414  /// Flag to determine if Block Jacobi is used instead of PILU
416 
417  /// Euclid level parameter for ILU(k) factorization
419 
420  /// Flag to set the level of printing from Euclid
421  /// when the Euclid destructor is called
422  /// 0: no printing (default)
423  /// 1: prints summary of runtime settings and timings
424  /// 2: as 1 plus prints memory usage
426 
427  /// Used to set the Hypre printing level for the Krylov
428  /// subspace solvers
430 
431  /// Flag to determine if non-zero values of the Hypre error flag
432  /// plus Hypre error messages are output to screen at various points
433  /// in the solve function, i.e. after:
434  /// 1. setting up the Hypre matrix
435  /// 2. setting up the preconditioner
436  /// 3. setting up the solver
437  /// 4. setting up the Hypre vectors
438  /// 5. solving
439  /// 6. getting the solution vector
440  /// 7. deallocation of solver data
442 
443  /// Internal flag which is true when hypre_setup or hypre_solve
444  /// can delete input matrix.
446 
447 #ifdef OOMPH_HAS_MPI
448  /// Internal flag which tell the solver if the rhs Vector is
449  /// distributed or not
451 
452  /// Internal flag which tell the solver if the solution Vector to
453  /// be returned is distributed or not
455 #endif
456 
457  private:
458  /// Hypre copies matrix data from oomph-lib's CRDoubleMatrix
459  /// into its own data structures, doubling the memory requirements.
460  /// As far as the Hypre solvers are concerned the oomph-lib matrix
461  /// is no longer required and could be deleted to save memory.
462  /// However, this will not be the expected behaviour and therefore
463  /// needs to be requested explicitly by the user by changing this
464  /// flag from false (its default) to true.
466 
467  /// The Hypre_IJMatrix version of the matrix used in solve(...),
468  /// resolve(...) or preconditioner_solve(...).
469  HYPRE_IJMatrix Matrix_ij;
470 
471  /// The Hypre_ParCSRMatrix version of the matrix used in solve(...),
472  /// resolve(...) or preconditioner_solve(...).
473  HYPRE_ParCSRMatrix Matrix_par;
474 
475  /// The Hypre solver used in solve(...), resolve(...) or
476  /// preconditioner_solve(...). [This is a C structure!]
477  HYPRE_Solver Solver;
478 
479  /// The internal Hypre preconditioner used in conjunction with
480  /// Solver. [This is a C structure!]
481  HYPRE_Solver Preconditioner;
482 
483  /// Used to keep track of which solver (if any) is currently stored.
484  unsigned Existing_solver;
485 
486  /// Used to keep track of which preconditioner (if any) is currently stored.
488 
489  /// the distribution for this helpers-
491  };
492 
493 
494  /// /////////////////////////////////////////////////////////////////////
495  /// /////////////////////////////////////////////////////////////////////
496  /// /////////////////////////////////////////////////////////////////////
497 
498  //=====================================================================
499  /// An LinearSolver class using the suite of Hypre solvers
500  /// to allow
501  ///
502  /// BoomerAMG (AMG), CG, GMRES or BiCGStab
503  ///
504  /// to be used for LinearSolver::solve(...) or LinearSolver::resolve(...).
505  ///
506  /// The Krylov subspace solvers (CG, GMRES and BiCGStab) may be
507  /// preconditioned using:
508  ///
509  /// BoomerAMG, Euclid or ParaSails
510  ///
511  /// By default GMRES without preconditioning is used.
512  //====================================================================
513  class HypreSolver : public LinearSolver, public HypreInterface
514  {
515  public:
516  /// Constructor
518  {
519  // Hypre copies matrix data from oomph-lib's CRDoubleMatrix
520  // and Distributed CRDoubleMatrix into its own data structures,
521  // doubling the memory requirements for the matrix.
522  // As far as the Hypre solvers are concerned the oomph-lib matrix
523  // is no longer required and could be deleted to save memory.
524  // However, this will not be the expected behaviour and therefore
525  // needs to be requested explicitly by the user by changing this
526  // flag from false (its default) to true.
527  Delete_matrix = false;
528 
529  // Do we want to output results of timings?
530  Doc_time = true;
531  }
532 
533  /// Empty destructor.
535 
536  /// Broken copy constructor.
537  HypreSolver(const HypreSolver&) = delete;
538 
539  /// Broken assignment operator.
540  void operator=(const HypreSolver&) = delete;
541 
542  /// Disable resolve function (overloads the LinearSolver
543  /// disable_resolve function).
545  {
546  Enable_resolve = false;
547  clean_up_memory();
548  }
549 
550  /// Access function to Max_iter
551  unsigned& max_iter()
552  {
553  return Max_iter;
554  }
555 
556  /// Access function to Tolerance
557  double& tolerance()
558  {
559  return Tolerance;
560  }
561 
562  /// Access function to Hypre_method flag -- specified via enumeration.
563  unsigned& hypre_method()
564  {
565  return Hypre_method;
566  }
567 
568  /// Access function to Internal_preconditioner flag -- specified
569  /// via enumeration.
571  {
573  }
574 
575  /// Function to select use of 'simple' AMG smoothers as controlled
576  /// by AMG_simple_smoother flag
578  {
580  }
581 
582  /// Access function to AMG_simple_smoother flag
584  {
585  return AMG_simple_smoother;
586  }
587 
588  /// Function to select use of 'complex' AMG smoothers as controlled
589  /// by AMG_complex_smoother flag
591  {
593  }
594 
595  /// Access function to AMG_complex_smoother flag
597  {
598  return AMG_complex_smoother;
599  }
600 
601  /// Access function to AMG_print_level
602  unsigned& amg_print_level()
603  {
604  return AMG_print_level;
605  }
606 
607  /// Access function to AMG_smoother_iterations
609  {
611  }
612 
613  /// Access function to AMG_coarsening flag
614  unsigned& amg_coarsening()
615  {
616  return AMG_coarsening;
617  }
618 
619  /// Access function to AMG_max_levels
620  unsigned& amg_max_levels()
621  {
622  return AMG_max_levels;
623  }
624 
625  /// Access function to AMG_damping parameter
626  double& amg_damping()
627  {
628  return AMG_damping;
629  }
630 
631  /// Access function to AMG_strength
632  double& amg_strength()
633  {
634  return AMG_strength;
635  }
636 
637  /// Access function to AMG_max_row_sum
638  double& amg_max_row_sum()
639  {
640  return AMG_max_row_sum;
641  }
642 
643  /// Access function to AMG_truncation
644  double& amg_truncation()
645  {
646  return AMG_truncation;
647  }
648 
649  /// Access function to ParaSails symmetry flag
651  {
652  return ParaSails_symmetry;
653  }
654 
655  /// Access function to ParaSails nlevel parameter
657  {
658  return ParaSails_nlevel;
659  }
660 
661  /// Access function to ParaSails thresh parameter
663  {
664  return ParaSails_thresh;
665  }
666 
667  /// Access function to ParaSails filter parameter
669  {
670  return ParaSails_filter;
671  }
672 
673  /// Access function to Euclid drop tolerance parameter
674  double& euclid_droptol()
675  {
676  return Euclid_droptol;
677  }
678 
679  /// Access function to Euclid level parameter
680  /// for ILU(k) factorization
682  {
683  return Euclid_level;
684  }
685 
686  /// Enable euclid rowScaling
688  {
689  Euclid_rowScale = true;
690  }
691 
692  /// Disable euclid row scaling
694  {
695  Euclid_rowScale = false;
696  }
697 
698  /// Enable use of Block Jacobi
699  /// as opposed to PILU
701  {
702  Euclid_using_BJ = true;
703  }
704 
705  /// Disable use of Block Jacobi,
706  /// so PILU will be used
708  {
709  Euclid_using_BJ = false;
710  }
711 
712  /// Function to switch on ILU(k) factorization for Euclid
713  /// (default is ILU(k) factorization)
715  {
716  Euclid_using_ILUT = false;
717  }
718 
719  /// Function to switch on ILUT factorization for Euclid
720  /// (default is ILU(k) factorization)
722  {
723  Euclid_using_ILUT = true;
724  }
725 
726  /// Function to set the level of printing from Euclid
727  /// when the Euclid destructor is called
728  /// 0: no printing (default)
729  /// 1: prints summary of runtime settings and timings
730  /// 2: as 1 plus prints memory usage
731  unsigned& euclid_print_level()
732  {
733  return Euclid_print_level;
734  }
735 
736  /// Access function to Krylov_print_level
737  unsigned& krylov_print_level()
738  {
739  return Krylov_print_level;
740  }
741 
742  /// Hypre copies matrix data from oomph-lib's CRDoubleMatrix
743  /// or DistributedCRDoubleMatrix into its own data structures,
744  /// doubling the memory requirements for the matrix.
745  /// As far as the Hypre solvers are concerned the oomph-lib matrix
746  /// is no longer required and could be deleted to save memory.
747  /// However, this will not be the expected behaviour and therefore
748  /// needs to be requested explicitly by the user by calling this function
749  /// which changes the
750  /// flag from false (its default) to true.
752  {
753  Delete_matrix = true;
754  }
755 
756  /// Hypre copies matrix data from oomph-lib's CRDoubleMatrix
757  /// or DistributedCRDoubleMatrix into its own data structures,
758  /// doubling the memory requirements for the matrix.
759  /// Calling this function ensures that the matrix is not deleted.
761  {
762  Delete_matrix = false;
763  }
764 
765 
766  /// Function which uses problem_pt's get_jacobian(...) function to
767  /// generate a linear system which is then solved. This function deletes
768  /// any existing internal data and then generates a new Hypre solver.
769  void solve(Problem* const& problem_pt, DoubleVector& solution);
770 
771  /// Function to solve the linear system defined by matrix_pt
772  /// and rhs. This function will delete any existing internal data
773  /// and generate a new Hypre solver.
774  /// \b Note: The matrix has to be of type CRDoubleMatrix or
775  /// Distributed CRDoubleMatrix.
776  /// \b Note: Hypre copies matrix data from oomph-lib's CRDoubleMatrix
777  /// or DistributedCRDoubleMatrix into its own data structures,
778  /// doubling the memory requirements for the matrix.
779  /// As far as the Hypre solvers are concerned, the oomph-lib matrix
780  /// is no longer required and could be deleted to save memory.
781  /// However, this will not be the expected behaviour and therefore
782  /// needs to be requested explicitly by the user by calling the
783  /// enable_delete_matrix() function.
784  void solve(DoubleMatrixBase* const& matrix_pt,
785  const DoubleVector& rhs,
786  DoubleVector& solution);
787 
788  /// Function to resolve a linear system using the existing solver
789  /// data, allowing a solve with a new right hand side vector. This
790  /// function must be used after a call to solve(...) with
791  /// enable_resolve set to true.
792  void resolve(const DoubleVector& rhs, DoubleVector& solution);
793 
794  /// Function deletes all solver data.
795  void clean_up_memory();
796 
797  private:
798  /// Hypre copies matrix data from oomph-lib's CRDoubleMatrix
799  /// or DistributedCRDoubleMatrix into its own data structures,
800  /// doubling the memory requirements for the matrix.
801  /// As far as the Hypre solvers are concerned the oomph-lib matrix
802  /// is no longer required and could be deleted to save memory.
803  /// However, this will not be the expected behaviour and therefore
804  /// needs to be requested explicitly by the user by changing this
805  /// flag from false (its default) to true.
807  };
808 
809 
810  /// /////////////////////////////////////////////////////////////////////
811  /// /////////////////////////////////////////////////////////////////////
812  /// /////////////////////////////////////////////////////////////////////
813 
814  //=====================================================================
815  /// An Preconditioner class using the suite of Hypre preconditioners
816  /// to allow
817  ///
818  /// BoomerAMG (Algebraic Multi Grid),
819  /// Euclid (ILU) or
820  /// ParaSails (Approximate inverse)
821  ///
822  /// to be used for Preconditioner::preconditioner_solve(...).
823  /// By default BoomerAMG is used.
824  //====================================================================
826  {
827  public:
828  /// Constructor. Provide optional string that is used
829  /// in annotation of performance
830  HyprePreconditioner(const std::string& context_string = "")
831  {
832  Context_string = context_string;
833 
834  // Hypre copies matrix data from oomph-lib's CRDoubleMatrix
835  // or DistributedCRDoubleMatrix into its own data structures,
836  // doubling the memory requirements for the matrix.
837  // As far as the Hypre solvers are concerned the oomph-lib matrix
838  // is no longer required and could be deleted to save memory.
839  // However, this will not be the expected behaviour and therefore
840  // needs to be requested explicitly by the user by changing this
841  // flag from false (its default) to true.
842  Delete_matrix = false;
843 
844  // Do we want to output results of timings?
845  Doc_time = true;
846 
847  // General control paramaters for using HYPRE as preconditioner
848  Tolerance = 0.0;
849  Max_iter = 1;
852 
853  // Initialise private double that accumulates the preconditioner
854  // solve time of thi instantiation of this class. Is reported
855  // in destructor if Report_my_cumulative_preconditioner_solve_time
856  // is set to true
859  }
860 
861  /// Destructor.
863  {
865  {
866  oomph_info << "~HyprePreconditioner() in context \" " << Context_string
867  << "\": My_cumulative_preconditioner_solve_time = "
869  }
870  }
871 
872  /// Broken copy constructor.
874 
875  /// Broken assignment operator.
876  void operator=(const HyprePreconditioner&) = delete;
877 
878  /// Static double that accumulates the preconditioner
879  /// solve time of all instantiations of this class. Reset
880  /// it manually, e.g. after every Newton solve, using
881  /// reset_cumulative_solve_times().
883 
884  /// Static double that accumulates the preconditioner
885  /// solve time of all instantiations of this class, labeled by
886  /// context string. Reset
887  /// it manually, e.g. after every Newton solve, using
888  /// reset_cumulative_solve_times().
889  static std::map<std::string, double> Context_based_cumulative_solve_time;
890 
891 
892  /// Static unsigned that accumulates the number of preconditioner
893  /// solves of all instantiations of this class. Reset
894  /// it manually, e.g. after every Newton solve, using
895  /// reset_cumulative_solve_times().
897 
898  /// Static unsigned that accumulates the number of preconditioner
899  /// solves of all instantiations of this class, labeled by
900  /// context string. Reset
901  /// it manually, e.g. after every Newton solve, using
902  /// reset_cumulative_solve_times().
903  static std::map<std::string, unsigned>
905 
906  /// Static unsigned that stores nrow for the most recent
907  /// instantiations of this class, labeled by
908  /// context string. Reset
909  /// it manually, e.g. after every Newton solve, using
910  /// reset_cumulative_solve_times().
911  static std::map<std::string, unsigned> Context_based_nrow;
912 
913  /// Report cumulative solve times of all instantiations of this
914  /// class
915  static void report_cumulative_solve_times();
916 
917  /// Reset cumulative solve times
918  static void reset_cumulative_solve_times();
919 
920  /// Enable reporting of cumulative solve time in destructor
922  {
924  }
925 
926  /// Disable reporting of cumulative solve time in destructor
928  {
930  }
931 
932  /// Enable documentation of preconditioner timings
934  {
935  Doc_time = true;
936  }
937 
938  /// Disable documentation of preconditioner timings
940  {
941  Doc_time = false;
942  }
943 
944  /// Access function to Hypre_method flag -- specified via enumeration.
945  unsigned& hypre_method()
946  {
947  return Hypre_method;
948  }
949 
950  /// Access function to Internal_preconditioner flag -- specified
951  /// via enumeration.
953  {
955  }
956 
957  /// Function to select BoomerAMG as the preconditioner
959  {
961  }
962 
963  /// Function to set the number of times to apply BoomerAMG
964  void set_amg_iterations(const unsigned& amg_iterations)
965  {
967  }
968 
969  /// Return function for Max_iter
970  unsigned& amg_iterations()
971  {
972  return Max_iter;
973  }
974 
975  /// Function to select use of 'simple' AMG smoothers as controlled
976  /// by the flag AMG_simple_smoother
978  {
980  }
981 
982  /// Access function to AMG_simple_smoother flag
984  {
985  return AMG_simple_smoother;
986  }
987 
988  /// Function to select use of 'complex' AMG smoothers as controlled
989  /// by the flag AMG_complex_smoother
991  {
993  }
994 
995  /// Access function to AMG_complex_smoother flag
997  {
998  return AMG_complex_smoother;
999  }
1000 
1001  /// Return function for the AMG_using_simple_smoothing_flag
1003  {
1005  }
1006 
1007  /// Access function to AMG_print_level
1008  unsigned& amg_print_level()
1009  {
1010  return AMG_print_level;
1011  }
1012 
1013  /// Access function to AMG_smoother_iterations
1015  {
1016  return AMG_smoother_iterations;
1017  }
1018 
1019  /// Access function to AMG_coarsening flag
1020  unsigned& amg_coarsening()
1021  {
1022  return AMG_coarsening;
1023  }
1024 
1025  /// Access function to AMG_max_levels
1026  unsigned& amg_max_levels()
1027  {
1028  return AMG_max_levels;
1029  }
1030 
1031  /// Access function to AMG_damping parameter
1032  double& amg_damping()
1033  {
1034  return AMG_damping;
1035  }
1036 
1037  /// Access function to AMG_strength
1038  double& amg_strength()
1039  {
1040  return AMG_strength;
1041  }
1042 
1043  /// Access function to AMG_max_row_sum
1045  {
1046  return AMG_max_row_sum;
1047  }
1048 
1049  /// Access function to AMG_truncation
1050  double& amg_truncation()
1051  {
1052  return AMG_truncation;
1053  }
1054 
1055  /// Function to select ParaSails as the preconditioner
1057  {
1059  }
1060 
1061  /// Access function to ParaSails symmetry flag
1063  {
1064  return ParaSails_symmetry;
1065  }
1066 
1067  /// Access function to ParaSails nlevel parameter
1069  {
1070  return ParaSails_nlevel;
1071  }
1072 
1073  /// Access function to ParaSails thresh parameter
1075  {
1076  return ParaSails_thresh;
1077  }
1078 
1079  /// Access function to ParaSails filter parameter
1081  {
1082  return ParaSails_filter;
1083  }
1084 
1085  /// Function to select use Euclid as the preconditioner
1086  void use_Euclid()
1087  {
1088  Hypre_method = Euclid;
1089  }
1090 
1091  /// Access function to Euclid drop tolerance parameter
1092  double& euclid_droptol()
1093  {
1094  return Euclid_droptol;
1095  }
1096 
1097  /// Access function to Euclid level parameter for ILU(k) factorization
1099  {
1100  return Euclid_level;
1101  }
1102 
1103  /// Enable euclid rowScaling
1105  {
1106  Euclid_rowScale = true;
1107  }
1108 
1109  /// Disable euclid row scaling
1111  {
1112  Euclid_rowScale = false;
1113  }
1114 
1115  /// Enable use of Block Jacobi
1116  /// as opposed to PILU
1118  {
1119  Euclid_using_BJ = true;
1120  }
1121 
1122  /// Disable use of Block Jacobi,
1123  /// so PILU will be used
1125  {
1126  Euclid_using_BJ = false;
1127  }
1128 
1129  /// Function to switch on ILU(k) factorization for Euclid
1130  /// (default is ILU(k) factorization)
1132  {
1133  Euclid_using_ILUT = false;
1134  }
1135 
1136  /// Function to switch on ILUT factorization for Euclid
1137  /// (default is ILU(k) factorization)
1139  {
1140  Euclid_using_ILUT = true;
1141  }
1142 
1143  /// Function to set the level of printing from Euclid
1144  /// when the Euclid destructor is called
1145  /// 0: no printing (default)
1146  /// 1: prints summary of runtime settings and timings
1147  /// 2: as 1 plus prints memory usage
1149  {
1150  return Euclid_print_level;
1151  }
1152 
1153  /// Hypre copies matrix data from oomph-lib's CRDoubleMatrix
1154  /// or DistributedCRDoubleMatrix into its own data structures,
1155  /// doubling the memory requirements for the matrix.
1156  /// As far as the Hypre solvers are concerned the oomph-lib matrix
1157  /// is no longer required and could be deleted to save memory.
1158  /// However, this will not be the expected behaviour and therefore
1159  /// needs to be requested explicitly by the user by calling this function
1160  /// which changes the
1161  /// flag from false (its default) to true.
1163  {
1164  Delete_matrix = true;
1165  }
1166 
1167  /// Hypre copies matrix data from oomph-lib's CRDoubleMatrix
1168  /// or DistributedCRDoubleMatrix into its own data structures,
1169  /// doubling the memory requirements for the matrix.
1170  /// Calling this function ensures that the matrix is not deleted.
1172  {
1173  Delete_matrix = false;
1174  }
1175 
1176  /// Function to set up a preconditioner for the linear
1177  /// system defined by matrix_pt. This function is required when
1178  /// preconditioning and must be called before using the
1179  /// preconditioner_solve(...) function. This interface allows
1180  /// HyprePreconditioner to be used as a Preconditioner object
1181  /// for oomph-lib's own IterativeLinearSolver class.
1182  /// \b Note: matrix_pt must point to an object of type
1183  /// CRDoubleMatrix or DistributedCRDoubleMatrix.
1184  /// \b Note: Hypre copies matrix data from oomph-lib's CRDoubleMatrix
1185  /// and DistributedCRDoubleMatrix into its own data structures,
1186  /// doubling the memory requirements for the matrix.
1187  /// As far as the Hypre solvers are concerned, the oomph-lib matrix
1188  /// is no longer required and could be deleted to save memory.
1189  /// However, this will not be the expected behaviour and therefore
1190  /// needs to be requested explicitly by the user by calling the
1191  /// enable_delete_matrix() function.
1192  void setup();
1193 
1194  /// Function applies solver to vector r for preconditioning.
1195  /// This requires a call to setup(...) first.
1196  /// \b Note: Hypre copies matrix data from oomph-lib's CRDoubleMatrix
1197  /// or DistributedCRDoubleMatrix into its own data structures,
1198  /// doubling the memory requirements for the matrix.
1199  /// As far as the Hypre solvers are concerned, the oomph-lib matrix
1200  /// is no longer required and could be deleted to save memory.
1201  /// However, this will not be the expected behaviour and therefore
1202  /// needs to be requested explicitly by the user with the
1203  /// delete_matrix() function.
1204  void preconditioner_solve(const DoubleVector& r, DoubleVector& z);
1205 
1206  /// Function deletes all solver data.
1207  void clean_up_memory();
1208 
1209  private:
1210  /// Hypre copies matrix data from oomph-lib's CRDoubleMatrix
1211  /// or DistributedCRDoubleMatrix into its own data structures,
1212  /// doubling the memory requirements for the matrix.
1213  /// As far as the Hypre solvers are concerned the oomph-lib matrix
1214  /// is no longer required and could be deleted to save memory.
1215  /// However, this will not be the expected behaviour and therefore
1216  /// needs to be requested explicitly by the user by changing this
1217  /// flag from false (its default) to true.
1219 
1220  // Flag is true to output results of timings
1221  bool Doc_time;
1222 
1223  /// Private double that accumulates the preconditioner
1224  /// solve time of thi instantiation of this class. Is reported
1225  /// in destructor if Report_my_cumulative_preconditioner_solve_time
1226  /// is set to true
1228 
1229  /// Bool to request report of My_cumulative_preconditioner_solve_time
1230  /// in destructor
1232 
1233  /// String can be used to provide context for annotation
1235  };
1236 
1237 
1238  /// /////////////////////////////////////////////////////////////////////
1239  /// /////////////////////////////////////////////////////////////////////
1240  /// /////////////////////////////////////////////////////////////////////
1241 
1242  //==================================================================
1243  /// Default settings for various uses of the Hypre solver
1244  //==================================================================
1245  namespace Hypre_default_settings
1246  {
1247  /// Set default parameters for use as preconditioner in
1248  /// for momentum block in Navier-Stokes problem
1250  HyprePreconditioner* hypre_preconditioner_pt);
1251 
1252  /// Set default parameters for use as preconditioner in
1253  /// 2D Poisson-type problem.
1255  HyprePreconditioner* hypre_preconditioner_pt);
1256 
1257  /// Set default parameters for use as preconditioner in
1258  /// 3D Poisson-type problem.
1260  HyprePreconditioner* hypre_preconditioner_pt);
1261 
1262  } // namespace Hypre_default_settings
1263 } // namespace oomph
1264 
1265 #endif
e
Definition: cfortran.h:571
A class for compressed row matrices. This is a distributable object.
Definition: matrices.h:888
Abstract base class for matrices of doubles – adds abstract interfaces for solving,...
Definition: matrices.h:261
A vector in the mathematical sense, initially developed for linear algebra type applications....
Definition: double_vector.h:58
An interface class to the suite of Hypre solvers and preconditioners to allow use of:
Definition: hypre_solver.h:134
bool Output_info
Flag is true to output info and results of timings.
Definition: hypre_solver.h:305
bool Delete_input_data
Internal flag which is true when hypre_setup or hypre_solve can delete input matrix.
Definition: hypre_solver.h:445
unsigned AMG_coarsening
AMG coarsening strategy. Coarsening types include: 0 = CLJP (parallel coarsening using independent se...
Definition: hypre_solver.h:386
int Euclid_level
Euclid level parameter for ILU(k) factorization.
Definition: hypre_solver.h:418
double AMGEuclidSmoother_drop_tol
Definition: hypre_solver.h:283
bool Euclid_using_BJ
Flag to determine if Block Jacobi is used instead of PILU.
Definition: hypre_solver.h:415
bool Euclid_rowScale
Flag to switch on Euclid row scaling.
Definition: hypre_solver.h:409
bool Returning_distributed_solution
Internal flag which tell the solver if the solution Vector to be returned is distributed or not.
Definition: hypre_solver.h:454
double AMG_strength
Connection strength threshold parameter for BoomerAMG.
Definition: hypre_solver.h:369
HypreInterface(const HypreInterface &)=delete
Broken copy constructor.
void enable_hypre_error_messages()
Turn on the hypre error messages.
Definition: hypre_solver.h:242
LinearAlgebraDistribution * Hypre_distribution_pt
the distribution for this helpers-
Definition: hypre_solver.h:490
int ParaSails_nlevel
ParaSails nlevel parameter.
Definition: hypre_solver.h:397
double ParaSails_filter
ParaSails filter parameter.
Definition: hypre_solver.h:403
bool AMGEuclidSmoother_use_row_scaling
Definition: hypre_solver.h:280
unsigned AMGEuclidSmoother_level
Definition: hypre_solver.h:282
HYPRE_Solver Solver
The Hypre solver used in solve(...), resolve(...) or preconditioner_solve(...). [This is a C structur...
Definition: hypre_solver.h:477
int ParaSails_symmetry
ParaSails symmetry flag, used to inform ParaSails of Symmetry of definitenss of problem and type of P...
Definition: hypre_solver.h:394
bool Using_distributed_rhs
Internal flag which tell the solver if the rhs Vector is distributed or not.
Definition: hypre_solver.h:450
HypreInterface()
Constructor.
Definition: hypre_solver.h:137
unsigned existing_preconditioner()
Function return value of which preconditioner (if any) is stored.
Definition: hypre_solver.h:273
HYPRE_ParCSRMatrix Matrix_par
The Hypre_ParCSRMatrix version of the matrix used in solve(...), resolve(...) or preconditioner_solve...
Definition: hypre_solver.h:473
void operator=(const HypreInterface &)=delete
Broken assignment operator.
void disable_hypre_error_messages()
Turn off hypre error messages.
Definition: hypre_solver.h:248
unsigned Internal_preconditioner
Preconditioner method flag used with Hypre's PCG, GMRES and BiCGStab in solve(...) or resolve(....
Definition: hypre_solver.h:320
void hypre_solve(const DoubleVector &rhs, DoubleVector &solution)
Helper function performs a solve if any solver exists.
Hypre_method_types
Enumerated flag to define which Hypre methods are used CAREFUL: DON'T CHANGE THE ORDER OF THESE!
Definition: hypre_solver.h:256
~HypreInterface()
Destructor.
Definition: hypre_solver.h:226
HYPRE_Solver Preconditioner
The internal Hypre preconditioner used in conjunction with Solver. [This is a C structure!...
Definition: hypre_solver.h:481
unsigned AMGEuclidSmoother_print_level
Definition: hypre_solver.h:284
bool AMGEuclidSmoother_use_block_jacobi
Definition: hypre_solver.h:279
unsigned AMG_complex_smoother
Complex smoothing methods used in BoomerAMG. Relaxation types are: 6 = Schwarz 7 = Pilut 8 = ParaSail...
Definition: hypre_solver.h:360
bool Delete_matrix
Hypre copies matrix data from oomph-lib's CRDoubleMatrix into its own data structures,...
Definition: hypre_solver.h:465
double ParaSails_thresh
ParaSails thresh parameter.
Definition: hypre_solver.h:400
double AMG_truncation
Interpolation truncation factor for BoomerAMG.
Definition: hypre_solver.h:372
unsigned Krylov_print_level
Used to set the Hypre printing level for the Krylov subspace solvers.
Definition: hypre_solver.h:429
unsigned AMG_smoother_iterations
The number of smoother iterations to apply.
Definition: hypre_solver.h:363
unsigned existing_solver()
Function to return value of which solver (if any) is currently stored.
Definition: hypre_solver.h:267
void hypre_clean_up_memory()
Function deletes all solver data.
bool Hypre_error_messages
Flag to determine if non-zero values of the Hypre error flag plus Hypre error messages are output to ...
Definition: hypre_solver.h:441
double Tolerance
Tolerance used to terminate solver.
Definition: hypre_solver.h:311
bool AMG_using_simple_smoothing
Flag to determine whether simple smoothers (determined by the AMG_simple_smoother flag) or complex sm...
Definition: hypre_solver.h:338
unsigned Existing_solver
Used to keep track of which solver (if any) is currently stored.
Definition: hypre_solver.h:484
unsigned Max_iter
Maximum number of iterations used in solver.
Definition: hypre_solver.h:308
double AMG_damping
Damping factor for BoomerAMG smoothed Jacobi or hybrid SOR.
Definition: hypre_solver.h:366
unsigned AMG_simple_smoother
Simple smoothing methods used in BoomerAMG. Relaxation types include: 0 = Jacobi 1 = Gauss-Seidel,...
Definition: hypre_solver.h:351
unsigned AMG_print_level
Used to set the Hypre printing level for AMG 0: no printout 1: print setup information 2: print solve...
Definition: hypre_solver.h:327
HYPRE_IJMatrix Matrix_ij
The Hypre_IJMatrix version of the matrix used in solve(...), resolve(...) or preconditioner_solve(....
Definition: hypre_solver.h:469
unsigned Euclid_print_level
Flag to set the level of printing from Euclid when the Euclid destructor is called 0: no printing (de...
Definition: hypre_solver.h:425
unsigned AMG_max_levels
Maximum number of levels used in AMG.
Definition: hypre_solver.h:330
unsigned Hypre_method
Hypre method flag. Valid values are specified in enumeration.
Definition: hypre_solver.h:314
void hypre_matrix_setup(CRDoubleMatrix *matrix_pt)
Function which sets values of First_global_row, Last_global_row and other partitioning data and creat...
double AMG_max_row_sum
Parameter to identify diagonally dominant parts of the matrix in AMG.
Definition: hypre_solver.h:333
double Euclid_droptol
Euclid drop tolerance for ILU(k) and ILUT factorization.
Definition: hypre_solver.h:406
unsigned Existing_preconditioner
Used to keep track of which preconditioner (if any) is currently stored.
Definition: hypre_solver.h:487
bool Euclid_using_ILUT
Flag to determine if ILUT (if true) or ILU(k) is used in Euclid.
Definition: hypre_solver.h:412
void hypre_solver_setup()
Sets up the data required for to use as an oomph-lib LinearSolver or Preconditioner....
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
Definition: hypre_solver.h:826
~HyprePreconditioner()
Destructor.
Definition: hypre_solver.h:862
unsigned & amg_smoother_iterations()
Access function to AMG_smoother_iterations.
unsigned & amg_simple_smoother()
Access function to AMG_simple_smoother flag.
Definition: hypre_solver.h:983
double & amg_max_row_sum()
Access function to AMG_max_row_sum.
static void reset_cumulative_solve_times()
Reset cumulative solve times.
double & amg_strength()
Access function to AMG_strength.
static std::map< std::string, unsigned > Context_based_cumulative_npreconditioner_solve
Static unsigned that accumulates the number of preconditioner solves of all instantiations of this cl...
Definition: hypre_solver.h:904
unsigned & amg_iterations()
Return function for Max_iter.
Definition: hypre_solver.h:970
void disable_euclid_using_BJ()
Disable use of Block Jacobi, so PILU will be used.
std::string Context_string
String can be used to provide context for annotation.
double & parasails_thresh()
Access function to ParaSails thresh parameter.
int & parasails_nlevel()
Access function to ParaSails nlevel parameter.
void enable_delete_matrix()
Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own da...
void disable_report_my_cumulative_preconditioner_solve_time()
Disable reporting of cumulative solve time in destructor.
Definition: hypre_solver.h:927
double & amg_damping()
Access function to AMG_damping parameter.
unsigned & amg_print_level()
Access function to AMG_print_level.
void set_amg_iterations(const unsigned &amg_iterations)
Function to set the number of times to apply BoomerAMG.
Definition: hypre_solver.h:964
void clean_up_memory()
Function deletes all solver data.
void disable_delete_matrix()
Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own da...
static double Cumulative_preconditioner_solve_time
Static double that accumulates the preconditioner solve time of all instantiations of this class....
Definition: hypre_solver.h:882
void euclid_using_ILUT()
Function to switch on ILUT factorization for Euclid (default is ILU(k) factorization)
void setup()
Function to set up a preconditioner for the linear system defined by matrix_pt. This function is requ...
void enable_euclid_rowScale()
Enable euclid rowScaling.
unsigned & hypre_method()
Access function to Hypre_method flag – specified via enumeration.
Definition: hypre_solver.h:945
void euclid_using_ILUK()
Function to switch on ILU(k) factorization for Euclid (default is ILU(k) factorization)
int & euclid_level()
Access function to Euclid level parameter for ILU(k) factorization.
bool Report_my_cumulative_preconditioner_solve_time
Bool to request report of My_cumulative_preconditioner_solve_time in destructor.
void amg_using_complex_smoothing()
Function to select use of 'complex' AMG smoothers as controlled by the flag AMG_complex_smoother.
Definition: hypre_solver.h:990
static std::map< std::string, double > Context_based_cumulative_solve_time
Static double that accumulates the preconditioner solve time of all instantiations of this class,...
Definition: hypre_solver.h:889
unsigned & amg_max_levels()
Access function to AMG_max_levels.
void enable_euclid_using_BJ()
Enable use of Block Jacobi as opposed to PILU.
double & amg_truncation()
Access function to AMG_truncation.
bool Delete_matrix
Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own da...
void enable_doc_time()
Enable documentation of preconditioner timings.
Definition: hypre_solver.h:933
void disable_doc_time()
Disable documentation of preconditioner timings.
Definition: hypre_solver.h:939
static std::map< std::string, unsigned > Context_based_nrow
Static unsigned that stores nrow for the most recent instantiations of this class,...
Definition: hypre_solver.h:911
void amg_using_simple_smoothing()
Function to select use of 'simple' AMG smoothers as controlled by the flag AMG_simple_smoother.
Definition: hypre_solver.h:977
bool & amg_using_simple_smoothing_flag()
Return function for the AMG_using_simple_smoothing_flag.
double My_cumulative_preconditioner_solve_time
Private double that accumulates the preconditioner solve time of thi instantiation of this class....
void use_Euclid()
Function to select use Euclid as the preconditioner.
double & euclid_droptol()
Access function to Euclid drop tolerance parameter.
unsigned & euclid_print_level()
Function to set the level of printing from Euclid when the Euclid destructor is called 0: no printing...
void disable_euclid_rowScale()
Disable euclid row scaling.
void operator=(const HyprePreconditioner &)=delete
Broken assignment operator.
void use_BoomerAMG()
Function to select BoomerAMG as the preconditioner.
Definition: hypre_solver.h:958
unsigned & internal_preconditioner()
Access function to Internal_preconditioner flag – specified via enumeration.
Definition: hypre_solver.h:952
static void report_cumulative_solve_times()
Report cumulative solve times of all instantiations of this class.
void use_ParaSails()
Function to select ParaSails as the preconditioner.
int & parasails_symmetry()
Access function to ParaSails symmetry flag.
unsigned & amg_complex_smoother()
Access function to AMG_complex_smoother flag.
Definition: hypre_solver.h:996
unsigned & amg_coarsening()
Access function to AMG_coarsening flag.
static unsigned Cumulative_npreconditioner_solve
Static unsigned that accumulates the number of preconditioner solves of all instantiations of this cl...
Definition: hypre_solver.h:896
void preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Function applies solver to vector r for preconditioning. This requires a call to setup(....
HyprePreconditioner(const HyprePreconditioner &)=delete
Broken copy constructor.
HyprePreconditioner(const std::string &context_string="")
Constructor. Provide optional string that is used in annotation of performance.
Definition: hypre_solver.h:830
double & parasails_filter()
Access function to ParaSails filter parameter.
void enable_report_my_cumulative_preconditioner_solve_time()
Enable reporting of cumulative solve time in destructor.
Definition: hypre_solver.h:921
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
Definition: hypre_solver.h:514
unsigned & amg_simple_smoother()
Access function to AMG_simple_smoother flag.
Definition: hypre_solver.h:583
void disable_euclid_using_BJ()
Disable use of Block Jacobi, so PILU will be used.
Definition: hypre_solver.h:707
int & euclid_level()
Access function to Euclid level parameter for ILU(k) factorization.
Definition: hypre_solver.h:681
void disable_euclid_rowScale()
Disable euclid row scaling.
Definition: hypre_solver.h:693
double & amg_damping()
Access function to AMG_damping parameter.
Definition: hypre_solver.h:626
unsigned & max_iter()
Access function to Max_iter.
Definition: hypre_solver.h:551
void enable_euclid_using_BJ()
Enable use of Block Jacobi as opposed to PILU.
Definition: hypre_solver.h:700
double & euclid_droptol()
Access function to Euclid drop tolerance parameter.
Definition: hypre_solver.h:674
unsigned & hypre_method()
Access function to Hypre_method flag – specified via enumeration.
Definition: hypre_solver.h:563
unsigned & euclid_print_level()
Function to set the level of printing from Euclid when the Euclid destructor is called 0: no printing...
Definition: hypre_solver.h:731
void solve(Problem *const &problem_pt, DoubleVector &solution)
Function which uses problem_pt's get_jacobian(...) function to generate a linear system which is then...
double & amg_max_row_sum()
Access function to AMG_max_row_sum.
Definition: hypre_solver.h:638
unsigned & amg_print_level()
Access function to AMG_print_level.
Definition: hypre_solver.h:602
HypreSolver(const HypreSolver &)=delete
Broken copy constructor.
~HypreSolver()
Empty destructor.
Definition: hypre_solver.h:534
unsigned & krylov_print_level()
Access function to Krylov_print_level.
Definition: hypre_solver.h:737
double & tolerance()
Access function to Tolerance.
Definition: hypre_solver.h:557
void clean_up_memory()
Function deletes all solver data.
unsigned & amg_coarsening()
Access function to AMG_coarsening flag.
Definition: hypre_solver.h:614
double & parasails_thresh()
Access function to ParaSails thresh parameter.
Definition: hypre_solver.h:662
HypreSolver()
Constructor.
Definition: hypre_solver.h:517
double & parasails_filter()
Access function to ParaSails filter parameter.
Definition: hypre_solver.h:668
unsigned & amg_complex_smoother()
Access function to AMG_complex_smoother flag.
Definition: hypre_solver.h:596
bool Delete_matrix
Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own da...
Definition: hypre_solver.h:806
void amg_using_complex_smoothing()
Function to select use of 'complex' AMG smoothers as controlled by AMG_complex_smoother flag.
Definition: hypre_solver.h:590
void euclid_using_ILUT()
Function to switch on ILUT factorization for Euclid (default is ILU(k) factorization)
Definition: hypre_solver.h:721
void amg_using_simple_smoothing()
Function to select use of 'simple' AMG smoothers as controlled by AMG_simple_smoother flag.
Definition: hypre_solver.h:577
void resolve(const DoubleVector &rhs, DoubleVector &solution)
Function to resolve a linear system using the existing solver data, allowing a solve with a new right...
void enable_euclid_rowScale()
Enable euclid rowScaling.
Definition: hypre_solver.h:687
unsigned & amg_max_levels()
Access function to AMG_max_levels.
Definition: hypre_solver.h:620
int & parasails_symmetry()
Access function to ParaSails symmetry flag.
Definition: hypre_solver.h:650
unsigned & amg_smoother_iterations()
Access function to AMG_smoother_iterations.
Definition: hypre_solver.h:608
unsigned & internal_preconditioner()
Access function to Internal_preconditioner flag – specified via enumeration.
Definition: hypre_solver.h:570
void euclid_using_ILUK()
Function to switch on ILU(k) factorization for Euclid (default is ILU(k) factorization)
Definition: hypre_solver.h:714
void operator=(const HypreSolver &)=delete
Broken assignment operator.
double & amg_truncation()
Access function to AMG_truncation.
Definition: hypre_solver.h:644
int & parasails_nlevel()
Access function to ParaSails nlevel parameter.
Definition: hypre_solver.h:656
void disable_resolve()
Disable resolve function (overloads the LinearSolver disable_resolve function).
Definition: hypre_solver.h:544
void enable_delete_matrix()
Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own da...
Definition: hypre_solver.h:751
void disable_delete_matrix()
Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own da...
Definition: hypre_solver.h:760
double & amg_strength()
Access function to AMG_strength.
Definition: hypre_solver.h:632
Describes the distribution of a distributable linear algebra type object. Typically this is a contain...
Base class for all linear solvers. This merely defines standard interfaces for linear solvers,...
Definition: linear_solver.h:68
bool Doc_time
Boolean flag that indicates whether the time taken.
Definition: linear_solver.h:77
bool Enable_resolve
Boolean that indicates whether the matrix (or its factors, in the case of direct solver) should be st...
Definition: linear_solver.h:73
static bool mpi_has_been_initialised()
return true if MPI has been initialised
static OomphCommunicator * communicator_pt()
access to global communicator. This is the oomph-lib equivalent of MPI_COMM_WORLD
An OomphLibError object which should be thrown when an run-time error is encountered....
Preconditioner base class. Gives an interface to call all other preconditioners through and stores th...
////////////////////////////////////////////////////////////////// //////////////////////////////////...
Definition: problem.h:151
int hypre__global_error
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
void create_HYPRE_Matrix(CRDoubleMatrix *oomph_matrix, HYPRE_IJMatrix &hypre_ij_matrix, HYPRE_ParCSRMatrix &hypre_par_matrix, LinearAlgebraDistribution *dist_pt)
Helper function to create a serial HYPRE_IJMatrix and HYPRE_ParCSRMatrix from a CRDoubleMatrix NOTE: ...
void create_HYPRE_Vector(const DoubleVector &oomph_vec, const LinearAlgebraDistribution *dist_pt, HYPRE_IJVector &hypre_ij_vector, HYPRE_ParVector &hypre_par_vector)
Helper function to create a HYPRE_IJVector and HYPRE_ParVector.
int check_HYPRE_error_flag(std::ostringstream &message)
Helper function to check the Hypre error flag, return the message associated with any error,...
unsigned AMG_coarsening
Default AMG coarsening strategy. Coarsening types include: 0 = CLJP (parallel coarsening using indepe...
double AMG_truncation
AMG interpolation truncation factor.
double AMG_strength
Default for AMG strength (0.25 recommended for 2D problems; larger (0.5-0.75, say) for 3D.
void euclid_settings_helper(const bool &use_block_jacobi, const bool &use_row_scaling, const bool &use_ilut, const int &level, const double &drop_tol, const int &print_level, HYPRE_Solver &euclid_object)
Helper function to set Euclid options using a command line like array.
void set_defaults_for_navier_stokes_momentum_block(HyprePreconditioner *hypre_preconditioner_pt)
Set default parameters for use as preconditioner in for momentum block in Navier-Stokes problem.
Definition: hypre_solver.cc:85
void set_defaults_for_3D_poisson_problem(HyprePreconditioner *hypre_preconditioner_pt)
Set default parameters for use as preconditioner in 3D Poisson-type problem.
Definition: hypre_solver.cc:72
void set_defaults_for_2D_poisson_problem(HyprePreconditioner *hypre_preconditioner_pt)
Set default parameters for use as preconditioner in 2D Poisson-type problem.
Definition: hypre_solver.cc:45
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
OomphInfo oomph_info
Single (global) instantiation of the OomphInfo object – this is used throughout the library as a "rep...