]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Introduce MGTransferBlockMatrixFreeBase 13608/head
authorPeter Munch <peterrmuench@gmail.com>
Sun, 10 Apr 2022 14:03:26 +0000 (16:03 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Sun, 10 Apr 2022 14:03:26 +0000 (16:03 +0200)
include/deal.II/multigrid/mg_transfer_matrix_free.h
source/multigrid/mg_transfer_matrix_free.cc
source/multigrid/mg_transfer_matrix_free.inst.in

index 3e558e0954bd9d1eaa7f26b7e3c7226b2b79ef65..84ea6ba8b010c729ae976d50210784b86eb1bdf3 100644 (file)
@@ -310,78 +310,17 @@ private:
 };
 
 
+
 /**
- * Implementation of the MGTransferBase interface for which the transfer
- * operations is implemented in a matrix-free way based on the interpolation
- * matrices of the underlying finite element. This requires considerably less
- * memory than MGTransferPrebuilt and can also be considerably faster than
- * that variant.
- *
- * This class works with LinearAlgebra::distributed::BlockVector and
- * performs exactly the same transfer operations for each block as
- * MGTransferMatrixFree.
- * Both the cases that the same DoFHandler is used for all the blocks
- * and that each block uses its own DoFHandler are supported.
+ * Base class of MGTransferBlockMatrixFree. While MGTransferBlockMatrixFree
+ * constains all the setup routines of the transfer operators for the blocks,
+ * this class simply applies them, e.g., for restricting and prolongating.
  */
-template <int dim, typename Number>
-class MGTransferBlockMatrixFree
+template <int dim, typename Number, typename TransferType>
+class MGTransferBlockMatrixFreeBase
   : public MGTransferBase<LinearAlgebra::distributed::BlockVector<Number>>
 {
 public:
-  /**
-   * Constructor without constraint matrices. Use this constructor only with
-   * discontinuous finite elements or with no local refinement.
-   */
-  MGTransferBlockMatrixFree() = default;
-
-  /**
-   * Constructor with constraints. Equivalent to the default constructor
-   * followed by initialize_constraints().
-   */
-  MGTransferBlockMatrixFree(const MGConstrainedDoFs &mg_constrained_dofs);
-
-  /**
-   * Same as above for the case that each block has its own DoFHandler.
-   */
-  MGTransferBlockMatrixFree(
-    const std::vector<MGConstrainedDoFs> &mg_constrained_dofs);
-
-  /**
-   * Destructor.
-   */
-  virtual ~MGTransferBlockMatrixFree() override = default;
-
-  /**
-   * Initialize the constraints to be used in build().
-   */
-  void
-  initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs);
-
-  /**
-   * Same as above for the case that each block has its own DoFHandler.
-   */
-  void
-  initialize_constraints(
-    const std::vector<MGConstrainedDoFs> &mg_constrained_dofs);
-
-  /**
-   * Reset the object to the state it had right after the default constructor.
-   */
-  void
-  clear();
-
-  /**
-   * Actually build the information for the prolongation for each level.
-   */
-  void
-  build(const DoFHandler<dim, dim> &dof_handler);
-
-  /**
-   * Same as above for the case that each block has its own DoFHandler.
-   */
-  void
-  build(const std::vector<const DoFHandler<dim, dim> *> &dof_handler);
-
   /**
    * Prolongate a vector from level <tt>to_level-1</tt> to level
    * <tt>to_level</tt> using the embedding matrices of the underlying finite
@@ -481,29 +420,106 @@ public:
     const MGLevelObject<LinearAlgebra::distributed::BlockVector<Number>> &src)
     const;
 
-  /**
-   * Memory used by this object.
-   */
-  std::size_t
-  memory_consumption() const;
-
   /**
    * This class can both be used with a single DoFHandler
    * or a separate DoFHandler for each block.
    */
   static const bool supports_dof_handler_vector = true;
 
-private:
+protected:
   /**
    * Non-block matrix-free versions of transfer operation.
    */
-  std::vector<MGTransferMatrixFree<dim, Number>> matrix_free_transfer_vector;
+  std::vector<TransferType> matrix_free_transfer_vector;
 
   /**
    * A flag to indicate whether the same DoFHandler is used for all
    * the components or if each block has its own DoFHandler.
    */
