]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add useful methods to TrilinosWrappers::PreconditionAMG::AdditionalData.
authorJean-Paul Pelteret <jppelteret@gmail.com>
Tue, 12 Mar 2019 08:35:07 +0000 (09:35 +0100)
committerJean-Paul Pelteret <jppelteret@gmail.com>
Tue, 12 Mar 2019 08:42:00 +0000 (09:42 +0100)
These include:
- the initialization of a parameter list, which can be later modified by
the
user, and
- the initialisation of the null space settings in an existing parameter
list.

doc/news/changes/minor/20190310Jean-PaulPelteret-1 [new file with mode: 0644]
include/deal.II/lac/trilinos_precondition.h
source/lac/trilinos_precondition_ml.cc

diff --git a/doc/news/changes/minor/20190310Jean-PaulPelteret-1 b/doc/news/changes/minor/20190310Jean-PaulPelteret-1
new file mode 100644 (file)
index 0000000..6030de9
--- /dev/null
@@ -0,0 +1,6 @@
+Changed: The TrilinosWrappers::PreconditionAMG::AdditionalData data structure
+is now able to return a parameter list, which can be adjusted and fine-tuned by 
+the user and later used to initialize the AMG preconditioner. It can also
+initialize the null space settings of an existing parameter list.
+<br>
+(Jean-Paul Pelteret, 2019/03/10)
index cb6d2ee11b6297f8131a03c3f729d85ead49f72c..aba329ce8a6ec1c6791d89e03c94f68bf694bafc 100644 (file)
@@ -34,6 +34,7 @@
 #      include <Epetra_SerialComm.h>
 #    endif
 #    include <Epetra_Map.h>
+#    include <Epetra_MultiVector.h>
 #    include <Epetra_RowMatrix.h>
 #    include <Epetra_Vector.h>
 #    include <Teuchos_ParameterList.hpp>
@@ -1382,6 +1383,77 @@ namespace TrilinosWrappers
                      const char *       smoother_type    = "Chebyshev",
                      const char *       coarse_type      = "Amesos-KLU");
 
