]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Remove MatrixType from MGTransferGlobalCoarsening 11708/head
authorPeter Munch <peterrmuench@gmail.com>
Sun, 7 Feb 2021 10:27:52 +0000 (11:27 +0100)
committerPeter Munch <peterrmuench@gmail.com>
Mon, 8 Feb 2021 10:17:07 +0000 (11:17 +0100)
include/deal.II/multigrid/mg_transfer_global_coarsening.h
tests/multigrid-global-coarsening/multigrid_a_01.cc
tests/multigrid-global-coarsening/multigrid_p_01.cc

index dc88c7e8fafc39a55bb0ae8066899d07944788ba..c69d210fb8588292b28fdf03cf15d52f01239ad1 100644 (file)
@@ -247,32 +247,26 @@ private:
  * one of these element, as well as, systems with different elements or other
  * elements are currently not implemented.
  */
-template <typename MatrixType, typename VectorType>
+template <int dim, typename VectorType>
 class MGTransferGlobalCoarsening : public dealii::MGTransferBase<VectorType>
 {
 public:
-  static_assert(std::is_same<typename MatrixType::value_type,
-                             typename VectorType::value_type>::value,
-                "Types do not match.");
-
-  /**
-   * Dimension.
-   */
-  static const int dim = MatrixType::dim;
-
   /**
    * Value type.
    */
-  using Number = typename MatrixType::value_type;
+  using Number = typename VectorType::value_type;
 
   /**
-   * Constructor taking an operator for each level (minimum requirement is
-   * that the operator provides the function initialize_dof_vector()) and
-   * transfer operators (with the coarsest level kept empty in @p transfer).
+   * Constructor taking a collection of transfer operators (with the coarsest
+   * level kept
+   * empty in @p transfer) and an optional function that initializes the
+   * internal level vectors within the function call copy_to_mg() if used in the
+   * context of PreconditionMG.
    */
   MGTransferGlobalCoarsening(
-    const MGLevelObject<MatrixType> &                         matrices,
-    const MGLevelObject<MGTwoLevelTransfer<dim, VectorType>> &transfer);
+    const MGLevelObject<MGTwoLevelTransfer<dim, VectorType>> &transfer,
+    const std::function<void(const unsigned int, VectorType &)>
+      &initialize_dof_vector = {});
 
   /**
    * Perform prolongation.
@@ -315,8 +309,16 @@ public:
                const MGLevelObject<VectorType> &src) const;
 
 private:
-  const MGLevelObject<MatrixType> &                         matrices;
+  /**
+   * Collection of the two-level transfer operators.
+   */
   const MGLevelObject<MGTwoLevelTransfer<dim, VectorType>> &transfer;
+
+  /**
+   * Function to initialize internal level vectors.
+   */
+  const std::function<void(const unsigned int, VectorType &)>
+    initialize_dof_vector;
 };
 
 
@@ -368,22 +370,20 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
 
 
 
-template <typename MatrixType, typename VectorType>
-MGTransferGlobalCoarsening<MatrixType, VectorType>::MGTransferGlobalCoarsening(
-  const MGLevelObject<MatrixType> &                         matrices,
-  const MGLevelObject<MGTwoLevelTransfer<dim, VectorType>> &transfer)
-  : matrices(matrices)
-  , transfer(transfer)
-{
-  AssertDimension(matrices.max_level() - matrices.min_level(),
-                  transfer.max_level() - transfer.min_level());
-}
+template <int dim, typename VectorType>
+MGTransferGlobalCoarsening<dim, VectorType>::MGTransferGlobalCoarsening(
+  const MGLevelObject<MGTwoLevelTransfer<dim, VectorType>> &transfer,
+  const std::function<void(const unsigned int, VectorType &)>
+    &initialize_dof_vector)
+  : transfer(transfer)
+  , initialize_dof_vector(initialize_dof_vector)
+{}
 
 
 
-template <typename MatrixType, typename VectorType>
+template <int dim, typename VectorType>
 void
-MGTransferGlobalCoarsening<MatrixType, VectorType>::prolongate(
+MGTransferGlobalCoarsening<dim, VectorType>::prolongate(
   const unsigned int to_level,
   VectorType &       dst,
   const VectorType & src) const
@@ -393,9 +393,9 @@ MGTransferGlobalCoarsening<MatrixType, VectorType>::prolongate(
 
 
 
-template <typename MatrixType, typename VectorType>
+template <int dim, typename VectorType>
 void
-MGTransferGlobalCoarsening<MatrixType, VectorType>::restrict_and_add(
+MGTransferGlobalCoarsening<dim, VectorType>::restrict_and_add(
   const unsigned int from_level,
   VectorType &       dst,
   const VectorType & src) const
@@ -405,10 +405,10 @@ MGTransferGlobalCoarsening<MatrixType, VectorType>::restrict_and_add(
 
 
 
-template <typename MatrixType, typename VectorType>
+template <int dim, typename VectorType>
 template <class InVector, int spacedim>
 void
-MGTransferGlobalCoarsening<MatrixType, VectorType>::copy_to_mg(
+MGTransferGlobalCoarsening<dim, VectorType>::copy_to_mg(
   const DoFHandler<dim, spacedim> &dof_handler,
   MGLevelObject<VectorType> &      dst,
   const InVector &                 src) const
@@ -416,17 +416,17 @@ MGTransferGlobalCoarsening<MatrixType, VectorType>::copy_to_mg(
   (void)dof_handler;
 
   for (unsigned int level = dst.min_level(); level <= dst.max_level(); ++level)
-    matrices[level].initialize_dof_vector(dst[level]);
+    initialize_dof_vector(level, dst[level]);
 
   dst[dst.max_level()].copy_locally_owned_data_from(src);
 }
 
 
 
-template <typename MatrixType, typename VectorType>
+template <int dim, typename VectorType>
 template <class OutVector, int spacedim>
 void
-MGTransferGlobalCoarsening<MatrixType, VectorType>::copy_from_mg(
+MGTransferGlobalCoarsening<dim, VectorType>::copy_from_mg(
   const DoFHandler<dim, spacedim> &dof_handler,
   OutVector &                      dst,
   const MGLevelObject<VectorType> &src) const
index 83470e4650f81616103a77cc8268ba49d92c322d..1e2053b142d4e7d086942a1dffbdc6c8c3ef4627 100644 (file)
@@ -100,8 +100,9 @@ test(const unsigned int n_refinements,
                                                constraints[l + 1],
                                                constraints[l]);
 
-  MGTransferGlobalCoarsening<Operator<dim, Number>, VectorType> transfer(
-    operators, transfers);
+  MGTransferGlobalCoarsening<dim, VectorType> transfer(
+    transfers,
+    [&](const auto l, auto &vec) { operators[l].initialize_dof_vector(vec); });
 
   GMGParameters mg_data; // TODO
 
index 623d35b62d87ea7eb984fa77893ffce5973bfc1f..863be260ed15f2f7bd5f557f04c9620213f0fcc6 100644 (file)
@@ -101,8 +101,9 @@ test(const unsigned int n_refinements,
                                                 constraints[l + 1],
                                                 constraints[l]);
 
-  MGTransferGlobalCoarsening<Operator<dim, Number>, VectorType> transfer(
-    operators, transfers);
+  MGTransferGlobalCoarsening<dim, VectorType> transfer(
+    transfers,
+    [&](const auto l, auto &vec) { operators[l].initialize_dof_vector(vec); });
 
   GMGParameters mg_data; // TODO
 

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.