]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add settings to BoomerAMG. 11991/head
authorMaximilian Bergbauer <bergbauer@lnm.mw.tum.de>
Tue, 30 Mar 2021 17:39:40 +0000 (19:39 +0200)
committerMaximilian Bergbauer <bergbauer@lnm.mw.tum.de>
Wed, 31 Mar 2021 17:06:13 +0000 (19:06 +0200)
doc/news/changes/minor/20210331MaximilianBergbauer [new file with mode: 0644]
include/deal.II/lac/petsc_precondition.h
source/lac/petsc_precondition.cc

diff --git a/doc/news/changes/minor/20210331MaximilianBergbauer b/doc/news/changes/minor/20210331MaximilianBergbauer
new file mode 100644 (file)
index 0000000..e493755
--- /dev/null
@@ -0,0 +1,4 @@
+Add: Adds new settings to PETScWrappers::PreconditionBoomerAMG::initialize()
+in the struct PETScWrappers::PreconditionBoomerAMG::AdditionalData.
+<br>
+(Maximilian Bergbauer, 2021/03/31)
index 6e4f36dc4ff4515e4fe859ba2617f384613071bd..e0a677b387d0b5d8a86d7d3307f8a424cf497e1a 100644 (file)
@@ -684,21 +684,51 @@ namespace PETScWrappers
      */
     struct AdditionalData
     {
+      /**
+       * Defines the available relaxation types for BoomerAMG.
+       */
+      enum class RelaxationType
+      {
+        Jacobi,
+        sequentialGaussSeidel,
+        seqboundaryGaussSeidel,
+        SORJacobi,
+        backwardSORJacobi,
+        symmetricSORJacobi,
+        l1scaledSORJacobi,
+        GaussianElimination,
+        l1GaussSeidel,
+        backwardl1GaussSeidel,
+        CG,
+        Chebyshev,
+        FCFJacobi,
+        l1scaledJacobi,
+        None
+      };
+
       /**
        * Constructor. Note that BoomerAMG offers a lot more options to set
        * than what is exposed here.
        */
-      AdditionalData(const bool         symmetric_operator = false,
-                     const double       strong_threshold   = 0.25,
-                     const double       max_row_sum        = 0.9,
-                     const unsigned int aggressive_coarsening_num_levels = 0,
-                     const bool         output_details = false);
+      AdditionalData(
+        const bool           symmetric_operator               = false,
+        const double         strong_threshold                 = 0.25,
+        const double         max_row_sum                      = 0.9,
+        const unsigned int   aggressive_coarsening_num_levels = 0,
+        const bool           output_details                   = false,
+        const RelaxationType relaxation_type_up   = RelaxationType::SORJacobi,
+        const RelaxationType relaxation_type_down = RelaxationType::SORJacobi,
+        const RelaxationType relaxation_type_coarse =
+          RelaxationType::GaussianElimination,
+        const unsigned int n_sweeps_coarse = 1,
+        const double       tol             = 0.0,
+        const unsigned int max_iter        = 1,
+        const bool         w_cycle         = false);
 
       /**
        * Set this flag to true if you have a symmetric system matrix and you
        * want to use a solver which assumes a symmetric preconditioner like
-       * CG. The relaxation is done with SSOR/Jacobi when set to true and with
-       * SOR/Jacobi otherwise.
+       * CG.
        */
       bool symmetric_operator;
 
@@ -734,6 +764,42 @@ namespace PETScWrappers
        * preconditioner is constructed.
        */
       bool output_details;
+
+      /**
+       * Choose relaxation type up.
+       */
+      RelaxationType relaxation_type_up;
+
+      /**
+       * Choose relaxation type down.
+       */
+      RelaxationType relaxation_type_down;
+
+      /**
+       * Choose relaxation type coarse.
+       */
+      RelaxationType relaxation_type_coarse;
+
+      /**
+       * Choose number of sweeps on coarse grid.
+       */
+      unsigned int n_sweeps_coarse;
+
+      /**
+       * Choose BommerAMG tolerance.
+       */
+      double tol;
+
+      /**
+       * Choose BommerAMG maximum number of cycles.
+       */
+      unsigned int max_iter;
+
+      /**
+       * Defines whether a w-cycle should be used instead of the standard
+       * setting of a v-cycle.
+       */
+      bool w_cycle;
     };
 
     /**
index 13aa14efb9134de4626fec4d3bdaa57c34cc1b12..510453baa7f751732e9dd1b576dcc93e0e6eda51 100644 (file)
@@ -435,20 +435,109 @@ namespace PETScWrappers
   /* ----------------- PreconditionBoomerAMG -------------------- */
 
   PreconditionBoomerAMG::AdditionalData::AdditionalData(
-    const bool         symmetric_operator,
-    const double       strong_threshold,
-    const double       max_row_sum,
-    const unsigned int aggressive_coarsening_num_levels,
-    const bool         output_details)
+    const bool           symmetric_operator,
+    const double         strong_threshold,
+    const double         max_row_sum,
+    const unsigned int   aggressive_coarsening_num_levels,
+    const bool           output_details,
+    const RelaxationType relaxation_type_up,
+    const RelaxationType relaxation_type_down,
+    const RelaxationType relaxation_type_coarse,
+    const unsigned int   n_sweeps_coarse,
+    const double         tol,
+    const unsigned int   max_iter,
+    const bool           w_cycle)
     : symmetric_operator(symmetric_operator)
     , strong_threshold(strong_threshold)
     , max_row_sum(max_row_sum)
     , aggressive_coarsening_num_levels(aggressive_coarsening_num_levels)
     , output_details(output_details)