+      /**
+       * Fill in a @p parameter_list that can be used to initialize the
+       * AMG preconditioner.
+       *
+       * The @p matrix is used in conjunction with the @p constant_modes to
+       * configure the null space settings for the preconditioner.
+       * The @p distributed_constant_modes are initialized by this function, and
+       * must remain in scope until PreconditionAMG::initialize() has been
+       * called.
+       *
+       * @note The set parameters reflect the current settings in this
+       * object, with various options being set both directly though the state
+       * of the member variables (e.g. the "smoother: type") as well as
+       * indirectly (e.g. the "aggregation: type"). If you wish to have
+       * fine-grained control over the configuration of the AMG preconditioner,
+       * then you can create the parameter list using this function (which
+       * conveniently sets the null space of the operator), change the relevant
+       * settings, and use the amended parameters list as an argument to
+       * PreconditionAMG::initialize(), instead of the AdditionalData object
+       * itself.
+       *
+       * See the documentation for the
+       * <a
+       * href="https://trilinos.org/docs/dev/packages/ml/doc/html/index.html">
+       * Trilinos ML package</a> for details on what options are available for
+       * modification.
+       *
+       * @note Any user-defined parameters that are not in conflict with those
+       * set by this data structure will be retained.
+       */
+      void
+      set_parameters(
+        Teuchos::ParameterList &             parameter_list,
+        std::unique_ptr<Epetra_MultiVector> &distributed_constant_modes,
+        const Epetra_RowMatrix &             matrix) const;
+
+      /**
+       * Fill in a parameter list that can be used to initialize the
+       * AMG preconditioner.
+       *
+       * @note Any user-defined parameters that are not in conflict with those
+       * set by this data structure will be retained.
+       */
+      void
+      set_parameters(
+        Teuchos::ParameterList &             parameter_list,
+        std::unique_ptr<Epetra_MultiVector> &distributed_constant_modes,
+        const SparseMatrix &                 matrix) const;
+
+      /**
+       * Configure the null space setting in the @p parameter_list for
+       * the input @p matrix based on the state of the @p constant_modes
+       * variable.
+       */
+      void
+      set_operator_null_space(
+        Teuchos::ParameterList &             parameter_list,
+        std::unique_ptr<Epetra_MultiVector> &distributed_constant_modes,
+        const Epetra_RowMatrix &             matrix) const;
+
+      /**
+       * Configure the null space setting in the @p parameter_list for
+       * the input @p matrix based on the state of the @p constant_modes
+       * variable.
+       */
+      void
+      set_operator_null_space(
+        Teuchos::ParameterList &             parameter_list,
+        std::unique_ptr<Epetra_MultiVector> &distributed_constant_modes,
+        const SparseMatrix &                 matrix) const;
+
       /**
        * Determines whether the AMG preconditioner should be optimized for
        * elliptic problems (ML option smoothed aggregation SA, using a
index 0375697c39203a688d191d75cb26f3568aa5bed9..848315efaa74bc2cec138ceac816099877e1b8c9 100644 (file)
@@ -63,31 +63,14 @@ namespace TrilinosWrappers
   {}
 
 
-  PreconditionAMG::~PreconditionAMG()
-  {
-    preconditioner.reset();
-    trilinos_matrix.reset();
-  }
-
-
 
   void
-  PreconditionAMG::initialize(const SparseMatrix &  matrix,
-                              const AdditionalData &additional_data)
-  {
-    initialize(matrix.trilinos_matrix(), additional_data);
-  }
-
-
-
-  void
-  PreconditionAMG::initialize(const Epetra_RowMatrix &matrix,
-                              const AdditionalData &  additional_data)
+  PreconditionAMG::AdditionalData::set_parameters(
+    Teuchos::ParameterList &             parameter_list,
+    std::unique_ptr<Epetra_MultiVector> &distributed_constant_modes,
+    const Epetra_RowMatrix &             matrix) const
   {
-    // Build the AMG preconditioner.
-    Teuchos::ParameterList parameter_list;
-
-    if (additional_data.elliptic == true)
+    if (elliptic == true)
       {
         ML_Epetra::SetDefaults("SA", parameter_list);
 
@@ -99,7 +82,7 @@ namespace TrilinosWrappers
         // standard choice uncoupled. if higher order, right now we also just
         // use Uncoupled, but we should be aware that maybe MIS might be
         // needed
-        if (additional_data.higher_order_elements)
+        if (higher_order_elements)
           parameter_list.set("aggregation: type", "Uncoupled");
       }
     else
@@ -109,8 +92,8 @@ namespace TrilinosWrappers
         parameter_list.set("aggregation: block scaling", true);
       }
 
-    parameter_list.set("smoother: type", additional_data.smoother_type);
-    parameter_list.set("coarse: type", additional_data.coarse_type);
+    parameter_list.set("smoother: type", smoother_type);
+    parameter_list.set("coarse: type", coarse_type);
 
     // Force re-initialization of the random seed to make ML deterministic
     // (only supported in trilinos >12.2):
@@ -118,34 +101,43 @@ namespace TrilinosWrappers
     parameter_list.set("initialize random seed", true);
 #  endif
 
-    parameter_list.set("smoother: sweeps",
-                       static_cast<int>(additional_data.smoother_sweeps));
-    parameter_list.set("cycle applications",
-                       static_cast<int>(additional_data.n_cycles));
-    if (additional_data.w_cycle == true)
+    parameter_list.set("smoother: sweeps", static_cast<int>(smoother_sweeps));
+    parameter_list.set("cycle applications", static_cast<int>(n_cycles));
+    if (w_cycle == true)
       parameter_list.set("prec type", "MGW");
     else
       parameter_list.set("prec type", "MGV");
 
     parameter_list.set("smoother: Chebyshev alpha", 10.);
     parameter_list.set("smoother: ifpack overlap",
-                       static_cast<int>(additional_data.smoother_overlap));
-    parameter_list.set("aggregation: threshold",
-                       additional_data.aggregation_threshold);
+                       static_cast<int>(smoother_overlap));
+    parameter_list.set("aggregation: threshold", aggregation_threshold);
     parameter_list.set("coarse: max size", 2000);
 
-    if (additional_data.output_details)
+    if (output_details)
       parameter_list.set("ML output", 10);
     else
       parameter_list.set("ML output", 0);
 
+    set_operator_null_space(parameter_list, distributed_constant_modes, matrix);
+  }
+
+
+
+  void
+  PreconditionAMG::AdditionalData::set_operator_null_space(
+    Teuchos::ParameterList &             parameter_list,
+    std::unique_ptr<Epetra_MultiVector> &ptr_distributed_constant_modes,
+    const Epetra_RowMatrix &             matrix) const
+  {
     const Epetra_Map &domain_map = matrix.OperatorDomainMap();
 
-    const size_type constant_modes_dimension =
-      additional_data.constant_modes.size();
-    Epetra_MultiVector distributed_constant_modes(
-      domain_map, constant_modes_dimension > 0 ? constant_modes_dimension : 1);
-    std::vector<double> dummy(constant_modes_dimension);
+    const size_type constant_modes_dimension = constant_modes.size();
+    ptr_distributed_constant_modes.reset(new Epetra_MultiVector(
+      domain_map, constant_modes_dimension > 0 ? constant_modes_dimension : 1));
+    Assert(ptr_distributed_constant_modes, ExcNotInitialized());
+    Epetra_MultiVector &distributed_constant_modes =
+      *ptr_distributed_constant_modes;
 
     if (constant_modes_dimension > 0)
       {
@@ -159,7 +151,7 @@ namespace TrilinosWrappers
                                     TrilinosWrappers::global_length(
                                       distributed_constant_modes)));
         const bool constant_modes_are_global =
-          additional_data.constant_modes[0].size() == global_size;
+          constant_modes[0].size() == global_size;
         const size_type my_size = domain_map.NumMyElements();
 
         // Reshape null space as a contiguous vector of doubles so that
@@ -168,10 +160,9 @@ namespace TrilinosWrappers
           constant_modes_are_global ? global_size : my_size;
         for (size_type d = 0; d < constant_modes_dimension; ++d)
           {
-            Assert(
-              additional_data.constant_modes[d].size() == expected_mode_size,
-              ExcDimensionMismatch(additional_data.constant_modes[d].size(),
-                                   expected_mode_size));
+            Assert(constant_modes[d].size() == expected_mode_size,
+                   ExcDimensionMismatch(constant_modes[d].size(),
+                                        expected_mode_size));
             for (size_type row = 0; row < my_size; ++row)
               {
                 const TrilinosWrappers::types::int_type mode_index =
@@ -179,7 +170,7 @@ namespace TrilinosWrappers
                     TrilinosWrappers::global_index(domain_map, row) :
                     row;
                 distributed_constant_modes[d][row] =
-                  additional_data.constant_modes[d][mode_index];
+                  constant_modes[d][mode_index];
               }
           }
         (void)expected_mode_size;
@@ -190,13 +181,75 @@ namespace TrilinosWrappers
         if (my_size > 0)
           parameter_list.set("null space: vectors",
                              distributed_constant_modes.Values());
-        // We need to set a valid pointer to data even if there is no data on
-        // the current processor. Therefore, pass a dummy in that case
         else
-          parameter_list.set("null space: vectors", dummy.data());
+          {
+            // We need to set a valid pointer to data even if there is no data
+            // on the current processor. Therefore, pass a dummy in that case
+            static std::vector<double> dummy;
+            if (dummy.size() != constant_modes_dimension)
+              dummy.resize(constant_modes_dimension);
+            parameter_list.set("null space: vectors", dummy.data());
+          }
       }
+  }
+
+
+
+  void
+  PreconditionAMG::AdditionalData::set_parameters(
+    Teuchos::ParameterList &             parameter_list,
+    std::unique_ptr<Epetra_MultiVector> &distributed_constant_modes,
+    const SparseMatrix &                 matrix) const
+  {
+    return set_parameters(parameter_list,
+                          distributed_constant_modes,
+                          matrix.trilinos_matrix());
+  }
+
+
+
+  void
+  PreconditionAMG::AdditionalData::set_operator_null_space(
+    Teuchos::ParameterList &             parameter_list,
+    std::unique_ptr<Epetra_MultiVector> &distributed_constant_modes,
+    const SparseMatrix &                 matrix) const
+  {
+    return set_operator_null_space(parameter_list,
+                                   distributed_constant_modes,
+                                   matrix.trilinos_matrix());
+  }
+
+
+
+  PreconditionAMG::~PreconditionAMG()
+  {
+    preconditioner.reset();
+    trilinos_matrix.reset();
+  }
+
+
+
+  void
+  PreconditionAMG::initialize(const SparseMatrix &  matrix,
+                              const AdditionalData &additional_data)
+  {
+    initialize(matrix.trilinos_matrix(), additional_data);
+  }
+
+
+
+  void
+  PreconditionAMG::initialize(const Epetra_RowMatrix &matrix,
+                              const AdditionalData &  additional_data)
+  {
+    // Build the AMG preconditioner.
+    Teuchos::ParameterList              ml_parameters;
+    std::unique_ptr<Epetra_MultiVector> distributed_constant_modes;
+    additional_data.set_parameters(ml_parameters,
+                                   distributed_constant_modes,
+                                   matrix);
 
-    initialize(matrix, parameter_list);
+    initialize(matrix, ml_parameters);
 
     if (additional_data.output_details)
       {

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.