]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Make a specialization the general template. 16723/head
authorWolfgang Bangerth <bangerth@colostate.edu>
Wed, 6 Mar 2024 23:45:55 +0000 (16:45 -0700)
committerWolfgang Bangerth <bangerth@colostate.edu>
Wed, 6 Mar 2024 23:45:55 +0000 (16:45 -0700)
We can do that because the general template does not work anyway, and can
not have been used anywhere.

include/deal.II/multigrid/mg_transfer_global_coarsening.h
include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h

index dc84d2a4a51cc397a2128af2236be820b4a3d3fe..f84b2ea795eaafcdfa4b58842a0f9da5df4587f2 100644 (file)
@@ -207,66 +207,27 @@ namespace MGTransferGlobalCoarseningTools
 
 
 /**
- * Abstract base class for transfer operators between two multigrid levels.
+ * An abstract base class for transfer operators between two multigrid levels.
+ * Specialization for LinearAlgebra::distributed::Vector. The implementation of
+ * restriction and prolongation between levels is delegated to derived classes,
+ * which implement prolongate_and_add_internal() and restrict_and_add_internal()
+ * accordingly.
  */
 template <typename VectorType>
 class MGTwoLevelTransferBase : public Subscriptor
 {
 public:
-  /**
-   * Perform prolongation on a solution vector.
-   */
-  virtual void
-  prolongate_and_add(VectorType &dst, const VectorType &src) const = 0;
-
-  /**
-   * Perform restriction on a residual vector.
-   */
-  virtual void
-  restrict_and_add(VectorType &dst, const VectorType &src) const = 0;
-
-  /**
-   * Perform interpolation of a solution vector from the fine level to the
-   * coarse level. This function is different from restriction, where a
-   * weighted residual is transferred to a coarser level (transposition of
-   * prolongation matrix). In other words, restriction acts on right hand
-   * side vectors, whereas interpolation acts on solution vectors.
-   */
-  virtual void
-  interpolate(VectorType &dst, const VectorType &src) const = 0;
-
-  /**
-   * Enable inplace vector operations if external and internal vectors
-   * are compatible.
-   */
-  virtual void
-  enable_inplace_operations_if_possible(
-    const std::shared_ptr<const Utilities::MPI::Partitioner>
-      &partitioner_coarse,
-    const std::shared_ptr<const Utilities::MPI::Partitioner>
-      &partitioner_fine) = 0;
+  static_assert(
+    std::is_same_v<
+      VectorType,
+      LinearAlgebra::distributed::Vector<typename VectorType::value_type>>,
+    "This class is currently only implemented for vectors of "
+    "type LinearAlgebra::distributed::Vector.");
 
   /**
-   * Return the memory consumption of the allocated memory in this class.
+   * The scalar type used by the vector-type template argument.
    */
-  virtual std::size_t
-  memory_consumption() const = 0;
-};
-
-
-/**
- * Base class for transfer operators between two multigrid levels.
- * Specialization for LinearAlgebra::distributed::Vector. The implementation of
- * restriction and prolongation between levels is delegated to derived classes,
- * which implement prolongate_and_add_internal() and restrict_and_add_internal()
- * accordingly.
- */
-template <typename Number>
-class MGTwoLevelTransferBase<LinearAlgebra::distributed::Vector<Number>>
-  : public Subscriptor
-{
-public:
-  using VectorType = LinearAlgebra::distributed::Vector<Number>;
+  using Number = typename VectorType::value_type;
 
   /**
    * Default constructor.
@@ -335,8 +296,7 @@ protected:
    * partitioners.
    */
   void
-  update_ghost_values(
-    const LinearAlgebra::distributed::Vector<Number> &vec) const;
+  update_ghost_values(const VectorType &vec) const;
 
   /**
    * A wrapper around compress() optimized in case the
@@ -344,8 +304,7 @@ protected:
    * partitioners.
    */
   void
-  compress(LinearAlgebra::distributed::Vector<Number> &vec,
-           const VectorOperation::values               op) const;
+  compress(VectorType &vec, const VectorOperation::values op) const;
 
   /**
    * A wrapper around zero_out_ghost_values() optimized in case the
@@ -353,8 +312,7 @@ protected:
    * partitioners.
    */
   void
-  zero_out_ghost_values(
-    const LinearAlgebra::distributed::Vector<Number> &vec) const;
+  zero_out_ghost_values(const VectorType &vec) const;
 
   /**
    * Enable inplace vector operations if external and internal vectors
@@ -447,57 +405,24 @@ template <int dim, typename VectorType>
 class MGTwoLevelTransfer : public MGTwoLevelTransferBase<VectorType>
 {
 public:
-  /**
-   * @copydoc MGTwoLevelTransferBase::prolongate_and_add
-   */
-  void
-  prolongate_and_add(VectorType &dst, const VectorType &src) const override;
+  static_assert(
+    std::is_same_v<
+      VectorType,
+      LinearAlgebra::distributed::Vector<typename VectorType::value_type>>,
+    "This class is currently only implemented for vectors of "
+    "type LinearAlgebra::distributed::Vector.");
 
   /**
-   * @copydoc MGTwoLevelTransferBase::restrict_and_add
+   * The scalar type used by the vector-type template argument.
    */
-  void
-  restrict_and_add(VectorType &dst, const VectorType &src) const override;
+  using Number = typename VectorType::value_type;
 
   /**
-   * @copydoc MGTwoLevelTransferBase::interpolate
-   */
-  void
-  interpolate(VectorType &dst, const VectorType &src) const override;
-
-  /**
-   * Enable inplace vector operations if external and internal vectors
-   * are compatible.
+   * A data type representing a vectorized array of the same kind of objects
+   * stored in the `VectorType`.
    */
-  void
-  enable_inplace_operations_if_possible(
-    const std::shared_ptr<const Utilities::MPI::Partitioner>
-      &partitioner_coarse,
-    const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner_fine)
-    override;
-
-  /**
-   * Return the memory consumption of the allocated memory in this class.
-   */
-  std::size_t
-  memory_consumption() const override;
-};
-
-
-
-/**
- * Class for transfer between two multigrid levels for p- or global coarsening.
- * Specialization for LinearAlgebra::distributed::Vector.
- *
- * The implementation of this class is explained in detail in @cite munch2022gc.
- */
-template <int dim, typename Number>
-class MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>
-  : public MGTwoLevelTransferBase<LinearAlgebra::distributed::Vector<Number>>
-{
   using VectorizedArrayType = VectorizedArray<Number>;
 
-public:
   /**
    * Set up global coarsening between the given DoFHandler objects (
    * @p dof_handler_fine and @p dof_handler_coarse). The transfer
@@ -573,9 +498,7 @@ public:
    * @copydoc MGTwoLevelTransferBase::interpolate
    */
   void
-  interpolate(
-    LinearAlgebra::distributed::Vector<Number>       &dst,
-    const LinearAlgebra::distributed::Vector<Number> &src) const override;
+  interpolate(VectorType &dst, const VectorType &src) const override;
 
   /**
    * Enable inplace vector operations if external and internal vectors
@@ -596,14 +519,12 @@ public:
 
 protected:
   void
-  prolongate_and_add_internal(
-    LinearAlgebra::distributed::Vector<Number>       &dst,
-    const LinearAlgebra::distributed::Vector<Number> &src) const override;
+  prolongate_and_add_internal(VectorType       &dst,
+                              const VectorType &src) const override;
 
   void
-  restrict_and_add_internal(
-    LinearAlgebra::distributed::Vector<Number>       &dst,
-    const LinearAlgebra::distributed::Vector<Number> &src) const override;
+  restrict_and_add_internal(VectorType       &dst,
+                            const VectorType &src) const override;
 
 private:
   /**
index f2df2c85f318c43aea01def00cd653fb2703343c..3025bd18a95675371c5fc13413cba20e2bc8db32 100644 (file)
@@ -2852,20 +2852,18 @@ namespace MGTransferGlobalCoarseningTools
 
 
 
-template <typename Number>
-MGTwoLevelTransferBase<
-  LinearAlgebra::distributed::Vector<Number>>::MGTwoLevelTransferBase()
+template <typename VectorType>
+MGTwoLevelTransferBase<VectorType>::MGTwoLevelTransferBase()
   : vec_fine_needs_ghost_update(true)
 {}
 
 
 
-template <typename Number>
+template <typename VectorType>
 void
-MGTwoLevelTransferBase<LinearAlgebra::distributed::Vector<Number>>::
-  prolongate_and_add(
-    LinearAlgebra::distributed::Vector<Number>       &dst,
-    const LinearAlgebra::distributed::Vector<Number> &src) const
+MGTwoLevelTransferBase<VectorType>::prolongate_and_add(
+  VectorType       &dst,
+  const VectorType &src) const
 {
   const bool  use_dst_inplace = this->vec_fine.size() == 0;
   auto *const vec_fine_ptr    = use_dst_inplace ? &dst : &this->vec_fine;
@@ -2903,12 +2901,11 @@ MGTwoLevelTransferBase<LinearAlgebra::distributed::Vector<Number>>::
 
 
 
-template <int dim, typename Number>
+template <int dim, typename VectorType>
 void
-MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
-  prolongate_and_add_internal(
-    LinearAlgebra::distributed::Vector<Number>       &dst,
-    const LinearAlgebra::distributed::Vector<Number> &src) const
+MGTwoLevelTransfer<dim, VectorType>::prolongate_and_add_internal(
+  VectorType       &dst,
+  const VectorType &src) const
 {
   const unsigned int n_lanes = VectorizedArrayType::size();
 
@@ -3033,11 +3030,11 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
 
 
 
-template <typename Number>
+template <typename VectorType>
 void
-MGTwoLevelTransferBase<LinearAlgebra::distributed::Vector<Number>>::
-  restrict_and_add(LinearAlgebra::distributed::Vector<Number>       &dst,
-                   const LinearAlgebra::distributed::Vector<Number> &src) const
+MGTwoLevelTransferBase<VectorType>::restrict_and_add(
+  VectorType       &dst,
+  const VectorType &src) const
 {
   const bool        use_src_inplace = this->vec_fine.size() == 0;
   const auto *const vec_fine_ptr    = use_src_inplace ? &src : &this->vec_fine;
@@ -3084,12 +3081,11 @@ MGTwoLevelTransferBase<LinearAlgebra::distributed::Vector<Number>>::
 
 
 
-template <int dim, typename Number>
+template <int dim, typename VectorType>
 void
-MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
-  restrict_and_add_internal(
-    LinearAlgebra::distributed::Vector<Number>       &dst,
-    const LinearAlgebra::distributed::Vector<Number> &src) const
+MGTwoLevelTransfer<dim, VectorType>::restrict_and_add_internal(
+  VectorType       &dst,
+  const VectorType &src) const
 {
   const unsigned int n_lanes = VectorizedArrayType::size();
 
@@ -3215,11 +3211,10 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
 
 
 
-template <int dim, typename Number>
+template <int dim, typename VectorType>
 void
-MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
-  interpolate(LinearAlgebra::distributed::Vector<Number>       &dst,
-              const LinearAlgebra::distributed::Vector<Number> &src) const
+MGTwoLevelTransfer<dim, VectorType>::interpolate(VectorType       &dst,
+                                                 const VectorType &src) const
 {
   const unsigned int n_lanes = VectorizedArrayType::size();
 
@@ -3397,10 +3392,12 @@ namespace internal
   } // namespace
 } // namespace internal
 
-template <typename Number>
+
+
+template <typename VectorType>
 template <int dim, std::size_t width, typename IndexType>
 void
-MGTwoLevelTransferBase<LinearAlgebra::distributed::Vector<Number>>::
+MGTwoLevelTransferBase<VectorType>::
   internal_enable_inplace_operations_if_possible(
     const std::shared_ptr<const Utilities::MPI::Partitioner>
       &external_partitioner_coarse,
@@ -3465,14 +3462,15 @@ MGTwoLevelTransferBase<LinearAlgebra::distributed::Vector<Number>>::
     }
 }
 
-template <int dim, typename Number>
+
+
+template <int dim, typename VectorType>
 void
-MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
-  enable_inplace_operations_if_possible(
-    const std::shared_ptr<const Utilities::MPI::Partitioner>
-      &external_partitioner_coarse,
-    const std::shared_ptr<const Utilities::MPI::Partitioner>
-      &external_partitioner_fine)
+MGTwoLevelTransfer<dim, VectorType>::enable_inplace_operations_if_possible(
+  const std::shared_ptr<const Utilities::MPI::Partitioner>
+    &external_partitioner_coarse,
+  const std::shared_ptr<const Utilities::MPI::Partitioner>
+    &external_partitioner_fine)
 {
   this->internal_enable_inplace_operations_if_possible(
     external_partitioner_coarse,
@@ -3484,15 +3482,15 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
 
 
 
-template <int dim, typename Number>
+template <int dim, typename VectorType>
 void
-MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
-  reinit_geometric_transfer(const DoFHandler<dim>           &dof_handler_fine,
-                            const DoFHandler<dim>           &dof_handler_coarse,
-                            const AffineConstraints<Number> &constraints_fine,
-                            const AffineConstraints<Number> &constraints_coarse,
-                            const unsigned int               mg_level_fine,
-                            const unsigned int               mg_level_coarse)
+MGTwoLevelTransfer<dim, VectorType>::reinit_geometric_transfer(
+  const DoFHandler<dim>           &dof_handler_fine,
+  const DoFHandler<dim>           &dof_handler_coarse,
+  const AffineConstraints<Number> &constraints_fine,
+  const AffineConstraints<Number> &constraints_coarse,
+  const unsigned int               mg_level_fine,
+  const unsigned int               mg_level_coarse)
 {
   internal::MGTwoLevelTransferImplementation::reinit_geometric_transfer(
     dof_handler_fine,
@@ -3506,16 +3504,15 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
 
 
 
-template <int dim, typename Number>
+template <int dim, typename VectorType>
 void
-MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
-  reinit_polynomial_transfer(
-    const DoFHandler<dim>           &dof_handler_fine,
-    const DoFHandler<dim>           &dof_handler_coarse,
-    const AffineConstraints<Number> &constraints_fine,
-    const AffineConstraints<Number> &constraints_coarse,
-    const unsigned int               mg_level_fine,
-    const unsigned int               mg_level_coarse)
+MGTwoLevelTransfer<dim, VectorType>::reinit_polynomial_transfer(
+  const DoFHandler<dim>           &dof_handler_fine,
+  const DoFHandler<dim>           &dof_handler_coarse,
+  const AffineConstraints<Number> &constraints_fine,
+  const AffineConstraints<Number> &constraints_coarse,
+  const unsigned int               mg_level_fine,
+  const unsigned int               mg_level_coarse)
 {
   internal::MGTwoLevelTransferImplementation::reinit_polynomial_transfer(
     dof_handler_fine,
@@ -3529,9 +3526,9 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
 
 
 
-template <int dim, typename Number>
+template <int dim, typename VectorType>
 void
-MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::reinit(
+MGTwoLevelTransfer<dim, VectorType>::reinit(
   const DoFHandler<dim>           &dof_handler_fine,
   const DoFHandler<dim>           &dof_handler_coarse,
   const AffineConstraints<Number> &constraints_fine,
@@ -3622,11 +3619,11 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::reinit(
 
 
 
-template <int dim, typename Number>
+template <int dim, typename VectorType>
 bool
-MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
-  fast_polynomial_transfer_supported(const unsigned int fe_degree_fine,
-                                     const unsigned int fe_degree_coarse)
+MGTwoLevelTransfer<dim, VectorType>::fast_polynomial_transfer_supported(
+  const unsigned int fe_degree_fine,
+  const unsigned int fe_degree_coarse)
 {
   CellTransferFactory cell_transfer(fe_degree_fine, fe_degree_coarse);
   CellProlongatorTest cell_transfer_test;
@@ -3636,10 +3633,9 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
 
 
 
-template <int dim, typename Number>
+template <int dim, typename VectorType>
 std::size_t
-MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
-  memory_consumption() const
+MGTwoLevelTransfer<dim, VectorType>::memory_consumption() const
 {
   std::size_t size = 0;
 
@@ -3664,11 +3660,10 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
 
 
 
-template <typename Number>
+template <typename VectorType>
 void
-MGTwoLevelTransferBase<LinearAlgebra::distributed::Vector<Number>>::
-  update_ghost_values(
-    const LinearAlgebra::distributed::Vector<Number> &vec) const
+MGTwoLevelTransferBase<VectorType>::update_ghost_values(
+  const VectorType &vec) const
 {
   if ((vec.get_partitioner().get() == this->partitioner_coarse.get()) &&
       (this->partitioner_coarse_embedded != nullptr))
@@ -3686,11 +3681,11 @@ MGTwoLevelTransferBase<LinearAlgebra::distributed::Vector<Number>>::
 
 
 
-template <typename Number>
+template <typename VectorType>
 void
-MGTwoLevelTransferBase<LinearAlgebra::distributed::Vector<Number>>::compress(
-  LinearAlgebra::distributed::Vector<Number> &vec,
-  const VectorOperation::values               op) const
+MGTwoLevelTransferBase<VectorType>::compress(
+  VectorType                   &vec,
+  const VectorOperation::values op) const
 {
   Assert(op == VectorOperation::add, ExcNotImplemented());
 
@@ -3710,11 +3705,10 @@ MGTwoLevelTransferBase<LinearAlgebra::distributed::Vector<Number>>::compress(
 
 
 
-template <typename Number>
+template <typename VectorType>
 void
-MGTwoLevelTransferBase<LinearAlgebra::distributed::Vector<Number>>::
-  zero_out_ghost_values(
-    const LinearAlgebra::distributed::Vector<Number> &vec) const
+MGTwoLevelTransferBase<VectorType>::zero_out_ghost_values(
+  const VectorType &vec) const
 {
   if ((vec.get_partitioner().get() == this->partitioner_coarse.get()) &&
       (this->partitioner_coarse_embedded != nullptr))

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.