+    , relaxation_type_up(relaxation_type_up)
+    , relaxation_type_down(relaxation_type_down)
+    , relaxation_type_coarse(relaxation_type_coarse)
+    , n_sweeps_coarse(n_sweeps_coarse)
+    , tol(tol)
+    , max_iter(max_iter)
+    , w_cycle(w_cycle)
   {}
 
 
 
+  namespace
+  {
+    /**
+     * Converts the enums for the different relaxation types to the respective
+     * strings for PETSc.
+     */
+    std::string
+    to_string(
+      PreconditionBoomerAMG::AdditionalData::RelaxationType relaxation_type)
+    {
+      std::string string_type;
+
+      switch (relaxation_type)
+        {
+          case PreconditionBoomerAMG::AdditionalData::RelaxationType::Jacobi:
+            string_type = "Jacobi";
+            break;
+          case PreconditionBoomerAMG::AdditionalData::RelaxationType::
+            sequentialGaussSeidel:
+            string_type = "sequential-Gauss-Seidel";
+            break;
+          case PreconditionBoomerAMG::AdditionalData::RelaxationType::
+            seqboundaryGaussSeidel:
+            string_type = "seqboundary-Gauss-Seidel";
+            break;
+          case PreconditionBoomerAMG::AdditionalData::RelaxationType::SORJacobi:
+            string_type = "SOR/Jacobi";
+            break;
+          case PreconditionBoomerAMG::AdditionalData::RelaxationType::
+            backwardSORJacobi:
+            string_type = "backward-SOR/Jacobi";
+            break;
+          case PreconditionBoomerAMG::AdditionalData::RelaxationType::
+            symmetricSORJacobi:
+            string_type = "symmetric-SOR/Jacobi";
+            break;
+          case PreconditionBoomerAMG::AdditionalData::RelaxationType::
+            l1scaledSORJacobi:
+            string_type = " l1scaled-SOR/Jacobi";
+            break;
+          case PreconditionBoomerAMG::AdditionalData::RelaxationType::
+            GaussianElimination:
+            string_type = "Gaussian-elimination";
+            break;
+          case PreconditionBoomerAMG::AdditionalData::RelaxationType::
+            l1GaussSeidel:
+            string_type = "l1-Gauss-Seidel";
+            break;
+          case PreconditionBoomerAMG::AdditionalData::RelaxationType::
+            backwardl1GaussSeidel:
+            string_type = "backward-l1-Gauss-Seidel";
+            break;
+          case PreconditionBoomerAMG::AdditionalData::RelaxationType::CG:
+            string_type = "CG";
+            break;
+          case PreconditionBoomerAMG::AdditionalData::RelaxationType::Chebyshev:
+            string_type = "Chebyshev";
+            break;
+          case PreconditionBoomerAMG::AdditionalData::RelaxationType::FCFJacobi:
+            string_type = "FCF-Jacobi";
+            break;
+          case PreconditionBoomerAMG::AdditionalData::RelaxationType::
+            l1scaledJacobi:
+            string_type = "l1scaled-Jacobi";
+            break;
+          case PreconditionBoomerAMG::AdditionalData::RelaxationType::None:
+            string_type = "None";
+            break;
+          default:
+            Assert(false, ExcNotImplemented());
+        }
+      return string_type;
+    }
+  } // namespace
+
   PreconditionBoomerAMG::PreconditionBoomerAMG(
     const MPI_Comm &      comm,
     const AdditionalData &additional_data_)
@@ -495,29 +584,84 @@ namespace PETScWrappers
                      std::to_string(
                        additional_data.aggressive_coarsening_num_levels));
 
