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-2025 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
43extern hypre_Error hypre__global_error;
44
45
46namespace oomph
47{
48 //==================================================================
49 /// Helper functions for use with the Hypre library
50 //==================================================================
51 namespace HypreHelpers
52 {
53 /// Number of active Hypre solvers (smart pointer like behaviuor
54 /// to make sure that the initialise/finalize functions are only
55 /// called the required number of times
56 extern unsigned Number_of_active_hypre_solvers;
57
58 /// Default for AMG strength (0.25 recommended for 2D problems;
59 /// larger (0.5-0.75, say) for 3D
60 extern double AMG_strength;
61
62 /// Default AMG coarsening strategy. Coarsening types include:
63 /// 0 = CLJP (parallel coarsening using independent sets)
64 /// 1 = classical RS with no boundary treatment (not recommended
65 /// in parallel)
66 /// 3 = modified RS with 3rd pass to add C points on the boundaries
67 /// 6 = Falgout (uses 1 then CLJP using interior coarse points as
68 /// first independent set)
69 /// 8 = PMIS (parallel coarsening using independent sets - lower
70 /// complexities than 0, maybe also slower convergence)
71 /// 10= HMIS (one pass RS on each processor then PMIS on interior
72 /// coarse points as first independent set)
73 /// 11= One pass RS on each processor (not recommended)
74 extern unsigned AMG_coarsening;
75
76 /// AMG interpolation truncation factor
77 extern double AMG_truncation; // 0.0;
78
79 /// Helper function to check the Hypre error flag, return the
80 /// message associated with any error, and reset the error flag to zero.
81 int check_HYPRE_error_flag(std::ostringstream& message);
82
83 /// Helper function to create a HYPRE_IJVector and HYPRE_ParVector.
84 /// + If no MPI then serial vectors are created
85 /// + If MPI and serial input vector then distributed hypre vectors are
86 /// created
87 /// + If MPI and distributed input vector the distributed output vectors
88 /// are created.
89 extern void create_HYPRE_Vector(const DoubleVector& oomph_vec,
90 const LinearAlgebraDistribution* dist_pt,
91 HYPRE_IJVector& hypre_ij_vector,
92 HYPRE_ParVector& hypre_par_vector);
93
94 /// Helper function to create an empty HYPRE_IJVector and
95 /// HYPRE_ParVector.
96 /// + If no MPI then serial vectors are created
97 /// + If MPI and serial distribution then distributed hypre vectors are
98 /// created
99 /// + If MPI and distributed input distribution the distributed output
100 /// vectors are created.
101 void create_HYPRE_Vector(const LinearAlgebraDistribution* oomph_vec,
102 HYPRE_IJVector& hypre_ij_vector,
103 HYPRE_ParVector& hypre_par_vector);
104
105 /// Helper function to create a serial HYPRE_IJMatrix and
106 /// HYPRE_ParCSRMatrix from a CRDoubleMatrix
107 void create_HYPRE_Matrix(CRDoubleMatrix* oomph_matrix,
108 HYPRE_IJMatrix& hypre_ij_matrix,
109 HYPRE_ParCSRMatrix& hypre_par_matrix,
110 LinearAlgebraDistribution* dist_pt);
111
112 /// Helper function to set Euclid options using a command line
113 /// like array.
114 void euclid_settings_helper(const bool& use_block_jacobi,
115 const bool& use_row_scaling,
116 const bool& use_ilut,
117 const int& level,
118 const double& drop_tol,
119 const int& print_level,
120 HYPRE_Solver& euclid_object);
121
122 } // namespace HypreHelpers
123
124
125 //=====================================================================
126 /// An interface class to the suite of Hypre solvers and preconditioners
127 /// to allow use of:
128 ///
129 /// BoomerAMG (AMG), CG, GMRES or BiCGStab, Euclid (ILU) or
130 /// ParaSails (Approximate inverse)
131 ///
132 /// Hypre's Krylov subspace solvers (CG, GMRES and BiCGStab) may be
133 /// preconditioned using:
134 ///
135 /// BoomerAMG, Euclid or ParaSails
136 ///
137 //====================================================================
139 {
140 public:
141 /// Constructor
143 {
144#ifdef PARANOID
145#ifndef HYPRE_SEQUENTIAL
146 // For the MPI version of Hypre, check MPI_Helpers::setup has been called
148 {
149 std::ostringstream error_message;
150 error_message << "When using the MPI version of Hypre please first\n"
151 << "call function MPI_Helpers::setup()\n";
152 throw OomphLibError(error_message.str(),
155 }
156#endif
157#endif
158
159 // Track number of instances of hypre solvers
161 {
163 }
165
166 // setup the distribution
168
169 // These keep track of which solver and preconditioner
170 // (if any) currently exist
173
174 // Do we want to output info and results of timings?
175 Output_info = true;
176
177 // General control paramaters
178 Tolerance = 1e-10;
179 Max_iter = 100;
182
183 // Default AMG control parameters -- these seem OK;
184 // see hypre documenation for details.
187 {
188 // Jacobi in parallel
190 }
191 else
192 {
193 // Gauss Seidel in serial
195 }
198 AMG_damping = 1.0;
202 AMG_max_levels = 100;
203 AMG_max_row_sum = 1.0;
204 AMG_print_level = 0;
205
206 // Parameters for using Euclicd as an AMG smoother (defaults copied
207 // from the normal defaults listed in the manual).
212 AMGEuclidSmoother_drop_tol = 0; // No dropping
214
215 // Print level for CG, GMRES and BiCGStab
217
218 // Set ParaSails default values
221 ParaSails_thresh = 0.1;
222 ParaSails_filter = 0.1;
223
224 // Set Euclid default values
225 Euclid_droptol = 0.0;
226 Euclid_rowScale = false;
227 Euclid_using_ILUT = false;
228 Euclid_using_BJ = false;
229 Euclid_level = 1;
231
232 // Set to true to periodically check the hypre error flag
233 // and output any messages
234 Hypre_error_messages = false;
235 }
236
237 /// Destructor.
239 {
240 // call function to delete solver data
242
243 // delete teh oomph-lib distribution
245
246 // Track number of instances of hypre solvers
248 {
250 }
251#ifdef PARANOID
253 {
254 // Just display the error; destructors aren't supposed to throw!
255 oomph_info << "ERROR in ~HypreInterface (where we shouldn't throw!):\n"
256 << "HypreHelpers::Number_of_active_hypre_solvers = 0 "
257 << std::endl;
258 }
259#endif
261 }
262
263 /// Broken copy constructor.
265
266 /// Broken assignment operator.
267 void operator=(const HypreInterface&) = delete;
268
269 /// Turn on the hypre error messages
271 {
273 }
274
275 /// Turn off hypre error messages
277 {
278 Hypre_error_messages = false;
279 }
280
281 /// Enumerated flag to define which Hypre methods are used
282 /// CAREFUL: DON'T CHANGE THE ORDER OF THESE!
293
294 /// Function to return value of which solver (if any) is currently stored.
296 {
297 return Existing_solver;
298 }
299
300 /// Function return value of which preconditioner (if any) is stored.
302 {
304 }
305
306
307 protected:
308 /// Function deletes all solver data.
310
311 /// Function which sets values of First_global_row,
312 /// Last_global_row and other partitioning data and creates the distributed
313 /// Hypre matrix (stored in Matrix_ij/Matrix_par) from the CRDoubleMatrix.
314 void hypre_matrix_setup(CRDoubleMatrix* matrix_pt);
315
316 /// Sets up the data required for to use as an oomph-lib
317 /// LinearSolver or Preconditioner. This must be called after
318 /// the Hypre matrix has been generated using hypre_matrix_setup(...).
319 void hypre_solver_setup();
320
321 /// Helper function performs a solve if any solver
322 /// exists.
324
325 /// Flag is true to output info and results of timings
327
328 /// Maximum number of iterations used in solver.
329 unsigned Max_iter;
330
331 /// Tolerance used to terminate solver.
332 double Tolerance;
333
334 /// Hypre method flag. Valid values are specified in enumeration
335 unsigned Hypre_method;
336
337 /// Preconditioner method flag used with Hypre's PCG,
338 /// GMRES and BiCGStab in solve(...) or resolve(...). Valid values
339 /// are BoomerAMG, Euclid, ParaSails or None (all enumerated above),
340 /// for any other value no preconditioner is set.
342
343 /// Used to set the Hypre printing level for AMG
344 /// 0: no printout
345 /// 1: print setup information
346 /// 2: print solve information
347 /// 3: print setup and solve information
349
350 /// Maximum number of levels used in AMG
352
353 /// Parameter to identify diagonally dominant parts of the matrix in AMG
355
356 /// Flag to determine whether simple smoothers (determined by the
357 /// AMG_simple_smoother flag) or complex smoothers (determined by the
358 /// AMG_complex_smoother flag are used in AMG
360
361 /// Simple smoothing methods used in BoomerAMG.
362 /// To use these methods set AMG_using_simple_smoothing to true.
363 /// Here are the current options (from hypre's HYPRE_parcsr_ls.h):
364 ///
365 /// (Optional) Defines the smoother to be used. It uses the given
366 /// smoother on the fine grid, the up and
367 /// the down cycle and sets the solver on the coarsest level to Gaussian
368 /// elimination (9). The default is \f$\ell_1\f$-Gauss-Seidel, forward solve (13)
369 /// on the down cycle and backward solve (14) on the up cycle.
370 ///
371 /// There are the following options for \e relax_type:
372 ///
373 /// - 0 : Jacobi
374 /// - 1 : Gauss-Seidel, sequential (very slow!)
375 /// - 2 : Gauss-Seidel, interior points in parallel, boundary
376 /// sequential (slow!)
377 /// - 3 : hybrid Gauss-Seidel or SOR, forward solve
378 /// - 4 : hybrid Gauss-Seidel or SOR, backward solve
379 /// - 5 : hybrid chaotic Gauss-Seidel (works only with OpenMP)
380 /// - 6 : hybrid symmetric Gauss-Seidel or SSOR
381 /// - 7 : Jacobi (uses Matvec)
382 /// - 8 : \f$\ell_1\f$-scaled hybrid symmetric Gauss-Seidel
383 /// - 9 : Gaussian elimination (only on coarsest level)
384 /// - 10 : On-processor direct forward solve for matrices with
385 /// triangular structure
386 /// - 11 : Two Stage approximation to GS. Uses the strict lower
387 /// part of the diagonal matrix
388 /// - 12 : Two Stage approximation to GS. Uses the strict lower
389 /// part of the diagonal matrix and a second iteration
390 /// for additional error approximation
391 /// - 13 : \f$\ell_1\f$ Gauss-Seidel, forward solve
392 /// - 14 : \f$\ell_1\f$ Gauss-Seidel, backward solve
393 /// - 15 : CG (warning - not a fixed smoother - may require FGMRES)
394 /// - 16 : Chebyshev
395 /// - 17 : FCF-Jacobi
396 /// - 18 : \f$\ell_1\f$-scaled jacobi
397 /// - 19 : Gaussian elimination (old version)
398 /// - 21 : The same as 8 except forcing serialization on CPU
399 /// ( # OMP-thread = 1)
400 /// - 29 : Direct solve: use Gaussian elimination & BLAS
401 /// (with pivoting) (old version)
402 /// - 30 : Kaczmarz
403 /// - 88: The same methods as 8 with a convergent l1-term
404 /// - 89: Symmetric l1-hybrid Gauss-Seidel (i.e., 13 followed by 14)
405 /// - 98 : LU with pivoting
406 /// - 99 : LU with pivoting
407 /// -199 : Matvec with the inverse
409
410 /// Complex smoothing methods used in BoomerAMG.
411 /// To use these methods set AMG_using_simple_smoothing to false
412 /// Here are the current options (from hypre's HYPRE_parcsr_ls.h):
413 ///
414 /// - 4 : FSAI (routines needed to set: HYPRE_BoomerAMGSetFSAIMaxSteps,
415 /// HYPRE_BoomerAMGSetFSAIMaxStepSize,
416 /// HYPRE_BoomerAMGSetFSAIEigMaxIters,
417 /// HYPRE_BoomerAMGSetFSAIKapTolerance)
418 /// - 5 : ParILUK (routines needed to set: HYPRE_ILUSetLevelOfFill,
419 /// HYPRE_ILUSetType)
420 /// - 6 : Schwarz (routines needed to set: HYPRE_BoomerAMGSetDomainType,
421 /// HYPRE_BoomerAMGSetOverlap, HYPRE_BoomerAMGSetVariant,
422 /// HYPRE_BoomerAMGSetSchwarzRlxWeight)
423 /// - 7 : Pilut (routines needed to set: HYPRE_BoomerAMGSetDropTol,
424 /// HYPRE_BoomerAMGSetMaxNzPerRow)
425 /// - 8 : ParaSails (routines needed to set: HYPRE_BoomerAMGSetSym,
426 /// HYPRE_BoomerAMGSetLevel, HYPRE_BoomerAMGSetFilter,
427 /// HYPRE_BoomerAMGSetThreshold)
428 /// - 9 : Euclid (routines needed to set: HYPRE_BoomerAMGSetEuclidFile)
429 ///
431
432 /// The number of smoother iterations to apply
434
435 /// Damping factor for BoomerAMG smoothed Jacobi or hybrid SOR
437
438 /// Connection strength threshold parameter for BoomerAMG
440
441 /// Interpolation truncation factor for BoomerAMG
443
444 /// AMG coarsening strategy. Coarsening types include:
445 /// 0 = CLJP (parallel coarsening using independent sets)
446 /// 1 = classical RS with no boundary treatment (not recommended
447 /// in parallel)
448 /// 3 = modified RS with 3rd pass to add C points on the boundaries
449 /// 6 = Falgout (uses 1 then CLJP using interior coarse points as
450 /// first independent set)
451 /// 8 = PMIS (parallel coarsening using independent sets - lower
452 /// complexities than 0, maybe also slower convergence)
453 /// 10= HMIS (one pass RS on each processor then PMIS on interior
454 /// coarse points as first independent set)
455 /// 11= One pass RS on each processor (not recommended)
457
458 /// ParaSails symmetry flag, used to inform ParaSails of
459 /// Symmetry of definitenss of problem and type of ParaSails
460 /// preconditioner:
461 /// 0 = nonsymmetric and/or indefinite problem, nonsymmetric preconditioner
462 /// 1 = SPD problem, and SPD (factored preconditioner)
463 /// 2 = nonsymmetric, definite problem and SDP (factored preconditoner)
465
466 /// ParaSails nlevel parameter
468
469 /// ParaSails thresh parameter
471
472 /// ParaSails filter parameter
474
475 /// Euclid drop tolerance for ILU(k) and ILUT factorization
477
478 /// Flag to switch on Euclid row scaling
480
481 /// Flag to determine if ILUT (if true) or ILU(k) is used in Euclid
483
484 /// Flag to determine if Block Jacobi is used instead of PILU
486
487 /// Euclid level parameter for ILU(k) factorization
489
490 /// Flag to set the level of printing from Euclid
491 /// when the Euclid destructor is called
492 /// 0: no printing (default)
493 /// 1: prints summary of runtime settings and timings
494 /// 2: as 1 plus prints memory usage
496
497
498 // If these are used comment and write access functions
505
506 /// Used to set the Hypre printing level for the Krylov
507 /// subspace solvers
509
510 /// Doc parameters used in hypre solver
512 {
513 oomph_info << " " << std::endl;
514 oomph_info << " " << std::endl;
515 oomph_info << "Hypre parameters:" << std::endl;
516 oomph_info << "=================" << std::endl;
517 oomph_info << "- Maximum number of iterations used in solver: Max_iter = "
518 << Max_iter << std::endl;
519
520 oomph_info << "- Tolerance used to terminate solver: Tolerance = "
521 << Tolerance << std::endl;
522
523 oomph_info << "- Hypre method flag. Valid values are specified in "
524 "enumeration: Hypre_method = "
525 << Hypre_method << std::endl;
526
527 oomph_info << "- Preconditioner method flag used with Hypre's PCG,"
528 << std::endl;
530 << " GMRES and BiCGStab in solve(...) or resolve(...). Valid values"
531 << std::endl;
533 << " are BoomerAMG, Euclid, ParaSails or None (all enumerated above),"
534 << std::endl;
535 oomph_info << " for any other value no preconditioner is set. "
536 << std::endl;
537 oomph_info << " Internal_preconditioner = " << Internal_preconditioner
538 << std::endl;
539
541 << "- Flag to set the Hypre printing level for AMG: AMG_print_level = "
542 << AMG_print_level << std::endl;
543
544 oomph_info << "- Maximum number of levels used in AMG: AMG_max_levels = "
545 << AMG_max_levels << std::endl;
546
547
548 oomph_info << "- Parameter to identify diagonally dominant parts of the "
549 "matrix in AMG:"
550 << std::endl;
551 oomph_info << " AMG_max_row_sum = " << AMG_max_row_sum << std::endl;
552
553
555 << "- Flag to determine whether simple smoothers (determined by the"
556 << std::endl;
558 << " AMG_simple_smoother flag) or complex smoothers (determined by the"
559 << std::endl;
560 oomph_info << " AMG_complex_smoother flag are used in AMG" << std::endl;
561 oomph_info << " AMG_using_simple_smoothing ="
562 << AMG_using_simple_smoothing << std::endl;
563
564 oomph_info << "- Simple smoothing method used in BoomerAMG. "
565 << std::endl;
566 oomph_info << " (only used if AMG_using_simple_smoothing is true"
567 << std::endl;
568 oomph_info << " AMG_simple_smoother = " << AMG_simple_smoother
569 << std::endl;
570
571
572 oomph_info << "- Complex smoothing method used in BoomerAMG. "
573 << std::endl;
574 oomph_info << " (only used if AMG_using_simple_smoothing is false"
575 << std::endl;
576 oomph_info << " AMG_complex_smoother = " << AMG_complex_smoother
577 << std::endl;
578
579 oomph_info << "- The number of smoother iterations to apply: "
580 "AMG_smoother_iterations = "
581 << AMG_smoother_iterations << std::endl;
582
584 << "- Damping factor for BoomerAMG smoothed Jacobi or hybrid SOR: "
585 << "AMG_damping = " << AMG_damping << std::endl;
586
587
588 oomph_info << "- Connection strength threshold parameter for BoomerAMG: "
589 "AMG_strength = "
590 << AMG_strength << std::endl;
591
592
594 << "- Interpolation truncation factor for BoomerAMG: AMG_truncation = "
595 << AMG_truncation << std::endl;
596
597 oomph_info << "- AMG coarsening strategy: AMG_coarsening = "
598 << AMG_coarsening << std::endl;
599
600
601 oomph_info << "NOTE: Skipping parameters for ParaSails and Euclid"
602 << std::endl;
603 oomph_info << " feel free to document these if you use thes solvers"
604 << std::endl;
605
606 oomph_info << "- Flag to control the Hypre printing level for the Krylov"
607 << std::endl;
608 oomph_info << " subspace solvers: Krylov_print_level = "
609 << Krylov_print_level << std::endl;
610
611
612 oomph_info << " " << std::endl;
613 oomph_info << " " << std::endl;
614 }
615
616
617 /// Flag to determine if non-zero values of the Hypre error flag
618 /// plus Hypre error messages are output to screen at various points
619 /// in the solve function, i.e. after:
620 /// 1. setting up the Hypre matrix
621 /// 2. setting up the preconditioner
622 /// 3. setting up the solver
623 /// 4. setting up the Hypre vectors
624 /// 5. solving
625 /// 6. getting the solution vector
626 /// 7. deallocation of solver data
628
629 /// Internal flag which is true when hypre_setup or hypre_solve
630 /// can delete input matrix.
632
633#ifdef OOMPH_HAS_MPI
634 /// Internal flag which tell the solver if the rhs Vector is
635 /// distributed or not
637
638 /// Internal flag which tell the solver if the solution Vector to
639 /// be returned is distributed or not
641#endif
642
643 private:
644 /// Hypre copies matrix data from oomph-lib's CRDoubleMatrix
645 /// into its own data structures, doubling the memory requirements.
646 /// As far as the Hypre solvers are concerned the oomph-lib matrix
647 /// is no longer required and could be deleted to save memory.
648 /// However, this will not be the expected behaviour and therefore
649 /// needs to be requested explicitly by the user by changing this
650 /// flag from false (its default) to true.
652
653 /// The Hypre_IJMatrix version of the matrix used in solve(...),
654 /// resolve(...) or preconditioner_solve(...).
656
657 /// The Hypre_ParCSRMatrix version of the matrix used in solve(...),
658 /// resolve(...) or preconditioner_solve(...).
660
661 /// The Hypre solver used in solve(...), resolve(...) or
662 /// preconditioner_solve(...). [This is a C structure!]
664
665 /// The internal Hypre preconditioner used in conjunction with
666 /// Solver. [This is a C structure!]
668
669 /// Used to keep track of which solver (if any) is currently stored.
671
672 /// Used to keep track of which preconditioner (if any) is currently stored.
674
675 /// the distribution for this helpers-
677 };
678
679
680 ////////////////////////////////////////////////////////////////////////
681 ////////////////////////////////////////////////////////////////////////
682 ////////////////////////////////////////////////////////////////////////
683
684 //=====================================================================
685 /// An LinearSolver class using the suite of Hypre solvers
686 /// to allow
687 ///
688 /// BoomerAMG (AMG), CG, GMRES or BiCGStab
689 ///
690 /// to be used for LinearSolver::solve(...) or LinearSolver::resolve(...).
691 ///
692 /// The Krylov subspace solvers (CG, GMRES and BiCGStab) may be
693 /// preconditioned using:
694 ///
695 /// BoomerAMG, Euclid or ParaSails
696 ///
697 /// By default GMRES without preconditioning is used.
698 //====================================================================
700 {
701 public:
702 /// Constructor
704 {
705 // Hypre copies matrix data from oomph-lib's CRDoubleMatrix
706 // and Distributed CRDoubleMatrix into its own data structures,
707 // doubling the memory requirements for the matrix.
708 // As far as the Hypre solvers are concerned the oomph-lib matrix
709 // is no longer required and could be deleted to save memory.
710 // However, this will not be the expected behaviour and therefore
711 // needs to be requested explicitly by the user by changing this
712 // flag from false (its default) to true.
713 Delete_matrix = false;
714
715 // Do we want to output results of timings?
716 Doc_time = true;
717 }
718
719 /// Empty destructor.
721
722 /// Broken copy constructor.
723 HypreSolver(const HypreSolver&) = delete;
724
725 /// Broken assignment operator.
726 void operator=(const HypreSolver&) = delete;
727
728 /// Disable resolve function (overloads the LinearSolver
729 /// disable_resolve function).
731 {
732 Enable_resolve = false;
734 }
735
736 /// Access function to Max_iter
737 unsigned& max_iter()
738 {
739 return Max_iter;
740 }
741
742 /// Access function to Tolerance
743 double& tolerance()
744 {
745 return Tolerance;
746 }
747
748 /// Access function to Hypre_method flag -- specified via enumeration.
749 unsigned& hypre_method()
750 {
751 return Hypre_method;
752 }
753
754 /// Access function to Internal_preconditioner flag -- specified
755 /// via enumeration.
757 {
759 }
760
761 /// Function to select use of 'simple' AMG smoothers as controlled
762 /// by AMG_simple_smoother flag
767
768 /// Access function to AMG_simple_smoother flag
770 {
771 return AMG_simple_smoother;
772 }
773
774 /// Function to select use of 'complex' AMG smoothers as controlled
775 /// by AMG_complex_smoother flag
780
781 /// Access function to AMG_complex_smoother flag
783 {
785 }
786
787 /// Access function to AMG_print_level
788 unsigned& amg_print_level()
789 {
790 return AMG_print_level;
791 }
792
793 /// Access function to AMG_smoother_iterations
795 {
797 }
798
799 /// Access function to AMG_coarsening flag
800 unsigned& amg_coarsening()
801 {
802 return AMG_coarsening;
803 }
804
805 /// Access function to AMG_max_levels
806 unsigned& amg_max_levels()
807 {
808 return AMG_max_levels;
809 }
810
811 /// Access function to AMG_damping parameter
812 double& amg_damping()
813 {
814 return AMG_damping;
815 }
816
817 /// Access function to AMG_strength
818 double& amg_strength()
819 {
820 return AMG_strength;
821 }
822
823 /// Access function to AMG_max_row_sum
825 {
826 return AMG_max_row_sum;
827 }
828
829 /// Access function to AMG_truncation
831 {
832 return AMG_truncation;
833 }
834
835 /// Access function to ParaSails symmetry flag
837 {
838 return ParaSails_symmetry;
839 }
840
841 /// Access function to ParaSails nlevel parameter
843 {
844 return ParaSails_nlevel;
845 }
846
847 /// Access function to ParaSails thresh parameter
849 {
850 return ParaSails_thresh;
851 }
852
853 /// Access function to ParaSails filter parameter
855 {
856 return ParaSails_filter;
857 }
858
859 /// Access function to Euclid drop tolerance parameter
861 {
862 return Euclid_droptol;
863 }
864
865 /// Access function to Euclid level parameter
866 /// for ILU(k) factorization
868 {
869 return Euclid_level;
870 }
871
872 /// Enable euclid rowScaling
874 {
875 Euclid_rowScale = true;
876 }
877
878 /// Disable euclid row scaling
880 {
881 Euclid_rowScale = false;
882 }
883
884 /// Enable use of Block Jacobi
885 /// as opposed to PILU
887 {
888 Euclid_using_BJ = true;
889 }
890
891 /// Disable use of Block Jacobi,
892 /// so PILU will be used
894 {
895 Euclid_using_BJ = false;
896 }
897
898 /// Function to switch on ILU(k) factorization for Euclid
899 /// (default is ILU(k) factorization)
901 {
902 Euclid_using_ILUT = false;
903 }
904
905 /// Function to switch on ILUT factorization for Euclid
906 /// (default is ILU(k) factorization)
908 {
909 Euclid_using_ILUT = true;
910 }
911
912 /// Function to set the level of printing from Euclid
913 /// when the Euclid destructor is called
914 /// 0: no printing (default)
915 /// 1: prints summary of runtime settings and timings
916 /// 2: as 1 plus prints memory usage
918 {
919 return Euclid_print_level;
920 }
921
922 /// Access function to Krylov_print_level
924 {
925 return Krylov_print_level;
926 }
927
928 /// Hypre copies matrix data from oomph-lib's CRDoubleMatrix
929 /// or DistributedCRDoubleMatrix into its own data structures,
930 /// doubling the memory requirements for the matrix.
931 /// As far as the Hypre solvers are concerned the oomph-lib matrix
932 /// is no longer required and could be deleted to save memory.
933 /// However, this will not be the expected behaviour and therefore
934 /// needs to be requested explicitly by the user by calling this function
935 /// which changes the
936 /// flag from false (its default) to true.
938 {
939 Delete_matrix = true;
940 }
941
942 /// Hypre copies matrix data from oomph-lib's CRDoubleMatrix
943 /// or DistributedCRDoubleMatrix into its own data structures,
944 /// doubling the memory requirements for the matrix.
945 /// Calling this function ensures that the matrix is not deleted.
947 {
948 Delete_matrix = false;
949 }
950
951
952 /// Function which uses problem_pt's get_jacobian(...) function to
953 /// generate a linear system which is then solved. This function deletes
954 /// any existing internal data and then generates a new Hypre solver.
955 void solve(Problem* const& problem_pt, DoubleVector& solution);
956
957 /// Function to solve the linear system defined by matrix_pt
958 /// and rhs. This function will delete any existing internal data
959 /// and generate a new Hypre solver.
960 /// \b Note: The matrix has to be of type CRDoubleMatrix or
961 /// Distributed CRDoubleMatrix.
962 /// \b Note: Hypre copies matrix data from oomph-lib's CRDoubleMatrix
963 /// or DistributedCRDoubleMatrix into its own data structures,
964 /// doubling the memory requirements for the matrix.
965 /// As far as the Hypre solvers are concerned, the oomph-lib matrix
966 /// is no longer required and could be deleted to save memory.
967 /// However, this will not be the expected behaviour and therefore
968 /// needs to be requested explicitly by the user by calling the
969 /// enable_delete_matrix() function.
970 void solve(DoubleMatrixBase* const& matrix_pt,
971 const DoubleVector& rhs,
973
974 /// Function to resolve a linear system using the existing solver
975 /// data, allowing a solve with a new right hand side vector. This
976 /// function must be used after a call to solve(...) with
977 /// enable_resolve set to true.
979
980 /// Function deletes all solver data.
981 void clean_up_memory();
982
983 private:
984 /// Hypre copies matrix data from oomph-lib's CRDoubleMatrix
985 /// or DistributedCRDoubleMatrix into its own data structures,
986 /// doubling the memory requirements for the matrix.
987 /// As far as the Hypre solvers are concerned the oomph-lib matrix
988 /// is no longer required and could be deleted to save memory.
989 /// However, this will not be the expected behaviour and therefore
990 /// needs to be requested explicitly by the user by changing this
991 /// flag from false (its default) to true.
993 };
994
995
996 ////////////////////////////////////////////////////////////////////////
997 ////////////////////////////////////////////////////////////////////////
998 ////////////////////////////////////////////////////////////////////////
999
1000 //=====================================================================
1001 /// An Preconditioner class using the suite of Hypre preconditioners
1002 /// to allow
1003 ///
1004 /// BoomerAMG (Algebraic Multi Grid),
1005 /// Euclid (ILU) or
1006 /// ParaSails (Approximate inverse)
1007 ///
1008 /// to be used for Preconditioner::preconditioner_solve(...).
1009 /// By default BoomerAMG is used.
1010 //====================================================================
1012 {
1013 public:
1014 /// Constructor. Provide optional string that is used
1015 /// in annotation of performance
1016 HyprePreconditioner(const std::string& context_string = "")
1017 {
1019
1020 // Hypre copies matrix data from oomph-lib's CRDoubleMatrix
1021 // or DistributedCRDoubleMatrix into its own data structures,
1022 // doubling the memory requirements for the matrix.
1023 // As far as the Hypre solvers are concerned the oomph-lib matrix
1024 // is no longer required and could be deleted to save memory.
1025 // However, this will not be the expected behaviour and therefore
1026 // needs to be requested explicitly by the user by changing this
1027 // flag from false (its default) to true.
1028 Delete_matrix = false;
1029
1030 // Do we want to output results of timings?
1031 Doc_time = true;
1032
1033 // General control paramaters for using HYPRE as preconditioner
1034 Tolerance = 0.0;
1035 Max_iter = 1;
1038
1039 // Initialise private double that accumulates the preconditioner
1040 // solve time of thi instantiation of this class. Is reported
1041 // in destructor if Report_my_cumulative_preconditioner_solve_time
1042 // is set to true
1045 }
1046
1047 /// Destructor.
1049 {
1051 {
1052 oomph_info << "~HyprePreconditioner() in context \" " << Context_string
1053 << "\": My_cumulative_preconditioner_solve_time = "
1055 }
1056 }
1057
1058 /// Broken copy constructor.
1060
1061 /// Broken assignment operator.
1062 void operator=(const HyprePreconditioner&) = delete;
1063
1064 /// Static double that accumulates the preconditioner
1065 /// solve time of all instantiations of this class. Reset
1066 /// it manually, e.g. after every Newton solve, using
1067 /// reset_cumulative_solve_times().
1069
1070 /// Static double that accumulates the preconditioner
1071 /// solve time of all instantiations of this class, labeled by
1072 /// context string. Reset
1073 /// it manually, e.g. after every Newton solve, using
1074 /// reset_cumulative_solve_times().
1075 static std::map<std::string, double> Context_based_cumulative_solve_time;
1076
1077
1078 /// Static unsigned that accumulates the number of preconditioner
1079 /// solves of all instantiations of this class. Reset
1080 /// it manually, e.g. after every Newton solve, using
1081 /// reset_cumulative_solve_times().
1083
1084 /// Static unsigned that accumulates the number of preconditioner
1085 /// solves of all instantiations of this class, labeled by
1086 /// context string. Reset
1087 /// it manually, e.g. after every Newton solve, using
1088 /// reset_cumulative_solve_times().
1089 static std::map<std::string, unsigned>
1091
1092 /// Static unsigned that stores nrow for the most recent
1093 /// instantiations of this class, labeled by
1094 /// context string. Reset
1095 /// it manually, e.g. after every Newton solve, using
1096 /// reset_cumulative_solve_times().
1097 static std::map<std::string, unsigned> Context_based_nrow;
1098
1099 /// Report cumulative solve times of all instantiations of this
1100 /// class
1101 static void report_cumulative_solve_times();
1102
1103 /// Reset cumulative solve times
1104 static void reset_cumulative_solve_times();
1105
1106 /// Enable reporting of cumulative solve time in destructor
1111
1112 /// Disable reporting of cumulative solve time in destructor
1117
1118 /// Enable documentation of preconditioner timings
1120 {
1121 Doc_time = true;
1122 }
1123
1124 /// Disable documentation of preconditioner timings
1126 {
1127 Doc_time = false;
1128 }
1129
1130 /// Access function to Hypre_method flag -- specified via enumeration.
1131 unsigned& hypre_method()
1132 {
1133 return Hypre_method;
1134 }
1135
1136 /// Access function to Internal_preconditioner flag -- specified
1137 /// via enumeration.
1139 {
1141 }
1142
1143 /// Function to select BoomerAMG as the preconditioner
1145 {
1147 }
1148
1149 /// Function to set the number of times to apply BoomerAMG
1151 {
1153 }
1154
1155 /// Return function for Max_iter
1156 unsigned& amg_iterations()
1157 {
1158 return Max_iter;
1159 }
1160
1161 /// Function to select use of 'simple' AMG smoothers as controlled
1162 /// by the flag AMG_simple_smoother
1167
1168 /// Access function to AMG_simple_smoother flag
1170 {
1171 return AMG_simple_smoother;
1172 }
1173
1174 /// Function to select use of 'complex' AMG smoothers as controlled
1175 /// by the flag AMG_complex_smoother
1180
1181 /// Access function to AMG_complex_smoother flag
1183 {
1184 return AMG_complex_smoother;
1185 }
1186
1187 /// Return function for the AMG_using_simple_smoothing_flag
1192
1193 /// Access function to AMG_print_level
1195 {
1196 return AMG_print_level;
1197 }
1198
1199 /// Access function to AMG_smoother_iterations
1201 {
1203 }
1204
1205 /// Access function to AMG_coarsening flag
1206 unsigned& amg_coarsening()
1207 {
1208 return AMG_coarsening;
1209 }
1210
1211 /// Access function to AMG_max_levels
1212 unsigned& amg_max_levels()
1213 {
1214 return AMG_max_levels;
1215 }
1216
1217 /// Access function to AMG_damping parameter
1218 double& amg_damping()
1219 {
1220 return AMG_damping;
1221 }
1222
1223 /// Access function to AMG_strength
1225 {
1226 return AMG_strength;
1227 }
1228
1229 /// Access function to AMG_max_row_sum
1231 {
1232 return AMG_max_row_sum;
1233 }
1234
1235 /// Access function to AMG_truncation
1237 {
1238 return AMG_truncation;
1239 }
1240
1241 /// Function to select ParaSails as the preconditioner
1243 {
1245 }
1246
1247 /// Access function to ParaSails symmetry flag
1249 {
1250 return ParaSails_symmetry;
1251 }
1252
1253 /// Access function to ParaSails nlevel parameter
1255 {
1256 return ParaSails_nlevel;
1257 }
1258
1259 /// Access function to ParaSails thresh parameter
1261 {
1262 return ParaSails_thresh;
1263 }
1264
1265 /// Access function to ParaSails filter parameter
1267 {
1268 return ParaSails_filter;
1269 }
1270
1271 /// Function to select use Euclid as the preconditioner
1273 {
1275 }
1276
1277 /// Access function to Euclid drop tolerance parameter
1279 {
1280 return Euclid_droptol;
1281 }
1282
1283 /// Access function to Euclid level parameter for ILU(k) factorization
1285 {
1286 return Euclid_level;
1287 }
1288
1289 /// Enable euclid rowScaling
1291 {
1292 Euclid_rowScale = true;
1293 }
1294
1295 /// Disable euclid row scaling
1297 {
1298 Euclid_rowScale = false;
1299 }
1300
1301 /// Enable use of Block Jacobi
1302 /// as opposed to PILU
1304 {
1305 Euclid_using_BJ = true;
1306 }
1307
1308 /// Disable use of Block Jacobi,
1309 /// so PILU will be used
1311 {
1312 Euclid_using_BJ = false;
1313 }
1314
1315 /// Function to switch on ILU(k) factorization for Euclid
1316 /// (default is ILU(k) factorization)
1318 {
1319 Euclid_using_ILUT = false;
1320 }
1321
1322 /// Function to switch on ILUT factorization for Euclid
1323 /// (default is ILU(k) factorization)
1325 {
1326 Euclid_using_ILUT = true;
1327 }
1328
1329 /// Function to set the level of printing from Euclid
1330 /// when the Euclid destructor is called
1331 /// 0: no printing (default)
1332 /// 1: prints summary of runtime settings and timings
1333 /// 2: as 1 plus prints memory usage
1335 {
1336 return Euclid_print_level;
1337 }
1338
1339 /// Hypre copies matrix data from oomph-lib's CRDoubleMatrix
1340 /// or DistributedCRDoubleMatrix into its own data structures,
1341 /// doubling the memory requirements for the matrix.
1342 /// As far as the Hypre solvers are concerned the oomph-lib matrix
1343 /// is no longer required and could be deleted to save memory.
1344 /// However, this will not be the expected behaviour and therefore
1345 /// needs to be requested explicitly by the user by calling this function
1346 /// which changes the
1347 /// flag from false (its default) to true.
1349 {
1350 Delete_matrix = true;
1351 }
1352
1353 /// Hypre copies matrix data from oomph-lib's CRDoubleMatrix
1354 /// or DistributedCRDoubleMatrix into its own data structures,
1355 /// doubling the memory requirements for the matrix.
1356 /// Calling this function ensures that the matrix is not deleted.
1358 {
1359 Delete_matrix = false;
1360 }
1361
1362 /// Function to set up a preconditioner for the linear
1363 /// system defined by matrix_pt. This function is required when
1364 /// preconditioning and must be called before using the
1365 /// preconditioner_solve(...) function. This interface allows
1366 /// HyprePreconditioner to be used as a Preconditioner object
1367 /// for oomph-lib's own IterativeLinearSolver class.
1368 /// \b Note: matrix_pt must point to an object of type
1369 /// CRDoubleMatrix or DistributedCRDoubleMatrix.
1370 /// \b Note: Hypre copies matrix data from oomph-lib's CRDoubleMatrix
1371 /// and DistributedCRDoubleMatrix into its own data structures,
1372 /// doubling the memory requirements for the matrix.
1373 /// As far as the Hypre solvers are concerned, the oomph-lib matrix
1374 /// is no longer required and could be deleted to save memory.
1375 /// However, this will not be the expected behaviour and therefore
1376 /// needs to be requested explicitly by the user by calling the
1377 /// enable_delete_matrix() function.
1378 void setup();
1379
1380 /// Function applies solver to vector r for preconditioning.
1381 /// This requires a call to setup(...) first.
1382 /// \b Note: Hypre copies matrix data from oomph-lib's CRDoubleMatrix
1383 /// or DistributedCRDoubleMatrix into its own data structures,
1384 /// doubling the memory requirements for the matrix.
1385 /// As far as the Hypre solvers are concerned, the oomph-lib matrix
1386 /// is no longer required and could be deleted to save memory.
1387 /// However, this will not be the expected behaviour and therefore
1388 /// needs to be requested explicitly by the user with the
1389 /// delete_matrix() function.
1391
1392 /// Function deletes all solver data.
1393 void clean_up_memory();
1394
1395 private:
1396 /// Hypre copies matrix data from oomph-lib's CRDoubleMatrix
1397 /// or DistributedCRDoubleMatrix into its own data structures,
1398 /// doubling the memory requirements for the matrix.
1399 /// As far as the Hypre solvers are concerned the oomph-lib matrix
1400 /// is no longer required and could be deleted to save memory.
1401 /// However, this will not be the expected behaviour and therefore
1402 /// needs to be requested explicitly by the user by changing this
1403 /// flag from false (its default) to true.
1405
1406 // Flag is true to output results of timings
1408
1409 /// Private double that accumulates the preconditioner
1410 /// solve time of thi instantiation of this class. Is reported
1411 /// in destructor if Report_my_cumulative_preconditioner_solve_time
1412 /// is set to true
1414
1415 /// Bool to request report of My_cumulative_preconditioner_solve_time
1416 /// in destructor
1418
1419 /// String can be used to provide context for annotation
1420 std::string Context_string;
1421 };
1422
1423
1424 ////////////////////////////////////////////////////////////////////////
1425 ////////////////////////////////////////////////////////////////////////
1426 ////////////////////////////////////////////////////////////////////////
1427
1428 //==================================================================
1429 /// Default settings for various uses of the Hypre solver
1430 //==================================================================
1431 namespace Hypre_default_settings
1432 {
1433 /// Set default parameters for use as preconditioner in
1434 /// for momentum block in Navier-Stokes problem
1437
1438 /// Set default parameters for use as preconditioner in
1439 /// 2D Poisson-type problem.
1442
1443 /// Set default parameters for use as preconditioner in
1444 /// 3D Poisson-type problem.
1447
1448 } // namespace Hypre_default_settings
1449} // namespace oomph
1450
1451#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....
An interface class to the suite of Hypre solvers and preconditioners to allow use of:
bool Output_info
Flag is true to output info and results of timings.
bool Delete_input_data
Internal flag which is true when hypre_setup or hypre_solve can delete input matrix.
unsigned AMG_coarsening
AMG coarsening strategy. Coarsening types include: 0 = CLJP (parallel coarsening using independent se...
int Euclid_level
Euclid level parameter for ILU(k) factorization.
bool Euclid_using_BJ
Flag to determine if Block Jacobi is used instead of PILU.
bool Euclid_rowScale
Flag to switch on Euclid row scaling.
bool Returning_distributed_solution
Internal flag which tell the solver if the solution Vector to be returned is distributed or not.
double AMG_strength
Connection strength threshold parameter for BoomerAMG.
HypreInterface(const HypreInterface &)=delete
Broken copy constructor.
void enable_hypre_error_messages()
Turn on the hypre error messages.
LinearAlgebraDistribution * Hypre_distribution_pt
the distribution for this helpers-
int ParaSails_nlevel
ParaSails nlevel parameter.
double ParaSails_filter
ParaSails filter parameter.
unsigned AMGEuclidSmoother_level
HYPRE_Solver Solver
The Hypre solver used in solve(...), resolve(...) or preconditioner_solve(...). [This is a C structur...
int ParaSails_symmetry
ParaSails symmetry flag, used to inform ParaSails of Symmetry of definitenss of problem and type of P...
bool Using_distributed_rhs
Internal flag which tell the solver if the rhs Vector is distributed or not.
HypreInterface()
Constructor.
unsigned existing_preconditioner()
Function return value of which preconditioner (if any) is stored.
HYPRE_ParCSRMatrix Matrix_par
The Hypre_ParCSRMatrix version of the matrix used in solve(...), resolve(...) or preconditioner_solve...
void operator=(const HypreInterface &)=delete
Broken assignment operator.
void disable_hypre_error_messages()
Turn off hypre error messages.
unsigned Internal_preconditioner
Preconditioner method flag used with Hypre's PCG, GMRES and BiCGStab in solve(...) or resolve(....
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!
~HypreInterface()
Destructor.
HYPRE_Solver Preconditioner
The internal Hypre preconditioner used in conjunction with Solver. [This is a C structure!...
unsigned AMGEuclidSmoother_print_level
bool AMGEuclidSmoother_use_block_jacobi
unsigned AMG_complex_smoother
Complex smoothing methods used in BoomerAMG. To use these methods set AMG_using_simple_smoothing to f...
bool Delete_matrix
Hypre copies matrix data from oomph-lib's CRDoubleMatrix into its own data structures,...
double ParaSails_thresh
ParaSails thresh parameter.
double AMG_truncation
Interpolation truncation factor for BoomerAMG.
unsigned Krylov_print_level
Used to set the Hypre printing level for the Krylov subspace solvers.
unsigned AMG_smoother_iterations
The number of smoother iterations to apply.
unsigned existing_solver()
Function to return value of which solver (if any) is currently stored.
void hypre_clean_up_memory()
Function deletes all solver data.
void doc_hypre_parameters()
Doc parameters used in hypre solver.
bool Hypre_error_messages
Flag to determine if non-zero values of the Hypre error flag plus Hypre error messages are output to ...
double Tolerance
Tolerance used to terminate solver.
bool AMG_using_simple_smoothing
Flag to determine whether simple smoothers (determined by the AMG_simple_smoother flag) or complex sm...
unsigned Existing_solver
Used to keep track of which solver (if any) is currently stored.
unsigned Max_iter
Maximum number of iterations used in solver.
double AMG_damping
Damping factor for BoomerAMG smoothed Jacobi or hybrid SOR.
unsigned AMG_simple_smoother
Simple smoothing methods used in BoomerAMG. To use these methods set AMG_using_simple_smoothing to tr...
unsigned AMG_print_level
Used to set the Hypre printing level for AMG 0: no printout 1: print setup information 2: print solve...
HYPRE_IJMatrix Matrix_ij
The Hypre_IJMatrix version of the matrix used in solve(...), resolve(...) or preconditioner_solve(....
unsigned Euclid_print_level
Flag to set the level of printing from Euclid when the Euclid destructor is called 0: no printing (de...
unsigned AMG_max_levels
Maximum number of levels used in AMG.
unsigned Hypre_method
Hypre method flag. Valid values are specified in enumeration.
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.
double Euclid_droptol
Euclid drop tolerance for ILU(k) and ILUT factorization.
unsigned Existing_preconditioner
Used to keep track of which preconditioner (if any) is currently stored.
bool Euclid_using_ILUT
Flag to determine if ILUT (if true) or ILU(k) is used in Euclid.
void hypre_solver_setup()
Sets up the data required for to use as an oomph-lib LinearSolver or Preconditioner....
An Preconditioner class using the suite of Hypre preconditioners to allow.
bool & amg_using_simple_smoothing_flag()
Return function for the AMG_using_simple_smoothing_flag.
unsigned & amg_simple_smoother()
Access function to AMG_simple_smoother flag.
static void reset_cumulative_solve_times()
Reset cumulative solve times.
unsigned & amg_print_level()
Access function to AMG_print_level.
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...
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.
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.
unsigned & internal_preconditioner()
Access function to Internal_preconditioner flag – specified via enumeration.
void set_amg_iterations(const unsigned &amg_iterations)
Function to set the number of times to apply BoomerAMG.
void clean_up_memory()
Function deletes all solver data.
unsigned & amg_max_levels()
Access function to AMG_max_levels.
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....
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.
double & parasails_filter()
Access function to ParaSails filter parameter.
void euclid_using_ILUK()
Function to switch on ILU(k) factorization for Euclid (default is 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.
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,...
void enable_euclid_using_BJ()
Enable use of Block Jacobi as opposed to PILU.
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.
double & amg_max_row_sum()
Access function to AMG_max_row_sum.
void disable_doc_time()
Disable documentation of preconditioner timings.
static std::map< std::string, unsigned > Context_based_nrow
Static unsigned that stores nrow for the most recent instantiations of this class,...
void amg_using_simple_smoothing()
Function to select use of 'simple' AMG smoothers as controlled by the flag AMG_simple_smoother.
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 & amg_damping()
Access function to AMG_damping parameter.
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.
unsigned & amg_coarsening()
Access function to AMG_coarsening flag.
unsigned & amg_smoother_iterations()
Access function to AMG_smoother_iterations.
double & parasails_thresh()
Access function to ParaSails thresh parameter.
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.
double & amg_strength()
Access function to AMG_strength.
double & amg_truncation()
Access function to AMG_truncation.
unsigned & euclid_print_level()
Function to set the level of printing from Euclid when the Euclid destructor is called 0: no printing...
unsigned & amg_complex_smoother()
Access function to AMG_complex_smoother flag.
static unsigned Cumulative_npreconditioner_solve
Static unsigned that accumulates the number of preconditioner solves of all instantiations of this cl...
void preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Function applies solver to vector r for preconditioning. This requires a call to setup(....
unsigned & hypre_method()
Access function to Hypre_method flag – specified via enumeration.
unsigned & amg_iterations()
Return function for Max_iter.
int & euclid_level()
Access function to Euclid level parameter for ILU(k) factorization.
int & parasails_nlevel()
Access function to ParaSails nlevel parameter.
HyprePreconditioner(const HyprePreconditioner &)=delete
Broken copy constructor.
double & euclid_droptol()
Access function to Euclid drop tolerance parameter.
HyprePreconditioner(const std::string &context_string="")
Constructor. Provide optional string that is used in annotation of performance.
int & parasails_symmetry()
Access function to ParaSails symmetry flag.
void enable_report_my_cumulative_preconditioner_solve_time()
Enable reporting of cumulative solve time in destructor.
An LinearSolver class using the suite of Hypre solvers to allow.
void disable_euclid_using_BJ()
Disable use of Block Jacobi, so PILU will be used.
void disable_euclid_rowScale()
Disable euclid row scaling.
unsigned & krylov_print_level()
Access function to Krylov_print_level.
void enable_euclid_using_BJ()
Enable use of Block Jacobi as opposed to PILU.
unsigned & hypre_method()
Access function to Hypre_method flag – specified via enumeration.
int & parasails_nlevel()
Access function to ParaSails nlevel parameter.
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...
HypreSolver(const HypreSolver &)=delete
Broken copy constructor.
~HypreSolver()
Empty destructor.
unsigned & amg_max_levels()
Access function to AMG_max_levels.
int & euclid_level()
Access function to Euclid level parameter for ILU(k) factorization.
double & euclid_droptol()
Access function to Euclid drop tolerance parameter.
void clean_up_memory()
Function deletes all solver data.
unsigned & amg_simple_smoother()
Access function to AMG_simple_smoother flag.
unsigned & internal_preconditioner()
Access function to Internal_preconditioner flag – specified via enumeration.
double & amg_truncation()
Access function to AMG_truncation.
double & amg_max_row_sum()
Access function to AMG_max_row_sum.
HypreSolver()
Constructor.
unsigned & max_iter()
Access function to Max_iter.
int & parasails_symmetry()
Access function to ParaSails symmetry flag.
bool Delete_matrix
Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own da...
void amg_using_complex_smoothing()
Function to select use of 'complex' AMG smoothers as controlled by AMG_complex_smoother flag.
void euclid_using_ILUT()
Function to switch on ILUT factorization for Euclid (default is ILU(k) factorization)
unsigned & amg_print_level()
Access function to AMG_print_level.
void amg_using_simple_smoothing()
Function to select use of 'simple' AMG smoothers as controlled by AMG_simple_smoother flag.
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.
unsigned & amg_complex_smoother()
Access function to AMG_complex_smoother flag.
unsigned & amg_coarsening()
Access function to AMG_coarsening flag.
unsigned & euclid_print_level()
Function to set the level of printing from Euclid when the Euclid destructor is called 0: no printing...
double & parasails_filter()
Access function to ParaSails filter parameter.
double & amg_strength()
Access function to AMG_strength.
void euclid_using_ILUK()
Function to switch on ILU(k) factorization for Euclid (default is ILU(k) factorization)
void operator=(const HypreSolver &)=delete
Broken assignment operator.
double & amg_damping()
Access function to AMG_damping parameter.
double & tolerance()
Access function to Tolerance.
unsigned & amg_smoother_iterations()
Access function to AMG_smoother_iterations.
void disable_resolve()
Disable resolve function (overloads the LinearSolver disable_resolve function).
void enable_delete_matrix()
Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own da...
void disable_delete_matrix()
Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own da...
double & parasails_thresh()
Access function to ParaSails thresh parameter.
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,...
bool Doc_time
Boolean flag that indicates whether the time taken.
bool Enable_resolve
Boolean that indicates whether the matrix (or its factors, in the case of direct solver) should be st...
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:154
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
hypre_Error hypre__global_error
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 Number_of_active_hypre_solvers
Number of active Hypre solvers (smart pointer like behaviuor to make sure that the initialise/finaliz...
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.
void set_defaults_for_3D_poisson_problem(HyprePreconditioner *hypre_preconditioner_pt)
Set default parameters for use as preconditioner in 3D Poisson-type problem.
void set_defaults_for_2D_poisson_problem(HyprePreconditioner *hypre_preconditioner_pt)
Set default parameters for use as preconditioner in 2D Poisson-type problem.
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
OomphInfo oomph_info
Single (global) instantiation of the OomphInfo object – this is used throughout the library as a "rep...