-  const bool same_for_all;
+  bool same_for_all;
+};
+
+
+
+/**
+ * Implementation of the MGTransferBase interface for which the transfer
+ * operations is implemented in a matrix-free way based on the interpolation
+ * matrices of the underlying finite element. This requires considerably less
+ * memory than MGTransferPrebuilt and can also be considerably faster than
+ * that variant.
+ *
+ * This class works with LinearAlgebra::distributed::BlockVector and
+ * performs exactly the same transfer operations for each block as
+ * MGTransferMatrixFree.
+ * Both the cases that the same DoFHandler is used for all the blocks
+ * and that each block uses its own DoFHandler are supported.
+ */
+template <int dim, typename Number>
+class MGTransferBlockMatrixFree
+  : public MGTransferBlockMatrixFreeBase<dim,
+                                         Number,
+                                         MGTransferMatrixFree<dim, Number>>
+{
+public:
+  /**
+   * Constructor without constraint matrices. Use this constructor only with
+   * discontinuous finite elements or with no local refinement.
+   */
+  MGTransferBlockMatrixFree() = default;
+
+  /**
+   * Constructor with constraints. Equivalent to the default constructor
+   * followed by initialize_constraints().
+   */
+  MGTransferBlockMatrixFree(const MGConstrainedDoFs &mg_constrained_dofs);
+
+  /**
+   * Same as above for the case that each block has its own DoFHandler.
+   */
+  MGTransferBlockMatrixFree(
+    const std::vector<MGConstrainedDoFs> &mg_constrained_dofs);
+
+  /**
+   * Destructor.
+   */
+  virtual ~MGTransferBlockMatrixFree() override = default;
+
+  /**
+   * Initialize the constraints to be used in build().
+   */
+  void
+  initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs);
+
+  /**
+   * Same as above for the case that each block has its own DoFHandler.
+   */
+  void
+  initialize_constraints(
+    const std::vector<MGConstrainedDoFs> &mg_constrained_dofs);
+
+  /**
+   * Reset the object to the state it had right after the default constructor.
+   */
+  void
+  clear();
+
+  /**
+   * Actually build the information for the prolongation for each level.
+   */
+  void
+  build(const DoFHandler<dim, dim> &dof_handler);
+
+  /**
+   * Same as above for the case that each block has its own DoFHandler.
+   */
+  void
+  build(const std::vector<const DoFHandler<dim, dim> *> &dof_handler);
+
+  /**
+   * Memory used by this object.
+   */
+  std::size_t
+  memory_consumption() const;
 };
 
 