-    std::stringstream ssStream;
-    ssStream << additional_data.max_row_sum;
-    set_option_value("-pc_hypre_boomeramg_max_row_sum", ssStream.str());
+    set_option_value("-pc_hypre_boomeramg_max_row_sum",
+                     std::to_string(additional_data.max_row_sum));
 
-    ssStream.str(""); // empty the stringstream
-    ssStream << additional_data.strong_threshold;
-    set_option_value("-pc_hypre_boomeramg_strong_threshold", ssStream.str());
+    set_option_value("-pc_hypre_boomeramg_strong_threshold",
+                     std::to_string(additional_data.strong_threshold));
 
-    if (additional_data.symmetric_operator)
+    // change to symmetric SOR/Jacobi when using a symmetric operator for
+    // backward compatibility
+    if (additional_data.relaxation_type_up ==
+          AdditionalData::RelaxationType::SORJacobi &&
+        additional_data.symmetric_operator)
       {
-        set_option_value("-pc_hypre_boomeramg_relax_type_up",
-                         "symmetric-SOR/Jacobi");
-        set_option_value("-pc_hypre_boomeramg_relax_type_down",
-                         "symmetric-SOR/Jacobi");
-        set_option_value("-pc_hypre_boomeramg_relax_type_coarse",
-                         "Gaussian-elimination");
+        additional_data.relaxation_type_up =
+          AdditionalData::RelaxationType::symmetricSORJacobi;
       }
-    else
+
+    if (additional_data.relaxation_type_down ==
+          AdditionalData::RelaxationType::SORJacobi &&
+        additional_data.symmetric_operator)
+      {
+        additional_data.relaxation_type_down =
+          AdditionalData::RelaxationType::symmetricSORJacobi;
+      }
+
+    if (additional_data.relaxation_type_coarse ==
+          AdditionalData::RelaxationType::SORJacobi &&
+        additional_data.symmetric_operator)
+      {
+        additional_data.relaxation_type_coarse =
+          AdditionalData::RelaxationType::symmetricSORJacobi;
+      }
+
+    auto relaxation_type_is_symmetric =
+      [](AdditionalData::RelaxationType relaxation_type) {
+        return relaxation_type == AdditionalData::RelaxationType::Jacobi ||
+               relaxation_type ==
+                 AdditionalData::RelaxationType::symmetricSORJacobi ||
+               relaxation_type ==
+                 AdditionalData::RelaxationType::GaussianElimination ||
+               relaxation_type == AdditionalData::RelaxationType::None ||
+               relaxation_type ==
+                 AdditionalData::RelaxationType::l1scaledJacobi ||
+               relaxation_type == AdditionalData::RelaxationType::CG ||
+               relaxation_type == AdditionalData::RelaxationType::Chebyshev;
+      };
+
+    if (additional_data.symmetric_operator &&
+        !relaxation_type_is_symmetric(additional_data.relaxation_type_up))
+      Assert(false,
+             ExcMessage("Use a symmetric smoother for relaxation_type_up"));
+
+    if (additional_data.symmetric_operator &&
+        !relaxation_type_is_symmetric(additional_data.relaxation_type_down))
+      Assert(false,
+             ExcMessage("Use a symmetric smoother for relaxation_type_down"));
+
+    if (additional_data.symmetric_operator &&
+        !relaxation_type_is_symmetric(additional_data.relaxation_type_coarse))
+      Assert(false,
+             ExcMessage("Use a symmetric smoother for relaxation_type_coarse"));
+
+    set_option_value("-pc_hypre_boomeramg_relax_type_up",
+                     to_string(additional_data.relaxation_type_up));
+    set_option_value("-pc_hypre_boomeramg_relax_type_down",
+                     to_string(additional_data.relaxation_type_down));
+    set_option_value("-pc_hypre_boomeramg_relax_type_coarse",
+                     to_string(additional_data.relaxation_type_coarse));
+    set_option_value("-pc_hypre_boomeramg_grid_sweeps_coarse",
+                     std::to_string(additional_data.n_sweeps_coarse));
+
+    set_option_value("-pc_hypre_boomeramg_tol",
+                     std::to_string(additional_data.tol));
+    set_option_value("-pc_hypre_boomeramg_max_iter",
+                     std::to_string(additional_data.max_iter));
+
+    if (additional_data.w_cycle)
       {
-        set_option_value("-pc_hypre_boomeramg_relax_type_up", "SOR/Jacobi");
-        set_option_value("-pc_hypre_boomeramg_relax_type_down", "SOR/Jacobi");
-        set_option_value("-pc_hypre_boomeramg_relax_type_coarse",
-                         "Gaussian-elimination");
+        set_option_value("-pc_hypre_boomeramg_cycle_type", "W");
       }
 
     ierr = PCSetFromOptions(pc);

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.