]> https://gitweb.dealii.org/ - dealii.git/commitdiff
MGTransferMatrixFree: remove unused template argument 15750/head
authorPeter Munch <peterrmuench@gmail.com>
Sat, 15 Jul 2023 10:02:02 +0000 (12:02 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Fri, 21 Jul 2023 09:41:27 +0000 (11:41 +0200)
doc/news/changes/incompatibilities/20230721Munch [new file with mode: 0644]
include/deal.II/multigrid/mg_transfer_global_coarsening.h
include/deal.II/multigrid/mg_transfer_matrix_free.h
source/multigrid/mg_transfer_matrix_free.cc

diff --git a/doc/news/changes/incompatibilities/20230721Munch b/doc/news/changes/incompatibilities/20230721Munch
new file mode 100644 (file)
index 0000000..31dde10
--- /dev/null
@@ -0,0 +1,8 @@
+Removed: The methods copy_to_mg, copy_from_mg,
+and interpolate_to_mg used of
+MGTransferGlobalCoarsening and MGTransferMatrixFree
+used to have spacedim as template argument. Since it
+is not used and the underlying functionality did not work,
+it has been removed.
+<br>
+(Peter Munch, 2023/07/21)
index 3d01d1420c78cfa7f1e2d2c831b3560382728152..0ef32204b449059acbd6a7863b6ad91e0bc535b8 100644 (file)
@@ -889,11 +889,11 @@ public:
    *
    * @note DoFHandler is not needed here, but is required by the interface.
    */
-  template <class InVector, int spacedim>
+  template <class InVector>
   void
-  copy_to_mg(const DoFHandler<dim, spacedim> &dof_handler,
-             MGLevelObject<VectorType> &      dst,
-             const InVector &                 src) const;
+  copy_to_mg(const DoFHandler<dim> &    dof_handler,
+             MGLevelObject<VectorType> &dst,
+             const InVector &           src) const;
 
   /**
    * Initialize internal vectors and copy the values on the finest
@@ -901,9 +901,9 @@ public:
    *
    * @note DoFHandler is not needed here, but is required by the interface.
    */
-  template <class OutVector, int spacedim>
+  template <class OutVector>
   void
-  copy_from_mg(const DoFHandler<dim, spacedim> &dof_handler,
+  copy_from_mg(const DoFHandler<dim> &          dof_handler,
                OutVector &                      dst,
                const MGLevelObject<VectorType> &src) const;
 
@@ -929,11 +929,11 @@ public:
    * is required to be able to use MGTransferGlobalCoarsening and
    * MGTransferMatrixFree as template argument.
    */
-  template <class InVector, int spacedim>
+  template <class InVector>
   void
-  interpolate_to_mg(const DoFHandler<dim, spacedim> &dof_handler,
-                    MGLevelObject<VectorType> &      dst,
-                    const InVector &                 src) const;
+  interpolate_to_mg(const DoFHandler<dim> &    dof_handler,
+                    MGLevelObject<VectorType> &dst,
+                    const InVector &           src) const;
 
   /**
    * Return the memory consumption of the allocated memory in this class.
@@ -1150,12 +1150,12 @@ MGTransferGlobalCoarsening<dim, VectorType>::restrict_and_add(
 
 
 template <int dim, typename VectorType>
-template <class InVector, int spacedim>
+template <class InVector>
 void
 MGTransferGlobalCoarsening<dim, VectorType>::copy_to_mg(
-  const DoFHandler<dim, spacedim> &dof_handler,
-  MGLevelObject<VectorType> &      dst,
-  const InVector &                 src) const
+  const DoFHandler<dim> &    dof_handler,
+  MGLevelObject<VectorType> &dst,
+  const InVector &           src) const
 {
   (void)dof_handler;
 
@@ -1173,10 +1173,10 @@ MGTransferGlobalCoarsening<dim, VectorType>::copy_to_mg(
 
 
 template <int dim, typename VectorType>
-template <class OutVector, int spacedim>
+template <class OutVector>
 void
 MGTransferGlobalCoarsening<dim, VectorType>::copy_from_mg(
-  const DoFHandler<dim, spacedim> &dof_handler,
+  const DoFHandler<dim> &          dof_handler,
   OutVector &                      dst,
   const MGLevelObject<VectorType> &src) const
 {
@@ -1212,12 +1212,12 @@ MGTransferGlobalCoarsening<dim, VectorType>::interpolate_to_mg(
 
 
 template <int dim, typename VectorType>
-template <class InVector, int spacedim>
+template <class InVector>
 void
 MGTransferGlobalCoarsening<dim, VectorType>::interpolate_to_mg(
-  const DoFHandler<dim, spacedim> &dof_handler,
-  MGLevelObject<VectorType> &      dst,
-  const InVector &                 src) const
+  const DoFHandler<dim> &    dof_handler,
+  MGLevelObject<VectorType> &dst,
+  const InVector &           src) const
 {
   (void)dof_handler;
 
index ddf6fa8d3030da0cc2bbdd4ca012d894a5a5e20c..cef72e3e466c66bab40f60b3646837c629701050 100644 (file)
@@ -102,7 +102,7 @@ public:
    * ignored and internal variants are used instead.
    */
   void
-  build(const DoFHandler<dim, dim> &dof_handler,
+  build(const DoFHandler<dim> &dof_handler,
         const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
           &external_partitioners =
             std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>());
@@ -112,7 +112,7 @@ public:
    * partitioners.
    */
   void
-  build(const DoFHandler<dim, dim> &dof_handler,
+  build(const DoFHandler<dim> &dof_handler,
         const std::function<void(const unsigned int,
                                  LinearAlgebra::distributed::Vector<Number> &)>
           &initialize_dof_vector);
@@ -181,10 +181,10 @@ public:
    *
    * The use of this function is demonstrated in step-66.
    */
-  template <typename BlockVectorType2, int spacedim>
+  template <typename BlockVectorType2>
   void
   interpolate_to_mg(
-    const DoFHandler<dim, spacedim> &                          dof_handler,
+    const DoFHandler<dim> &                                    dof_handler,
     MGLevelObject<LinearAlgebra::distributed::Vector<Number>> &dst,
     const BlockVectorType2 &                                   src) const;
 
@@ -387,42 +387,42 @@ public:
    * This function will initialize @p dst accordingly if needed as required by
    * the Multigrid class.
    */
-  template <typename BlockVectorType2, int spacedim>
+  template <typename BlockVectorType2>
   void
   copy_to_mg(
-    const DoFHandler<dim, spacedim> &                               dof_handler,
+    const DoFHandler<dim> &                                         dof_handler,
     MGLevelObject<LinearAlgebra::distributed::BlockVector<Number>> &dst,
     const BlockVectorType2 &                                        src) const;
 
   /**
    * Same as above for the case that each block has its own DoFHandler.
    */
-  template <typename BlockVectorType2, int spacedim>
+  template <typename BlockVectorType2>
   void
   copy_to_mg(
-    const std::vector<const DoFHandler<dim, spacedim> *> &          dof_handler,
+    const std::vector<const DoFHandler<dim> *> &                    dof_handler,
     MGLevelObject<LinearAlgebra::distributed::BlockVector<Number>> &dst,
     const BlockVectorType2 &                                        src) const;
 
   /**
    * Transfer from multi-level block-vector to normal vector.
    */
-  template <typename BlockVectorType2, int spacedim>
+  template <typename BlockVectorType2>
   void
   copy_from_mg(
-    const DoFHandler<dim, spacedim> &dof_handler,
-    BlockVectorType2 &               dst,
+    const DoFHandler<dim> &dof_handler,
+    BlockVectorType2 &     dst,
     const MGLevelObject<LinearAlgebra::distributed::BlockVector<Number>> &src)
     const;
 
   /**
    * Same as above for the case that each block has its own DoFHandler.
    */
-  template <typename BlockVectorType2, int spacedim>
+  template <typename BlockVectorType2>
   void
   copy_from_mg(
-    const std::vector<const DoFHandler<dim, spacedim> *> &dof_handler,
-    BlockVectorType2 &                                    dst,
+    const std::vector<const DoFHandler<dim> *> &dof_handler,
+    BlockVectorType2 &                          dst,
     const MGLevelObject<LinearAlgebra::distributed::BlockVector<Number>> &src)
     const;
 
@@ -515,13 +515,13 @@ public:
    * Actually build the information for the prolongation for each level.
    */
   void
-  build(const DoFHandler<dim, dim> &dof_handler);
+  build(const DoFHandler<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);
+  build(const std::vector<const DoFHandler<dim> *> &dof_handler);
 
   /**
    * Memory used by this object.
@@ -549,10 +549,10 @@ private:
 
 
 template <int dim, typename Number>
-template <typename BlockVectorType2, int spacedim>
+template <typename BlockVectorType2>
 void
 MGTransferMatrixFree<dim, Number>::interpolate_to_mg(
-  const DoFHandler<dim, spacedim> &                          dof_handler,
+  const DoFHandler<dim> &                                    dof_handler,
   MGLevelObject<LinearAlgebra::distributed::Vector<Number>> &dst,
   const BlockVectorType2 &                                   src) const
 {
@@ -563,7 +563,7 @@ MGTransferMatrixFree<dim, Number>::interpolate_to_mg(
          ExcDimensionMismatch(
            max_level, dof_handler.get_triangulation().n_global_levels() - 1));
 
-  const FiniteElement<dim, spacedim> &fe = dof_handler.get_fe();
+  const auto &fe = dof_handler.get_fe();
 
   for (unsigned int level = min_level; level <= max_level; ++level)
     if (dst[level].size() != dof_handler.n_dofs(level) ||
@@ -642,10 +642,10 @@ MGTransferMatrixFree<dim, Number>::interpolate_to_mg(
 
 
 template <int dim, typename Number, typename TransferType>
-template <typename BlockVectorType2, int spacedim>
+template <typename BlockVectorType2>
 void
 MGTransferBlockMatrixFreeBase<dim, Number, TransferType>::copy_to_mg(
-  const DoFHandler<dim, spacedim> &                               dof_handler,
+  const DoFHandler<dim> &                                         dof_handler,
   MGLevelObject<LinearAlgebra::distributed::BlockVector<Number>> &dst,
   const BlockVectorType2 &                                        src) const
 {
@@ -654,8 +654,8 @@ MGTransferBlockMatrixFreeBase<dim, Number, TransferType>::copy_to_mg(
            "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!"));
-  const std::vector<const DoFHandler<dim, spacedim> *> mg_dofs(src.n_blocks(),
-                                                               &dof_handler);
+  const std::vector<const DoFHandler<dim> *> mg_dofs(src.n_blocks(),
+                                                     &dof_handler);
 
   copy_to_mg(mg_dofs, dst, src);
 }
@@ -663,10 +663,10 @@ MGTransferBlockMatrixFreeBase<dim, Number, TransferType>::copy_to_mg(
 
 
 template <int dim, typename Number, typename TransferType>
-template <typename BlockVectorType2, int spacedim>
+template <typename BlockVectorType2>
 void
 MGTransferBlockMatrixFreeBase<dim, Number, TransferType>::copy_to_mg(
-  const std::vector<const DoFHandler<dim, spacedim> *> &          dof_handler,
+  const std::vector<const DoFHandler<dim> *> &                    dof_handler,
   MGLevelObject<LinearAlgebra::distributed::BlockVector<Number>> &dst,
   const BlockVectorType2 &                                        src) const
 {
@@ -702,11 +702,11 @@ MGTransferBlockMatrixFreeBase<dim, Number, TransferType>::copy_to_mg(
 }
 
 template <int dim, typename Number, typename TransferType>
-template <typename BlockVectorType2, int spacedim>
+template <typename BlockVectorType2>
 void
 MGTransferBlockMatrixFreeBase<dim, Number, TransferType>::copy_from_mg(
-  const DoFHandler<dim, spacedim> &dof_handler,
-  BlockVectorType2 &               dst,
+  const DoFHandler<dim> &dof_handler,
+  BlockVectorType2 &     dst,
   const MGLevelObject<LinearAlgebra::distributed::BlockVector<Number>> &src)
   const
 {
@@ -715,18 +715,18 @@ MGTransferBlockMatrixFreeBase<dim, Number, TransferType>::copy_from_mg(
            "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!"));
-  const std::vector<const DoFHandler<dim, spacedim> *> mg_dofs(dst.n_blocks(),
-                                                               &dof_handler);
+  const std::vector<const DoFHandler<dim> *> mg_dofs(dst.n_blocks(),
+                                                     &dof_handler);
 
   copy_from_mg(mg_dofs, dst, src);
 }
 
 template <int dim, typename Number, typename TransferType>
-template <typename BlockVectorType2, int spacedim>
+template <typename BlockVectorType2>
 void
 MGTransferBlockMatrixFreeBase<dim, Number, TransferType>::copy_from_mg(
-  const std::vector<const DoFHandler<dim, spacedim> *> &dof_handler,
-  BlockVectorType2 &                                    dst,
+  const std::vector<const DoFHandler<dim> *> &dof_handler,
+  BlockVectorType2 &                          dst,
   const MGLevelObject<LinearAlgebra::distributed::BlockVector<Number>> &src)
   const
 {
index da32097f733bb4cf958c736f5ee7cbbba2163c92..acdaed6a18c760e7a7476a9acccf89e3c216a013 100644 (file)
@@ -97,7 +97,7 @@ MGTransferMatrixFree<dim, Number>::clear()
 template <int dim, typename Number>
 void
 MGTransferMatrixFree<dim, Number>::build(
-  const DoFHandler<dim, dim> &dof_handler,
+  const DoFHandler<dim> &dof_handler,
   const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
     &external_partitioners)
 {
@@ -203,7 +203,7 @@ MGTransferMatrixFree<dim, Number>::build(
 template <int dim, typename Number>
 void
 MGTransferMatrixFree<dim, Number>::build(
-  const DoFHandler<dim, dim> &dof_handler,
+  const DoFHandler<dim> &dof_handler,
   const std::function<void(const unsigned int,
                            LinearAlgebra::distributed::Vector<Number> &)>
     &initialize_dof_vector)
@@ -756,7 +756,7 @@ MGTransferBlockMatrixFree<dim, Number>::initialize_constraints(
 template <int dim, typename Number>
 void
 MGTransferBlockMatrixFree<dim, Number>::build(
-  const DoFHandler<dim, dim> &dof_handler)
+  const DoFHandler<dim> &dof_handler)
 {
   AssertDimension(this->matrix_free_transfer_vector.size(), 1);
   this->matrix_free_transfer_vector[0].build(dof_handler);
@@ -767,7 +767,7 @@ MGTransferBlockMatrixFree<dim, Number>::build(
 template <int dim, typename Number>
 void
 MGTransferBlockMatrixFree<dim, Number>::build(
-  const std::vector<const DoFHandler<dim, dim> *> &dof_handler)
+  const std::vector<const DoFHandler<dim> *> &dof_handler)
 {
   AssertDimension(this->matrix_free_transfer_vector.size(), dof_handler.size());
   for (unsigned int i = 0; i < dof_handler.size(); ++i)

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.