@@ -607,10 +623,10 @@ MGTransferMatrixFree<dim, Number>::interpolate_to_mg(
 
 
 
-template <int dim, typename Number>
+template <int dim, typename Number, typename TransferType>
 template <typename Number2, int spacedim>
 void
-MGTransferBlockMatrixFree<dim, Number>::copy_to_mg(
+MGTransferBlockMatrixFreeBase<dim, Number, TransferType>::copy_to_mg(
   const DoFHandler<dim, spacedim> &                               dof_handler,
   MGLevelObject<LinearAlgebra::distributed::BlockVector<Number>> &dst,
   const LinearAlgebra::distributed::BlockVector<Number2> &        src) const
@@ -629,10 +645,10 @@ MGTransferBlockMatrixFree<dim, Number>::copy_to_mg(
 
 
 
-template <int dim, typename Number>
+template <int dim, typename Number, typename TransferType>
 template <typename Number2, int spacedim>
 void
-MGTransferBlockMatrixFree<dim, Number>::copy_to_mg(
+MGTransferBlockMatrixFreeBase<dim, Number, TransferType>::copy_to_mg(
   const std::vector<const DoFHandler<dim, spacedim> *> &          dof_handler,
   MGLevelObject<LinearAlgebra::distributed::BlockVector<Number>> &dst,
   const LinearAlgebra::distributed::BlockVector<Number2> &        src) const
@@ -720,10 +736,10 @@ MGTransferBlockMatrixFree<dim, Number>::copy_to_mg(
     }
 }
 
-template <int dim, typename Number>
+template <int dim, typename Number, typename TransferType>
 template <typename Number2, int spacedim>
 void
-MGTransferBlockMatrixFree<dim, Number>::copy_from_mg(
+MGTransferBlockMatrixFreeBase<dim, Number, TransferType>::copy_from_mg(
   const DoFHandler<dim, spacedim> &                 dof_handler,
   LinearAlgebra::distributed::BlockVector<Number2> &dst,
   const MGLevelObject<LinearAlgebra::distributed::BlockVector<Number>> &src)
@@ -736,10 +752,10 @@ MGTransferBlockMatrixFree<dim, Number>::copy_from_mg(
   copy_from_mg(mg_dofs, dst, src);
 }
 
-template <int dim, typename Number>
+template <int dim, typename Number, typename TransferType>
 template <typename Number2, int spacedim>
 void
-MGTransferBlockMatrixFree<dim, Number>::copy_from_mg(
+MGTransferBlockMatrixFreeBase<dim, Number, TransferType>::copy_from_mg(
   const std::vector<const DoFHandler<dim, spacedim> *> &dof_handler,
   LinearAlgebra::distributed::BlockVector<Number2> &    dst,
   const MGLevelObject<LinearAlgebra::distributed::BlockVector<Number>> &src)
index 6cb316961776f00e2a76ebc1abf1bd1c770e6966..19fcb4295fadcd9491ebcfaf58adc6b4ca81962e 100644 (file)
@@ -707,9 +707,9 @@ MGTransferMatrixFree<dim, Number>::memory_consumption() const
 template <int dim, typename Number>
 MGTransferBlockMatrixFree<dim, Number>::MGTransferBlockMatrixFree(
   const MGConstrainedDoFs &mg_c)
-  : same_for_all(true)
 {
-  matrix_free_transfer_vector.emplace_back(mg_c);
+  this->same_for_all = true;
+  this->matrix_free_transfer_vector.emplace_back(mg_c);
 }
 
 
@@ -717,10 +717,10 @@ MGTransferBlockMatrixFree<dim, Number>::MGTransferBlockMatrixFree(
 template <int dim, typename Number>
 MGTransferBlockMatrixFree<dim, Number>::MGTransferBlockMatrixFree(
   const std::vector<MGConstrainedDoFs> &mg_c)
-  : same_for_all(false)
 {
+  this->same_for_all = false;
   for (const auto &constrained_block_dofs : mg_c)
-    matrix_free_transfer_vector.emplace_back(constrained_block_dofs);
+    this->matrix_free_transfer_vector.emplace_back(constrained_block_dofs);
 }
 
 
@@ -730,13 +730,13 @@ void
 MGTransferBlockMatrixFree<dim, Number>::initialize_constraints(
   const MGConstrainedDoFs &mg_c)
 {
-  Assert(same_for_all,
+  Assert(this->same_for_all,
          ExcMessage("This object was initialized with support for usage with "
                     "one DoFHandler for each block, but this method assumes "
                     "that the same DoFHandler is used for all the blocks!"));
-  AssertDimension(matrix_free_transfer_vector.size(), 1);
+  AssertDimension(this->matrix_free_transfer_vector.size(), 1);
 
-  matrix_free_transfer_vector[0].initialize_constraints(mg_c);
+  this->matrix_free_transfer_vector[0].initialize_constraints(mg_c);
 }
 
 
@@ -746,24 +746,26 @@ void
 MGTransferBlockMatrixFree<dim, Number>::initialize_constraints(
   const std::vector<MGConstrainedDoFs> &mg_c)
 {
-  Assert(!same_for_all,
+  Assert(!this->same_for_all,
          ExcMessage("This object was initialized with support for using "
                     "the same DoFHandler for all the blocks, but this "
                     "method assumes that there is a separate DoFHandler "
                     "for each block!"));
-  AssertDimension(matrix_free_transfer_vector.size(), mg_c.size());
+  AssertDimension(this->matrix_free_transfer_vector.size(), mg_c.size());
 
   for (unsigned int i = 0; i < mg_c.size(); ++i)
-    matrix_free_transfer_vector[i].initialize_constraints(mg_c[i]);
+    this->matrix_free_transfer_vector[i].initialize_constraints(mg_c[i]);
 }
 
 
 
 template <int dim, typename Number>
 void
-MGTransferBlockMatrixFree<dim, Number>::clear()
+MGTransferBlockMatrixFree<dim, Number>::build(
+  const DoFHandler<dim, dim> &dof_handler)
 {
-  matrix_free_transfer_vector.clear();
+  AssertDimension(this->matrix_free_transfer_vector.size(), 1);
+  this->matrix_free_transfer_vector[0].build(dof_handler);
 }
 
 
@@ -771,29 +773,39 @@ MGTransferBlockMatrixFree<dim, Number>::clear()
 template <int dim, typename Number>
 void
 MGTransferBlockMatrixFree<dim, Number>::build(
-  const DoFHandler<dim, dim> &dof_handler)
+  const std::vector<const DoFHandler<dim, dim> *> &dof_handler)
 {
-  AssertDimension(matrix_free_transfer_vector.size(), 1);
-  matrix_free_transfer_vector[0].build(dof_handler);
+  AssertDimension(this->matrix_free_transfer_vector.size(), dof_handler.size());
+  for (unsigned int i = 0; i < dof_handler.size(); ++i)
+    this->matrix_free_transfer_vector[i].build(*dof_handler[i]);
 }
 
 
 
 template <int dim, typename Number>
-void
-MGTransferBlockMatrixFree<dim, Number>::build(
-  const std::vector<const DoFHandler<dim, dim> *> &dof_handler)
+std::size_t
+MGTransferBlockMatrixFree<dim, Number>::memory_consumption() const
 {
-  AssertDimension(matrix_free_transfer_vector.size(), dof_handler.size());
-  for (unsigned int i = 0; i < dof_handler.size(); ++i)
-    matrix_free_transfer_vector[i].build(*dof_handler[i]);
+  std::size_t total_memory_consumption = 0;
+  for (const auto &el : this->matrix_free_transfer_vector)
+    total_memory_consumption += el.memory_consumption();
+  return total_memory_consumption;
 }
 
 
 
 template <int dim, typename Number>
 void
-MGTransferBlockMatrixFree<dim, Number>::prolongate(
+MGTransferBlockMatrixFree<dim, Number>::clear()
+{
+  this->matrix_free_transfer_vector.clear();
+}
+
+
+
+template <int dim, typename Number, typename TransferType>
+void
+MGTransferBlockMatrixFreeBase<dim, Number, TransferType>::prolongate(
   const unsigned int                                     to_level,
   LinearAlgebra::distributed::BlockVector<Number> &      dst,
   const LinearAlgebra::distributed::BlockVector<Number> &src) const
@@ -815,9 +827,9 @@ MGTransferBlockMatrixFree<dim, Number>::prolongate(
 
 
 
-template <int dim, typename Number>
+template <int dim, typename Number, typename TransferType>
 void
-MGTransferBlockMatrixFree<dim, Number>::prolongate_and_add(
+MGTransferBlockMatrixFreeBase<dim, Number, TransferType>::prolongate_and_add(
   const unsigned int                                     to_level,
   LinearAlgebra::distributed::BlockVector<Number> &      dst,
   const LinearAlgebra::distributed::BlockVector<Number> &src) const
@@ -839,9 +851,9 @@ MGTransferBlockMatrixFree<dim, Number>::prolongate_and_add(
 
 
 
-template <int dim, typename Number>
+template <int dim, typename Number, typename TransferType>
 void
-MGTransferBlockMatrixFree<dim, Number>::restrict_and_add(
+MGTransferBlockMatrixFreeBase<dim, Number, TransferType>::restrict_and_add(
   const unsigned int                                     from_level,
   LinearAlgebra::distributed::BlockVector<Number> &      dst,
   const LinearAlgebra::distributed::BlockVector<Number> &src) const
@@ -863,18 +875,6 @@ MGTransferBlockMatrixFree<dim, Number>::restrict_and_add(
 
 
 
-template <int dim, typename Number>
-std::size_t
-MGTransferBlockMatrixFree<dim, Number>::memory_consumption() const
-{
-  std::size_t total_memory_consumption = 0;
-  for (const auto &el : matrix_free_transfer_vector)
-    total_memory_consumption += el.memory_consumption();
-  return total_memory_consumption;
-}
-
-
-
 // explicit instantiation
 #include "mg_transfer_matrix_free.inst"
 
index 2967b4c71a8f4adb46ca2a853c4ee6ff3514eaf0..868b82d7998e2a488b0901088fa6ce89a6f57fab 100644 (file)
@@ -18,5 +18,9 @@
 for (deal_II_dimension : DIMENSIONS; S1 : REAL_SCALARS)
   {
     template class MGTransferMatrixFree<deal_II_dimension, S1>;
+    template class MGTransferBlockMatrixFreeBase<
+      deal_II_dimension,
+      S1,
+      MGTransferMatrixFree<deal_II_dimension, S1>>;
     template class MGTransferBlockMatrixFree<deal_II_dimension, S1>;
   }

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.