hypre_solver.cc
Go to the documentation of this file.
1// LIC// ====================================================================
2// LIC// This file forms part of oomph-lib, the object-oriented,
3// LIC// multi-physics finite-element library, available
4// LIC// at http://www.oomph-lib.org.
5// LIC//
6// LIC// 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#include "hypre_solver.h"
27
28// For problem->get_jacobian(...)
29#include "problem.h"
30
31
32namespace oomph
33{
34 ////////////////////////////////////////////////////////////////////////
35 ////////////////////////////////////////////////////////////////////////
36 ////////////////////////////////////////////////////////////////////////
37
38 //==================================================================
39 /// Default settings for various uses of the HYPRE solver
40 //==================================================================
41 namespace Hypre_default_settings
42 {
43 /// Set default parameters for use as preconditioner in
44 /// 2D Poisson-type problem.
47 {
48 // Set iterations to 1
49 hypre_preconditioner_pt->set_amg_iterations(1);
50
51 // Use simple smoother
52 hypre_preconditioner_pt->amg_using_simple_smoothing();
53
54 // Smoother types:
55 // 0=Jacobi
56 // 1=Gauss-Seidel
57 hypre_preconditioner_pt->amg_simple_smoother() = 1;
58
59 // AMG preconditioner
61
62 // Choose strength parameter for amg
63 hypre_preconditioner_pt->amg_strength() = 0.25;
64
65 // Coarsening type
66 hypre_preconditioner_pt->amg_coarsening() = 0;
67 }
68
69
70 /// Set default parameters for use as preconditioner in
71 /// 3D Poisson-type problem.
74 {
75 // Set default settings as for 2D Poisson
77
78 // Change strength parameter for amg
79 hypre_preconditioner_pt->amg_strength() = 0.7;
80 }
81
82
83 /// Set default parameters for use as preconditioner in
84 /// for momentum block in Navier-Stokes problem
87 {
88 // Set default settings as for 2D Poisson
90
91 // Change smoother type:
92 // 0=Jacobi
93 // 1=Gauss-Seidel
94 hypre_preconditioner_pt->amg_simple_smoother() = 0;
95
96 // Set smoother damping
97 hypre_preconditioner_pt->amg_damping() = 0.5;
98
99 // Change strength parameter for amg
100 hypre_preconditioner_pt->amg_strength() = 0.75;
101 }
102
103 } // namespace Hypre_default_settings
104
105
106 ////////////////////////////////////////////////////////////////////////
107 ////////////////////////////////////////////////////////////////////////
108 // functions for HypreHelpers namespace
109 ////////////////////////////////////////////////////////////////////////
110 ////////////////////////////////////////////////////////////////////////
111
112 namespace HypreHelpers
113 {
114
115 //========================================================================
116 /// Number of active Hypre solvers (smart pointer like behaviuor
117 /// to make sure that the initialise/finalize functions are only
118 /// called the required number of times
119 //========================================================================
121
122 //========================================================================
123 /// Default for AMG strength (0.25 recommended for 2D problems;
124 /// larger (0.5-0.75, say) for 3D
125 //========================================================================
126 double AMG_strength = 0.25;
127
128 //========================================================================
129 /// Default AMG coarsening strategy. Coarsening types include:
130 /// 0 = CLJP (parallel coarsening using independent sets)
131 /// 1 = classical RS with no boundary treatment (not recommended
132 /// in parallel)
133 /// 3 = modified RS with 3rd pass to add C points on the boundaries
134 /// 6 = Falgout (uses 1 then CLJP using interior coarse points as
135 /// first independent set)
136 /// 8 = PMIS (parallel coarsening using independent sets - lower
137 /// complexities than 0, maybe also slower convergence)
138 /// 10= HMIS (one pass RS on each processor then PMIS on interior
139 /// coarse points as first independent set)
140 /// 11= One pass RS on each processor (not recommended)
141 //========================================================================
142 unsigned AMG_coarsening = 6;
143
144
145 //========================================================================
146 /// AMG interpolation truncation factor
147 //========================================================================
148 double AMG_truncation = 0.0;
149
150 //========================================================================
151 /// Helper function to check the Hypre error flag, return the message
152 /// associated with any error, and reset the global error flag to zero.
153 /// This function also returns the error value.
154 //========================================================================
155 int check_HYPRE_error_flag(std::ostringstream& message)
156 {
157 // get the Hypre error flag
158 int err = HYPRE_GetError();
159
160 // Tell us all about it...
161 if (err)
162 {
163 oomph_info << "Hypre error flag=" << err << std::endl;
164 char* error_message = new char[128];
165 HYPRE_DescribeError(err, error_message);
166 message << "WARNING: " << std::endl
167 << "HYPRE error message: " << error_message << std::endl;
168 delete[] error_message;
169 }
170
171 // reset Hypre's global error flag
172 hypre__global_error.error_flag = 0;
173
174 return err;
175 }
176
177
178 //========================================================================
179 /// Helper function to create a HYPRE_IJVector and HYPRE_ParVector.
180 /// + If no MPI then serial vectors are created
181 /// + If MPI and serial input vector then distributed hypre vectors are
182 /// created
183 /// + If MPI and distributed input vector the distributed output vectors
184 /// are created.
185 //========================================================================
190 {
191 // the lower and upper row of the vector on this processor
192 unsigned lower = dist_pt->first_row();
193 unsigned upper = dist_pt->first_row() + dist_pt->nrow_local() - 1;
194
195 // number of local rows
196 unsigned nrow_local = dist_pt->nrow_local();
197
198 // initialize Hypre vector
199#ifdef OOMPH_HAS_MPI
201 oomph_vec.distribution_pt()->communicator_pt()->mpi_comm(),
202 lower,
203 upper,
205#else
207#endif
210
211 // set up array containing indices
212 int* indices = new int[nrow_local];
213 double* values = new double[nrow_local];
214 const unsigned hypre_first_row = dist_pt->first_row();
215 unsigned j = 0;
216 if (!oomph_vec.distributed() && dist_pt->distributed())
217 {
219 }
220 const double* o_pt = oomph_vec.values_pt();
221 for (unsigned i = 0; i < nrow_local; i++)
222 {
224 values[i] = o_pt[j + i];
225 }
226
227 // insert values
228 HYPRE_IJVectorSetValues(hypre_ij_vector, nrow_local, indices, values);
229
230 // assemble vectors
233
234 // clean up
235 delete[] indices;
236 delete[] values;
237 }
238
239
240 //========================================================================
241 /// Helper function to create a HYPRE_IJVector and HYPRE_ParVector.
242 /// + If no MPI then serial vectors are created
243 /// + If MPI and serial input vector then distributed hypre vectors are
244 /// created
245 /// + If MPI and distributed input vector the distributed output vectors
246 /// are created.
247 //========================================================================
251 {
252 // the lower and upper row of the vector on this processor
253 unsigned lower = dist_pt->first_row();
254 unsigned upper = dist_pt->first_row() + dist_pt->nrow_local() - 1;
255
256 // initialize Hypre vector
257#ifdef OOMPH_HAS_MPI
259 dist_pt->communicator_pt()->mpi_comm(), lower, upper, &hypre_ij_vector);
260#else
262#endif
265
266 // assemble vectors
269 }
270
271
272 //========================================================================
273 /// Helper function to create a serial HYPRE_IJMatrix and
274 /// HYPRE_ParCSRMatrix from a CRDoubleMatrix
275 /// NOTE: dist_pt is rebuilt to match the distribution of the hypre solver
276 /// which is not necassarily the same as the oomph lib matrix
277 //========================================================================
282 {
283#ifdef PARANOID
284 // check that the matrix is built
285 if (!oomph_matrix->built())
286 {
287 std::ostringstream error_message;
288 error_message << "The matrix has not been built";
289 throw OomphLibError(error_message.str(),
292 }
293 // check the matrix is square
294 if (oomph_matrix->nrow() != oomph_matrix->ncol())
295 {
296 std::ostringstream error_message;
297 error_message << "create_HYPRE_Matrix require a square matrix. "
298 << "Matrix is " << oomph_matrix->nrow() << " by "
299 << oomph_matrix->ncol() << std::endl;
300 throw OomphLibError(error_message.str(),
303 }
304#endif
305
306#ifdef OOMPH_HAS_MPI
307 // Trap the case when we have compiled with MPI,
308 // but are running in serial
310 {
311 std::ostringstream error_stream;
313 << "Oomph-lib has been compiled with MPI support and "
314 << "you are using HYPRE.\n"
315 << "For this combination of flags, MPI must be initialised.\n"
316 << "Call MPI_Helpers::init() in the "
317 << "main() function of your driver code\n";
318 throw OomphLibError(
320 }
321#endif
322
323 // find number of rows/columns
324 const unsigned nrow = int(oomph_matrix->nrow());
325
326 // get pointers to the matrix
327
328 // column indices of matrix
329 const int* matrix_cols = oomph_matrix->column_index();
330
331 // entries of matrix
332 const double* matrix_vals = oomph_matrix->value();
333
334 // row starts
335 const int* matrix_row_start = oomph_matrix->row_start();
336
337 // build the distribution
338 if (oomph_matrix->distribution_pt()->distributed())
339 {
340 dist_pt->build(oomph_matrix->distribution_pt());
341 }
342 else
343 {
344 bool distributed = true;
345 if (oomph_matrix->distribution_pt()->communicator_pt()->nproc() == 1)
346 {
347 distributed = false;
348 }
349 dist_pt->build(oomph_matrix->distribution_pt()->communicator_pt(),
350 nrow,
351 distributed);
352 }
353
354 // initialize hypre matrix
355 unsigned lower = dist_pt->first_row();
356 unsigned upper = lower + dist_pt->nrow_local() - 1;
357
358#ifdef OOMPH_HAS_MPI
359 int i_err = 0;
360 i_err = HYPRE_IJMatrixCreate(dist_pt->communicator_pt()->mpi_comm(),
361 lower,
362 upper,
363 lower,
364 upper,
366 if (i_err)
367 {
368 throw OomphLibError("Error in HYPRE_IJMatrixCreate()",
371 }
372#else
373 int i_err = 0;
376 if (i_err)
377 {
378 throw OomphLibError("Error in HYPRE_IJMatrixCreate()",
381 }
382#endif
383
385 if (i_err)
386 {
387 throw OomphLibError("Error in HYPRE_IJMatrixSetObjectType()",
390 }
391
393 if (i_err)
394 {
395 throw OomphLibError("Error in HYPRE_IJMatrixInitialize()",
398 }
399 // set up a row map
400 // and first row / nrow_local
401 const unsigned hypre_nrow_local = dist_pt->nrow_local();
402 const unsigned hypre_first_row = dist_pt->first_row();
403 int* ncols_per_row = new int[hypre_nrow_local];
404 int* row_map = new int[hypre_nrow_local];
405 for (unsigned i = 0; i < hypre_nrow_local; i++)
406 {
407 unsigned j = i;
408 if (!oomph_matrix->distributed() && dist_pt->distributed())
409 {
411 }
414 }
415
416 // put values in HYPRE matrix
417 int local_start = 0;
418 if (!oomph_matrix->distributed() && dist_pt->distributed())
419 {
421 }
422
423
427 row_map,
430
431 // assemble matrix
434
435 // tidy up memory
436 delete[] ncols_per_row;
437 delete[] row_map;
438 }
439
440 //====================================================================
441 /// Helper function to set Euclid options using a command line
442 /// like array.
443 //=====================================================================
445 const bool& use_row_scaling,
446 const bool& use_ilut,
447 const int& level,
448 const double& drop_tol,
449 const int& print_level,
451 {
452 // Easier to use C-arrays rather than std::strings because Euclid takes
453 // char** as argument and because C++ doesn't provide decent number to
454 // std::string conversion functions.
455
456 int n_args = 0;
457 const char* args[22];
458
459 // first argument is empty string
460 args[n_args++] = "";
461
462 // switch on/off Block Jacobi ILU
463 args[n_args++] = "-bj";
465 {
466 args[n_args++] = "1";
467 }
468 else
469 {
470 args[n_args++] = "0";
471 }
472
473 // switch on/off row scaling
474 args[n_args++] = "-rowScale";
475 if (use_row_scaling)
476 {
477 args[n_args++] = "1";
478 }
479 else
480 {
481 args[n_args++] = "0";
482 }
483
484 // set level for ILU(k)
485 args[n_args++] = "-level";
486 char level_value[10];
487 snprintf(level_value, sizeof(level_value), "%d", level);
489
490 // // set drop tol for ILU(k) factorization
491 // args[n_args++] = "-sparseA";
492 // char droptol[20];
493 // snprintf(droptol, sizeof(droptol), "%f",drop_tol);
494 // args[n_args++] = droptol;
495
496 // // set ILUT factorization if required
497 // if (use_ilut)
498 // {
499 // args[n_args++] = "-ilut";
500 // args[n_args++] = droptol;
501 // }
502
503 // set printing of Euclid data
504 if (print_level == 0)
505 {
506 args[n_args++] = "-eu_stats";
507 args[n_args++] = "0";
508 args[n_args++] = "-eu_mem";
509 args[n_args++] = "0";
510 }
511 if (print_level == 1)
512 {
513 args[n_args++] = "-eu_stats";
514 args[n_args++] = "1";
515 args[n_args++] = "-eu_mem";
516 args[n_args++] = "0";
517 }
518 if (print_level == 2)
519 {
520 args[n_args++] = "-eu_stats";
521 args[n_args++] = "1";
522 args[n_args++] = "-eu_mem";
523 args[n_args++] = "1";
524 }
525
526 // set next entry in array to null
527 args[n_args] = 0;
528
529 // Send the parameters
530 HYPRE_EuclidSetParams(euclid_object, n_args, const_cast<char**>(args));
531 }
532
533 } // namespace HypreHelpers
534
535
536 ///////////////////////////////////////////////////////////////////////
537 ///////////////////////////////////////////////////////////////////////
538 // functions for HypreInterface class
539 ///////////////////////////////////////////////////////////////////////
540 ///////////////////////////////////////////////////////////////////////
541
542
543 //=============================================================================
544 /// Helper function which creates a Hypre matrix from a CRDoubleMatrix
545 /// If OOMPH-LIB has been set up for MPI use, the Hypre matrix is
546 /// distributed over the available processors.
547 //=============================================================================
549 {
550 // reset Hypre's global error flag
551 hypre__global_error.error_flag = 0;
552
553 // issue warning if the matrix is small compared to the number of processors
554 if (unsigned(2 * matrix_pt->distribution_pt()->communicator_pt()->nproc()) >
555 matrix_pt->nrow())
556 {
558 << "Warning: HYPRE based solvers may fail if 2*number of processors "
559 << "is greater than the number of unknowns!" << std::endl;
560 }
561
562 // store the distribution
563 // generate the Hypre matrix
566
567 // Output error messages if required
569 {
570 std::ostringstream message;
572 if (err)
573 {
574 throw OomphLibError(
576 }
577 }
578
579 // delete CRDoubleMatrix if required
581 {
582 matrix_pt->clear();
583 }
584 }
585
586
587 //=============================================================================
588 /// Sets up the solver data required for use in an oomph-lib
589 /// LinearSolver or Preconditioner, once the Hypre matrix has been
590 /// generated using hypre_matrix_setup(...).
591 //=============================================================================
593 {
594 // Store time
595 double t_start = TimingHelpers::timer();
596 double t_end = 0;
597
598
599 // reset Hypre's global error flag
600 hypre__global_error.error_flag = 0;
601
602 // create dummy Hypre vectors which are required for setup
611
612 // Set up internal preconditioner for CG, GMRES or BiCGSTAB
613 // --------------------------------------------------------
614 if ((Hypre_method >= CG) && (Hypre_method <= BiCGStab))
615 {
616 // AMG preconditioner
618 {
619 // set up BoomerAMG
630
632 {
634
635 // New parameter. From HYPRE_parcsr_ls.h:
636 // (Optional) Defines in which order the points are relaxed. There
637 // are the following options for \e relax_order:
638 //
639 // - 0 : the points are relaxed in natural or lexicographic order
640 // on each processor
641 // - 1 : CF-relaxation is used, i.e on the fine grid and the down
642 // cycle the
643 // coarse points are relaxed first, followed by the fine
644 // points; on the up cycle the F-points are relaxed first,
645 // followed by the C-points. On the coarsest level, if an
646 // iterative scheme is used, the points are relaxed in
647 // lexicographic order.
648 //
649 // The default is 0.
650 // NOTE: Using option causes non-convergence in parallel mode.
652
653 // Number of smoother iterations
655
656
657 // Relax weight/damping for all levels
659 }
660 else
661 {
666
667 // If we are using Euclid then set up additional Euclid only options
668 if (AMG_complex_smoother == 9)
669 {
678 }
679 }
680
682 }
683
684 // Euclid preconditioner
686 {
687#ifdef OOMPH_HAS_MPI
690#else
692#endif
693
694 // Set parameters
702
704 }
705
706 // ParaSails preconditioner
708 {
709#ifdef OOMPH_HAS_MPI
713#else
715#endif
721 }
722
723 // check error flag
725 {
726 std::ostringstream message;
728 if (err)
729 {
731 "HypreSolver::hypre_setup()",
733 }
734 }
735 } // end of setting up internal preconditioner
736
737
738 // set up solver
739 // -------------
741
742 // AMG solver
743 if (Hypre_method == BoomerAMG)
744 {
745 if (Output_info)
746 {
747 oomph_info << "Setting up BoomerAMG, ";
748 }
749
750 // set up BoomerAMG
761
763 {
765
766 // New parameter. From HYPRE_parcsr_ls.h:
767 // (Optional) Defines in which order the points are relaxed. There are
768 // the following options for \e relax_order:
769 //
770 // - 0 : the points are relaxed in natural or lexicographic order on
771 // each processor
772 // - 1 : CF-relaxation is used, i.e on the fine grid and the down
773 // cycle the
774 // coarse points are relaxed first, followed by the fine
775 // points; on the up cycle the F-points are relaxed first,
776 // followed by the C-points. On the coarsest level, if an
777 // iterative scheme is used, the points are relaxed in
778 // lexicographic order.
779 //
780 // The default is 0.
781 // NOTE: Using option causes non-convergence in parallel mode.
783
784 // Number of smoother iterations
786
787
788 // Relax weight/damping for all levels
790 }
791 else
792 {
796
797 /* Other settings
798 * 6 & Schwarz smoothers & HYPRE_BoomerAMGSetDomainType,
799 * HYPRE_BoomerAMGSetOverlap, \\
800 * & & HYPRE_BoomerAMGSetVariant, HYPRE_BoomerAMGSetSchwarzRlxWeight
801 * \\
802 * 7 & Pilut & HYPRE_BoomerAMGSetDropTol, HYPRE_BoomerAMGSetMaxNzPerRow
803 * \\
804 * 8 & ParaSails & HYPRE_BoomerAMGSetSym, HYPRE_BoomerAMGSetLevel, \\
805 * & & HYPRE_BoomerAMGSetFilter, HYPRE_BoomerAMGSetThreshold \\
806 * 9 & Euclid & HYPRE_BoomerAMGSetEuclidFile \\
807 */
808
809 // If we are using Euclid then set up additional Euclid only options
810 if (AMG_complex_smoother == 9)
811 {
820 }
821
822 // Add any others here as required...
823 }
824
825 // MemoryUsage::doc_memory_usage("before amg setup [solver]");
826 // MemoryUsage::insert_comment_to_continous_top("BEFORE AMG SETUP
827 // [SOLVER]");
828
830
831 // MemoryUsage::doc_memory_usage("after amg setup [solver]");
832 // MemoryUsage::insert_comment_to_continous_top("AFTER AMG SETUP
833 // [SOLVER]");
834
836 }
837
838 // Euclid solver
839 else if (Hypre_method == Euclid)
840 {
841 if (Output_info)
842 {
843 oomph_info << "Setting up Euclid, ";
844 }
845#ifdef OOMPH_HAS_MPI
847 &Solver);
848#else
850#endif
851
852 // Set parameters
860
863 }
864
865 // ParaSails preconditioner
866 else if (Hypre_method == ParaSails)
867 {
868 if (Output_info)
869 {
870 oomph_info << "Setting up ParaSails, ";
871 }
872#ifdef OOMPH_HAS_MPI
875#else
877#endif
881
884 }
885
886 // CG solver
887 else if (Hypre_method == CG)
888 {
889 if (Output_info)
890 {
891 oomph_info << "Setting up CG, ";
892 }
893
894#ifdef OOMPH_HAS_MPI
897#else
899#endif
904
905 // set preconditioner
907 {
908 if (Output_info)
909 {
910 oomph_info << " with BoomerAMG preconditioner, ";
911 }
912
917 }
918 else if (Internal_preconditioner == Euclid) // Euclid
919 {
920 if (Output_info)
921 {
922 oomph_info << " with Euclid ILU preconditioner, ";
923 }
924
929 }
930 else if (Internal_preconditioner == ParaSails) // ParaSails
931 {
932 if (Output_info)
933 {
934 oomph_info << " with ParaSails approximate inverse preconditioner, ";
935 }
936
941 }
942 else
943 {
944 if (Output_info)
945 {
946 oomph_info << " with no preconditioner";
947 }
948 }
949
950
955
956
958 }
959
960 // GMRES solver
961 else if (Hypre_method == GMRES)
962 {
963 if (Output_info)
964 {
965 oomph_info << "Setting up GMRES";
966 }
967
968#ifdef OOMPH_HAS_MPI
971#else
973#endif
979
980 // set preconditioner
982 {
983 if (Output_info)
984 {
985 oomph_info << " with BoomerAMG preconditioner, ";
986 }
987
992 }
993 else if (Internal_preconditioner == Euclid) // Euclid
994 {
995 if (Output_info)
996 {
997 oomph_info << " with Euclid ILU preconditioner, ";
998 }
999
1004 }
1005 else if (Internal_preconditioner == ParaSails) // ParaSails
1006 {
1007 if (Output_info)
1008 {
1009 oomph_info << " with ParaSails approximate inverse preconditioner, ";
1010 }
1011
1016 }
1017 else
1018 {
1019 if (Output_info)
1020 {
1021 oomph_info << " with no preconditioner";
1022 }
1023 }
1024
1029
1031 }
1032
1033 // BiCGStab solver
1034 else if (Hypre_method == BiCGStab)
1035 {
1036 if (Output_info)
1037 {
1038 oomph_info << "Setting up BiCGStab";
1039 }
1040#ifdef OOMPH_HAS_MPI
1043#else
1045#endif
1050
1051 // set preconditioner
1052 if (Internal_preconditioner == BoomerAMG) // AMG
1053 {
1054 if (Output_info)
1055 {
1056 oomph_info << " with BoomerAMG preconditioner, ";
1057 }
1058
1063 }
1064 else if (Internal_preconditioner == Euclid) // Euclid
1065 {
1066 if (Output_info)
1067 {
1068 oomph_info << " with Euclid ILU preconditioner, ";
1069 }
1070
1075 }
1076 else if (Internal_preconditioner == ParaSails) // ParaSails
1077 {
1078 if (Output_info)
1079 {
1080 oomph_info << " with ParaSails approximate inverse preconditioner, ";
1081 }
1082
1087 }
1088 else
1089 {
1090 if (Output_info)
1091 {
1092 oomph_info << " with no preconditioner, ";
1093 }
1094 }
1095
1100
1102 }
1103
1104 // no solver exists for this value of Solver flag
1105 else
1106 {
1107 std::ostringstream error_message;
1108 error_message << "Solver has been set to an invalid value. "
1109 << "current value=" << Solver;
1110 throw OomphLibError(
1111 error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1112 }
1113
1115 double solver_setup_time = t_end - t_start;
1116
1117 // destroy dummy hypre vectors
1120
1121 // check error flag
1123 {
1124 std::ostringstream message;
1126 if (err)
1127 {
1129 "HypreSolver::hypre_solver_setup()",
1131 }
1132 }
1133
1134 // output times
1135 if (Output_info)
1136 {
1137 oomph_info << "time for setup [s] : " << solver_setup_time << std::endl;
1138 }
1139
1140 // doc parameters
1141 // doc_hypre_parameters();
1142 }
1143
1144
1145 //===================================================================
1146 /// Helper function performs a solve if solver data has been set
1147 /// up using hypre_solver_setup(...).
1148 //====================================================================
1151 {
1152 // Record time
1153 double t_start = TimingHelpers::timer();
1154
1155 // Set up hypre vectors
1156 // --------------------
1157
1158 // Hypre vector for rhs values
1161
1162 // Hypre vector for solution
1165
1166 // set up rhs_values and vec_indices
1169
1172
1173 // check error flag
1175 {
1176 std::ostringstream message;
1178 if (err)
1179 {
1181 "HypreSolver::hypre_solve()",
1183 }
1184 }
1185
1186 // solve
1187 // -----
1188
1189 // for solver stats
1190 int iterations = 0;
1191 double norm = 0;
1192
1193 // Get the norm of rhs
1194 const double rhs_norm = rhs.norm();
1195 bool do_solving = false;
1196 if (rhs_norm > 0.0)
1197 {
1198 do_solving = true;
1199 }
1200
1201#ifdef OOMPH_HAS_MPI
1202 // We need to check whether any processor requires to solve, if that
1203 // is the case then do the solving
1205 {
1207 {
1208 unsigned this_processor_do_solving = 0;
1209 unsigned all_processors_do_solving = 0;
1210 if (do_solving)
1211 {
1213 }
1214 // Get the communicator
1216 // Communicate with all procesoors
1219 1,
1221 MPI_SUM,
1222 comm_pt->mpi_comm());
1224 {
1225 do_solving = true;
1226 }
1227 }
1228 }
1229#endif
1230
1231 if (do_solving)
1232 {
1234 {
1238 }
1239 else if (Existing_solver == CG)
1240 {
1245 HYPRE_PCGGetNumIterations(Solver, &iterations);
1247 }
1248 else if (Existing_solver == GMRES)
1249 {
1254 HYPRE_GMRESGetNumIterations(Solver, &iterations);
1256 }
1257 else if (Existing_solver == BiCGStab)
1258 {
1265 }
1266 else if (Existing_solver == Euclid)
1267 {
1269 }
1270 else if (Existing_solver == ParaSails)
1271 {
1273 }
1274
1275 // output any error message
1277 {
1278 std::ostringstream message;
1280 if (err)
1281 {
1283 "HypreSolver::hypre_solve()",
1285 }
1286 }
1287
1288 } // if (do_solving)
1289
1290 // Copy result to solution
1291 unsigned nrow_local = Hypre_distribution_pt->nrow_local();
1292 unsigned first_row = Hypre_distribution_pt->first_row();
1293 int* indices = new int[nrow_local];
1294 for (unsigned i = 0; i < nrow_local; i++)
1295 {
1296 indices[i] = first_row + i;
1297 }
1299 if (solution.built())
1300 {
1301 soln_dist_pt = new LinearAlgebraDistribution(solution.distribution_pt());
1302 }
1303 else
1304 {
1305 soln_dist_pt = new LinearAlgebraDistribution(rhs.distribution_pt());
1306 }
1307 solution.build(Hypre_distribution_pt, 0.0);
1309 solution_ij, nrow_local, indices, solution.values_pt());
1310 solution.redistribute(soln_dist_pt);
1311 delete[] indices;
1312 delete soln_dist_pt;
1313
1314 // output any error message
1316 {
1317 std::ostringstream message;
1319 if (err)
1320 {
1322 "HypreSolver::hypre_solve()",
1324 }
1325 }
1326
1327 // deallocation
1330
1331 // Record time
1332 double solve_time = 0;
1333 if (Output_info)
1334 {
1335 double t_end = TimingHelpers::timer();
1337 }
1338
1339 // output timings and info
1340 if (Output_info)
1341 {
1342 oomph_info << "Time for HYPRE solve [s] : " << solve_time << std::endl;
1343 }
1344
1345 // for iterative solvers output iterations and final norm
1346 if ((Hypre_method >= CG) && (Hypre_method <= BoomerAMG))
1347 {
1348 if (iterations > 1)
1349 {
1350 if (Output_info)
1351 oomph_info << "Number of iterations : " << iterations
1352 << std::endl;
1353 if (Output_info)
1354 oomph_info << "Final Relative Residual Norm : " << norm << std::endl;
1355 }
1356 }
1357 }
1358
1359
1360 //===================================================================
1361 /// hypre_clean_up_memory() deletes any existing Hypre solver and
1362 /// Hypre matrix
1363 //====================================================================
1365 {
1366 // is there an existing solver
1367 if (Existing_solver != None)
1368 {
1369 // delete matrix
1371
1372 // delete solver
1374 {
1376 }
1377 else if (Existing_solver == CG)
1378 {
1380 }
1381 else if (Existing_solver == GMRES)
1382 {
1384 }
1385 else if (Existing_solver == BiCGStab)
1386 {
1388 }
1389 else if (Existing_solver == Euclid)
1390 {
1392 }
1393 else if (Existing_solver == ParaSails)
1394 {
1396 }
1398
1399 // delete preconditioner
1401 {
1403 }
1404 else if (Existing_preconditioner == Euclid)
1405 {
1407 }
1409 {
1411 }
1413
1414 // check error flag
1416 {
1417 std::ostringstream message;
1419 if (err)
1420 {
1422 "HypreSolver::clean_up_memory()",
1424 }
1425 }
1426 }
1427 }
1428
1429
1430 ///////////////////////////////////////////////////////////////////////
1431 ///////////////////////////////////////////////////////////////////////
1432 // functions for HypreSolver class
1433 ///////////////////////////////////////////////////////////////////////
1434 ///////////////////////////////////////////////////////////////////////
1435
1436
1437 //===================================================================
1438 /// Problem-based solve function to generate the Jacobian matrix and
1439 /// residual vector and use HypreInterface::hypre_solver_setup(...)
1440 /// and HypreInterface::hypre_solve(...) to solve the linear system.
1441 /// This function will delete any existing data.
1442 /// Note: the returned solution vector is NOT distributed, i.e. all
1443 /// processors hold all values because this is what the Newton solver
1444 /// requires.
1445 //====================================================================
1447 {
1448 double t_start = TimingHelpers::timer();
1449
1450 // Set Output_time flag for HypreInterface
1452
1453 // Delete any existing solver data
1455
1456 // Set flag to allow deletion of the oomphlib Jacobian matrix
1457 // (we're in control)
1458 Delete_input_data = true;
1459
1460 // Get oomph-lib Jacobian matrix and residual vector
1461 DoubleVector residual;
1463 problem_pt->get_jacobian(residual, *matrix);
1464
1465 // Output times
1466 if (Doc_time)
1467 {
1468 oomph_info << "Time to generate Jacobian and residual [s] : ";
1469 double t_end = TimingHelpers::timer();
1470 oomph_info << t_end - t_start << std::endl;
1471 }
1472
1473 // generate hypre matrix
1475
1476 // call hypre_solver_setup to generate linear solver data
1478
1479 // perform hypre_solve
1480 hypre_solve(residual, solution);
1481
1482 // delete solver data if required
1483 if (!Enable_resolve)
1484 {
1486 }
1487 }
1488
1489
1490 //===================================================================
1491 /// Uses HypreInterface::hypre_solve(...) to solve the linear system
1492 /// for a CRDoubleMatrix or a DistributedCRDoubleMatrix.
1493 /// In the latter case, the rhs needs to be distributed too,
1494 /// i.e. the length of the rhs vector must be equal to the number of
1495 /// rows stored locally.
1496 /// Note: the returned solution vector is never distributed, i.e. all
1497 /// processors hold all values
1498 //====================================================================
1499 void HypreSolver::solve(DoubleMatrixBase* const& matrix_pt,
1500 const DoubleVector& rhs,
1502 {
1503#ifdef PARANOID
1504 // check the matrix is square
1505 if (matrix_pt->nrow() != matrix_pt->ncol())
1506 {
1507 std::ostringstream error_message;
1508 error_message << "HypreSolver require a square matrix. "
1509 << "Matrix is " << matrix_pt->nrow() << " by "
1510 << matrix_pt->ncol() << std::endl;
1511 throw OomphLibError(
1512 error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1513 }
1514#endif
1515
1516 // Set Output_time flag for HypreInterface
1518
1519 // Clean up existing solver data
1521
1522 // Set flag to decide if oomphlib matrix can be deleted
1523 // (Recall that Delete_matrix defaults to false).
1525
1526 // Try cast to a CRDoubleMatrix
1527 CRDoubleMatrix* cr_matrix_pt = dynamic_cast<CRDoubleMatrix*>(matrix_pt);
1528
1529 // If cast successful set things up for a serial solve
1530 if (cr_matrix_pt)
1531 {
1532 // rebuild the distribution
1533 this->build_distribution(cr_matrix_pt->distribution_pt());
1534
1535#ifdef PARANOID
1536 // check that rhs has the same distribution as the matrix (now stored
1537 // in Distribution_pt)
1538 if (*this->distribution_pt() != *rhs.distribution_pt())
1539 {
1540 std::ostringstream error_message;
1541 error_message << "The distribution of the rhs vector and the matrix "
1542 << " must be the same" << std::endl;
1543 throw OomphLibError(error_message.str(),
1546 }
1547 // if the solution is setup make sure it has the same distribution as
1548 // the matrix as well
1549 if (solution.built())
1550 {
1551 if (*this->distribution_pt() != *solution.distribution_pt())
1552 {
1553 std::ostringstream error_message;
1554 error_message << "The distribution of the solution vector is setup "
1555 << "there it must be the same as the matrix."
1556 << std::endl;
1557 throw OomphLibError(error_message.str(),
1560 }
1561 }
1562#endif
1563
1565 }
1566 else
1567 {
1568#ifdef PARANOID
1569 std::ostringstream error_message;
1570 error_message << "HypreSolver only work with "
1571 << "CRDoubleMatrix matrices" << std::endl;
1572 throw OomphLibError(
1573 error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1574#endif
1575 }
1576
1577 // call hypre_setup to generate Hypre matrix and linear solver data
1579
1580 // perform hypre_solve
1582
1583 // delete solver data if required
1584 if (!Enable_resolve)
1585 {
1587 }
1588 }
1589
1590
1591 //===================================================================
1592 /// Resolve performs a linear solve using current solver data (if
1593 /// such data exists).
1594 //====================================================================
1596 {
1597#ifdef PARANOID
1598 // check solver data exists
1599 if (existing_solver() == None)
1600 {
1601 std::ostringstream error_message;
1602 error_message << "resolve(...) requires that solver data has been "
1603 << "set up by a previous call to solve(...) after "
1604 << "a call to enable_resolve()" << std::endl;
1605 throw OomphLibError(
1606 error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1607 }
1608 // check that rhs has the same distribution as the matrix (now stored
1609 // in Distribution_pt)
1610 if (*this->distribution_pt() != *rhs.distribution_pt())
1611 {
1612 std::ostringstream error_message;
1613 error_message << "The distribution of the rhs vector and the matrix "
1614 << " must be the same" << std::endl;
1615 throw OomphLibError(
1616 error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1617 }
1618 // if the solution is setup make sure it has the same distribution as
1619 // the matrix as well
1620 if (solution.built())
1621 {
1622 if (*this->distribution_pt() != *solution.distribution_pt())
1623 {
1624 std::ostringstream error_message;
1625 error_message << "The distribution of the solution vector is setup "
1626 << "there it must be the same as the matrix."
1627 << std::endl;
1628 throw OomphLibError(error_message.str(),
1631 }
1632 }
1633#endif
1634
1635 // Set Output_info flag for HypreInterface
1637
1638 // solve
1640
1641 // Note: never delete solver data as the preconditioner is typically
1642 // called repeatedly.
1643 }
1644
1645
1646 //===================================================================
1647 /// clean_up_memory() deletes any existing solver data
1648 //====================================================================
1653
1654
1655 ///////////////////////////////////////////////////////////////////////
1656 ///////////////////////////////////////////////////////////////////////
1657 // functions for HyprePreconditioner class
1658 ///////////////////////////////////////////////////////////////////////
1659 ///////////////////////////////////////////////////////////////////////
1660
1661
1662 //=============================================================================
1663 /// Static double that accumulates the preconditioner
1664 /// solve time of all instantiations of this class. Reset
1665 /// it manually, e.g. after every Newton solve, using
1666 /// reset_cumulative_solve_times().
1667 //=============================================================================
1669
1670 //=============================================================================
1671 /// map of static doubles that accumulates the preconditioner
1672 /// solve time of all instantiations of this class, labeled by
1673 /// context string. Reset
1674 /// it manually, e.g. after every Newton solve, using
1675 /// reset_cumulative_solve_times().
1676 //=============================================================================
1677 std::map<std::string, double>
1679
1680 //=============================================================================
1681 /// Static unsigned that accumulates the number of preconditioner
1682 /// solves of all instantiations of this class. Reset
1683 /// it manually, e.g. after every Newton solve, using
1684 /// reset_cumulative_solve_times().
1685 //=============================================================================
1687
1688 //=============================================================================
1689 /// Static unsigned that accumulates the number of preconditioner
1690 /// solves of all instantiations of this class, labeled by
1691 /// context string. Reset
1692 /// it manually, e.g. after every Newton solve, using
1693 /// reset_cumulative_solve_times().
1694 //=============================================================================
1695 std::map<std::string, unsigned>
1697
1698 //=============================================================================
1699 /// Static unsigned that stores nrow for the most recent
1700 /// instantiations of this class, labeled by
1701 /// context string. Reset
1702 /// it manually, e.g. after every Newton solve, using
1703 /// reset_cumulative_solve_times().
1704 //=============================================================================
1705 std::map<std::string, unsigned> HyprePreconditioner::Context_based_nrow;
1706
1707 //=============================================================================
1708 /// Report cumulative solve times of all instantiations of this
1709 /// class
1710 //=============================================================================
1712 {
1713 oomph_info << "\n\n=====================================================\n";
1714 oomph_info << "Cumulative HyprePreconditioner solve time "
1716 << " for " << Cumulative_npreconditioner_solve << " solves";
1718 {
1719 oomph_info << " ( "
1722 << " per solve )";
1723 }
1724 oomph_info << std::endl << std::endl;
1726 {
1727 oomph_info << "Breakdown by context: " << std::endl;
1728 for (std::map<std::string, double>::iterator it =
1731 it++)
1732 {
1734 << (*it).first << " " << (*it).second << " for "
1736 << " solves";
1738 {
1740 << " ( "
1741 << (*it).second /
1742 double(
1744 << " per solve; "
1745 << (*it).second /
1746 double(
1748 .first]) /
1749 double(Context_based_nrow[(*it).first])
1750 << " per solve per dof )";
1751 }
1752 oomph_info << std::endl;
1753 }
1754 }
1755 oomph_info << "\n=====================================================\n";
1756 oomph_info << std::endl;
1757 }
1758
1759 //=============================================================================
1760 /// Reset cumulative solve times
1761 //=============================================================================
1770
1771
1772 //=============================================================================
1773 /// An interface to allow HypreSolver to be used as a Preconditioner
1774 /// for the oomph-lib IterativeLinearSolver class.
1775 /// Matrix has to be of type CRDoubleMatrix or DistributedCRDoubleMatrix.
1776 //=============================================================================
1778 {
1779 // Set Output_info flag for HypreInterface
1781
1782#ifdef PARANOID
1783 // check the matrix is square
1784 if (matrix_pt()->nrow() != matrix_pt()->ncol())
1785 {
1786 std::ostringstream error_message;
1787 error_message << "HyprePreconditioner require a square matrix. "
1788 << "Matrix is " << matrix_pt()->nrow() << " by "
1789 << matrix_pt()->ncol() << std::endl;
1790 throw OomphLibError(
1791 error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1792 }
1793#endif
1794
1795 // clean up any previous solver data
1797
1798 // set flag to decide if oomphlib matrix can be deleted
1799 // (Recall that Delete_matrix defaults to false).
1801
1802 // Try casting to a CRDoubleMatrix
1804
1805 // If cast successful set things up for a serial solve
1806 if (cr_matrix_pt)
1807 {
1808 this->build_distribution(cr_matrix_pt->distribution_pt());
1810 }
1811 else
1812 {
1813#ifdef PARANOID
1814 std::ostringstream error_message;
1815 error_message << "HyprePreconditioner only work with "
1816 << "CRDoubleMatrix matrices" << std::endl;
1817 throw OomphLibError(
1818 error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1819#endif
1820 }
1821
1822 if (Context_string != "")
1823 {
1824 oomph_info << "Setup of HyprePreconditioner in context \" "
1825 << Context_string << "\": nrow, nrow_local, nnz "
1826 << cr_matrix_pt->nrow() << " " << cr_matrix_pt->nrow_local()
1827 << " " << cr_matrix_pt->nnz() << std::endl;
1828 }
1830
1831 // call hypre_solver_setup
1833 }
1834
1835 //===================================================================
1836 /// Preconditioner_solve uses a hypre solver to precondition vector r
1837 //====================================================================
1839 DoubleVector& z)
1840 {
1841 // Store time
1842 double t_start = TimingHelpers::timer();
1843
1844#ifdef PARANOID
1845 // check solver data exists
1846 if (existing_solver() == None)
1847 {
1848 std::ostringstream error_message;
1849 error_message << "preconditioner_solve(...) requires that data has "
1850 << "been set up using the function setup(...)" << std::endl;
1851 throw OomphLibError(
1852 error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1853 }
1854 // check that rhs has the same distribution as the matrix (now stored
1855 // in Distribution_pt)
1856 if (*this->distribution_pt() != *r.distribution_pt())
1857 {
1858 std::ostringstream error_message;
1859 error_message << "The distribution of the rhs vector and the matrix "
1860 << " must be the same" << std::endl;
1861 throw OomphLibError(
1862 error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1863 }
1864
1865 // if the solution is setup make sure it has the same distribution as
1866 // the matrix as well
1867 if (z.built())
1868 {
1869 if (*this->distribution_pt() != *z.distribution_pt())
1870 {
1871 std::ostringstream error_message;
1872 error_message << "The distribution of the solution vector is setup "
1873 << "there it must be the same as the matrix."
1874 << std::endl;
1875 throw OomphLibError(error_message.str(),
1878 }
1879 }
1880#endif
1881
1882 // Switch off any timings for the solve
1883 Output_info = false;
1884
1885 // perform hypre_solve
1886 hypre_solve(r, z);
1887
1888 // Add to cumulative solve time
1889 double t_end = TimingHelpers::timer();
1893 if (Context_string != "")
1894 {
1897 }
1898 }
1899
1900
1901 //===================================================================
1902 /// clean_up_memory() deletes any existing Hypre solver and
1903 /// Hypre matrix
1904 //====================================================================
1909
1910} // namespace oomph
cstr elem_len * i
Definition cfortran.h:603
The conjugate gradient method.
The conjugate gradient method.
A class for compressed row matrices. This is a distributable object.
Definition matrices.h:888
unsigned long nrow() const
Return the number of rows of the matrix.
Definition matrices.h:1002
unsigned nrow() const
access function to the number of global rows.
LinearAlgebraDistribution * distribution_pt() const
access to the LinearAlgebraDistribution
void build_distribution(const LinearAlgebraDistribution *const dist_pt)
setup the distribution of this distributable linear algebra object
Abstract base class for matrices of doubles – adds abstract interfaces for solving,...
Definition matrices.h:261
virtual unsigned long ncol() const =0
Return the number of columns of the matrix.
virtual unsigned long nrow() const =0
Return the number of rows of the matrix.
A vector in the mathematical sense, initially developed for linear algebra type applications....
double size() const
Calculate the size of the element (length, area, volume,...) in Eulerian computational coordinates....
Definition elements.cc:4320
The GMRES method.
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.
double AMG_strength
Connection strength threshold parameter for BoomerAMG.
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...
HYPRE_ParCSRMatrix Matrix_par
The Hypre_ParCSRMatrix version of the matrix used in solve(...), resolve(...) or preconditioner_solve...
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.
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...
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.
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.
static void reset_cumulative_solve_times()
Reset cumulative solve times.
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...
std::string Context_string
String can be used to provide context for annotation.
void clean_up_memory()
Function deletes all solver data.
static double Cumulative_preconditioner_solve_time
Static double that accumulates the preconditioner solve time of all instantiations of this class....
void setup()
Function to set up a preconditioner for the linear system defined by matrix_pt. This function is requ...
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,...
bool Delete_matrix
Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own da...
static std::map< std::string, unsigned > Context_based_nrow
Static unsigned that stores nrow for the most recent instantiations of this class,...
double My_cumulative_preconditioner_solve_time
Private double that accumulates the preconditioner solve time of thi instantiation of this class....
static void report_cumulative_solve_times()
Report cumulative solve times of all instantiations of this class.
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(....
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...
void clean_up_memory()
Function deletes all solver data.
bool Delete_matrix
Hypre copies matrix data from oomph-lib's CRDoubleMatrix or DistributedCRDoubleMatrix into its own da...
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...
Describes the distribution of a distributable linear algebra type object. Typically this is a contain...
unsigned first_row() const
access function for the first row on this processor. If not distributed then this is just zero.
OomphCommunicator * communicator_pt() const
const access to the communicator pointer
unsigned nrow_local() const
access function for the num of local rows on this processor. If no MPI then Nrow is returned.
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 oomph-lib wrapper to the MPI_Comm communicator object. Just contains an MPI_Comm object (which is ...
An OomphLibError object which should be thrown when an run-time error is encountered....
An OomphLibWarning object which should be created as a temporary object to issue a warning....
Preconditioner base class. Gives an interface to call all other preconditioners through and stores th...
virtual DoubleMatrixBase * matrix_pt() const
Get function for matrix pointer.
////////////////////////////////////////////////////////////////// //////////////////////////////////...
Definition problem.h:154
virtual void get_jacobian(DoubleVector &residuals, DenseDoubleMatrix &jacobian)
Return the fully-assembled Jacobian and residuals for the problem Interface for the case when the Jac...
Definition problem.cc:3935
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.
double timer()
returns the time in seconds after some point in past
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
std::string lower(const std::string &text)
Definition fpdiff.cc:123
OomphInfo oomph_info
Single (global) instantiation of the OomphInfo object – this is used throughout the library as a "rep...