]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Introduce MGTwoLevelTransfer::reinit() 12484/head
authorPeter Munch <peterrmuench@gmail.com>
Tue, 22 Jun 2021 07:40:05 +0000 (09:40 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Tue, 29 Jun 2021 08:45:34 +0000 (10:45 +0200)
15 files changed:
doc/news/changes/minor/20210623Munch-1 [new file with mode: 0644]
examples/step-75/step-75.cc
include/deal.II/multigrid/mg_transfer_global_coarsening.h
include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h
tests/multigrid-global-coarsening/mg_transfer_a_01.cc
tests/multigrid-global-coarsening/mg_transfer_a_02.cc
tests/multigrid-global-coarsening/mg_transfer_a_03.cc
tests/multigrid-global-coarsening/mg_transfer_p_01.cc
tests/multigrid-global-coarsening/mg_transfer_p_02.cc
tests/multigrid-global-coarsening/mg_transfer_p_03.cc
tests/multigrid-global-coarsening/mg_transfer_p_04.cc
tests/multigrid-global-coarsening/mg_transfer_p_05.cc
tests/multigrid-global-coarsening/multigrid_a_01.cc
tests/multigrid-global-coarsening/multigrid_p_01.cc
tests/multigrid-global-coarsening/multigrid_p_02.cc

diff --git a/doc/news/changes/minor/20210623Munch-1 b/doc/news/changes/minor/20210623Munch-1
new file mode 100644 (file)
index 0000000..2acd580
--- /dev/null
@@ -0,0 +1,5 @@
+New: The new unified function MGTwoLevelTransfer::reinit() selects automatically
+if MGTwoLevelTransfer::reinit_geometric_transfer() or 
+MGTwoLevelTransfer::reinit_polynomial_transfer() is needed.
+<br>
+(Peter Munch, 2021/06/23)
index b946d3455a093d54466ce19e50898460260fb13a..951d173b0139cd542909835ff727b26186da94ad 100644 (file)
@@ -700,9 +700,8 @@ namespace Step75
         fe_index_for_degree[degree] = i;
       }
 
-    unsigned int minlevel   = 0;
-    unsigned int minlevel_p = n_h_levels;
-    unsigned int maxlevel   = n_h_levels + n_p_levels - 1;
+    unsigned int minlevel = 0;
+    unsigned int maxlevel = n_h_levels + n_p_levels - 1;
 
     dof_handlers.resize(minlevel, maxlevel);
     operators.resize(minlevel, maxlevel);
@@ -802,17 +801,11 @@ namespace Step75
 
     // Set up intergrid operators and collect transfer operators within a single
     // operator as needed by the Multigrid solver class.
-    for (unsigned int level = minlevel; level < minlevel_p; ++level)
-      transfers[level + 1].reinit_geometric_transfer(dof_handlers[level + 1],
-                                                     dof_handlers[level],
-                                                     constraints[level + 1],
-                                                     constraints[level]);
-
-    for (unsigned int level = minlevel_p; level < maxlevel; ++level)
-      transfers[level + 1].reinit_polynomial_transfer(dof_handlers[level + 1],
-                                                      dof_handlers[level],
-                                                      constraints[level + 1],
-                                                      constraints[level]);
+    for (unsigned int level = minlevel; level < maxlevel; ++level)
+      transfers[level + 1].reinit(dof_handlers[level + 1],
+                                  dof_handlers[level],
+                                  constraints[level + 1],
+                                  constraints[level]);
 
     MGTransferGlobalCoarsening<dim, VectorType> transfer(
       transfers, [&](const auto l, auto &vec) {
index 9187abc51589f0668bc27908a857ddc5166527c2..38790e66e5835339ecece076dbc93acd248141bd 100644 (file)
@@ -216,6 +216,27 @@ public:
     const unsigned int mg_level_fine   = numbers::invalid_unsigned_int,
     const unsigned int mg_level_coarse = numbers::invalid_unsigned_int);
 
+  /**
+   * Set up transfer operator between the given DoFHandler objects (
+   * @p dof_handler_fine and @p dof_handler_coarse). Depending on the
+   * underlying Triangulation objects polynomial or geometrical global
+   * coarsening is performed.
+   *
+   * @note While geometric transfer can be only performed on active levels
+   *   (`numbers::invalid_unsigned_int`), polynomial transfers can also be
+   *   performed on coarse-grid levels.
+   *
+   * @note The function polynomial_transfer_supported() can be used to
+   *   check if the given polynomial coarsening strategy is supported.
+   */
+  void
+  reinit(const DoFHandler<dim> &          dof_handler_fine,
+         const DoFHandler<dim> &          dof_handler_coarse,
+         const AffineConstraints<Number> &constraint_fine,
+         const AffineConstraints<Number> &constraint_coarse,
+         const unsigned int mg_level_fine   = numbers::invalid_unsigned_int,
+         const unsigned int mg_level_coarse = numbers::invalid_unsigned_int);
+
   /**
    * Check if a fast templated version of the polynomial transfer between
    * @p fe_degree_fine and @p fe_degree_coarse is available.
index 83f55a237453f9bd637e65916f067033c4a55137..3e1e5e5e4fa8c51fe8747763cb12311e1b18e25c 100644 (file)
@@ -411,6 +411,8 @@ namespace internal
     const std::function<unsigned int()> active_fe_index_function;
   };
 
+
+
   template <int dim>
   class FineDoFHandlerView
   {
@@ -2804,6 +2806,96 @@ MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
 
 
 
+template <int dim, typename Number>
+void
+MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::reinit(
+  const DoFHandler<dim> &          dof_handler_fine,
+  const DoFHandler<dim> &          dof_handler_coarse,
+  const AffineConstraints<Number> &constraint_fine,
+  const AffineConstraints<Number> &constraint_coarse,
+  const unsigned int               mg_level_fine,
+  const unsigned int               mg_level_coarse)
+{
+  // determine if polynomial transfer can be performed via the following two
+  // criteria:
+  // 1) multigrid levels can be only used with polynomial transfer
+  bool do_polynomial_transfer =
+    (mg_level_fine != numbers::invalid_unsigned_int) ||
+    (mg_level_coarse != numbers::invalid_unsigned_int);
+
+  // 2) the meshes are identical
+  if (do_polynomial_transfer == false)
+    {
+      const internal::CellIDTranslator<dim> cell_id_translator(
+        dof_handler_fine.get_triangulation().n_global_coarse_cells(),
+        dof_handler_fine.get_triangulation().n_global_levels());
+
+      AssertDimension(
+        dof_handler_fine.get_triangulation().n_global_coarse_cells(),
+        dof_handler_coarse.get_triangulation().n_global_coarse_cells());
+      AssertIndexRange(dof_handler_coarse.get_triangulation().n_global_levels(),
+                       dof_handler_fine.get_triangulation().n_global_levels() +
+                         1);
+
+      IndexSet is_locally_owned_fine(cell_id_translator.size());
+      IndexSet is_locally_owned_coarse(cell_id_translator.size());
+
+      for (const auto &cell : dof_handler_fine.active_cell_iterators())
+        if (!cell->is_artificial() && cell->is_locally_owned())
+          is_locally_owned_fine.add_index(cell_id_translator.translate(cell));
+
+      for (const auto &cell : dof_handler_coarse.active_cell_iterators())
+        if (!cell->is_artificial() && cell->is_locally_owned())
+          is_locally_owned_coarse.add_index(cell_id_translator.translate(cell));
+
+      const MPI_Comm communicator = dof_handler_fine.get_communicator();
+
+      std::vector<unsigned int> owning_ranks(
+        is_locally_owned_coarse.n_elements());
+
+      Utilities::MPI::internal::ComputeIndexOwner::ConsensusAlgorithmsPayload
+        process(is_locally_owned_fine,
+                is_locally_owned_coarse,
+                communicator,
+                owning_ranks,
+                false);
+
+      Utilities::MPI::ConsensusAlgorithms::Selector<
+        std::pair<types::global_cell_index, types::global_cell_index>,
+        unsigned int>
+        consensus_algorithm(process, communicator);
+      consensus_algorithm.run();
+
+      bool all_cells_found = true;
+
+      for (unsigned i = 0; i < is_locally_owned_coarse.n_elements(); ++i)
+        all_cells_found &= (owning_ranks[i] != numbers::invalid_unsigned_int);
+
+      do_polynomial_transfer =
+        Utilities::MPI::min(static_cast<unsigned int>(all_cells_found),
+                            communicator) == 1;
+    }
+
+  if (do_polynomial_transfer)
+    internal::MGTwoLevelTransferImplementation::reinit_polynomial_transfer(
+      dof_handler_fine,
+      dof_handler_coarse,
+      constraint_fine,
+      constraint_coarse,
+      mg_level_fine,
+      mg_level_coarse,
+      *this);
+  else
+    internal::MGTwoLevelTransferImplementation::reinit_geometric_transfer(
+      dof_handler_fine,
+      dof_handler_coarse,
+      constraint_fine,
+      constraint_coarse,
+      *this);
+}
+
+
+
 template <int dim, typename Number>
 bool
 MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>::
index 621f8e8b8d5e379778a6421c6c931cef3f257b67..5d0a2a48570910e8a9113e1d42bcb81dbe0b8f44 100644 (file)
@@ -107,10 +107,10 @@ do_test(const FiniteElement<dim> &fe_fine, const FiniteElement<dim> &fe_coarse)
 
   // setup transfer operator
   MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>> transfer;
-  transfer.reinit_geometric_transfer(dof_handler_fine,
-                                     dof_handler_coarse,
-                                     constraint_fine,
-                                     constraint_coarse);
+  transfer.reinit(dof_handler_fine,
+                  dof_handler_coarse,
+                  constraint_fine,
+                  constraint_coarse);
 
   test_transfer_operator(transfer, dof_handler_fine, dof_handler_coarse);
 }
index f9ba6d81af375ae5e27a3dac17b9ef10c26f0c85..7f624d96ab40891092a5defdb75b926a2fe564d1 100644 (file)
@@ -94,10 +94,10 @@ do_test(const FiniteElement<dim> &fe_fine, const FiniteElement<dim> &fe_coarse)
 
   // setup transfer operator
   MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>> transfer;
-  transfer.reinit_geometric_transfer(dof_handler_fine,
-                                     dof_handler_coarse,
-                                     constraint_fine,
-                                     constraint_coarse);
+  transfer.reinit(dof_handler_fine,
+                  dof_handler_coarse,
+                  constraint_fine,
+                  constraint_coarse);
 
   test_transfer_operator(transfer, dof_handler_fine, dof_handler_coarse);
 }
index c9566506837e3a1d064c9193afdc7d9a80347f88..230c2ef27fd6a65024fa80fed10e759c5a92cba7 100644 (file)
@@ -102,10 +102,10 @@ do_test(const FiniteElement<dim> &fe_fine, const FiniteElement<dim> &fe_coarse)
 
   // setup transfer operator
   MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>> transfer;
-  transfer.reinit_geometric_transfer(dof_handler_fine,
-                                     dof_handler_coarse,
-                                     constraint_fine,
-                                     constraint_coarse);
+  transfer.reinit(dof_handler_fine,
+                  dof_handler_coarse,
+                  constraint_fine,
+                  constraint_coarse);
 
   test_transfer_operator(transfer, dof_handler_fine, dof_handler_coarse);
 }
index 91561c532799fb8c68c6891acf55e68276244306..ada0702cfe369f6002bf598da88443aaffe6ffa1 100644 (file)
@@ -92,10 +92,10 @@ do_test(const FiniteElement<dim> &fe_fine, const FiniteElement<dim> &fe_coarse)
 
   // setup transfer operator
   MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>> transfer;
-  transfer.reinit_polynomial_transfer(dof_handler_fine,
-                                      dof_handler_coarse,
-                                      constraint_fine,
-                                      constraint_coarse);
+  transfer.reinit(dof_handler_fine,
+                  dof_handler_coarse,
+                  constraint_fine,
+                  constraint_coarse);
 
   test_transfer_operator(transfer, dof_handler_fine, dof_handler_coarse);
 }
index 0a32b9630c83c5b932643b2e43befebebb6c3ec9..0acd948bfc267fe0d292cf275e17c0bc7f2cce65 100644 (file)
@@ -140,10 +140,10 @@ do_test()
       // setup transfer operator
       MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>
         transfer;
-      transfer.reinit_polynomial_transfer(dof_handler_fine,
-                                          dof_handler_coarse,
-                                          constraint_fine,
-                                          constraint_coarse);
+      transfer.reinit(dof_handler_fine,
+                      dof_handler_coarse,
+                      constraint_fine,
+                      constraint_coarse);
 
       test_transfer_operator(transfer, dof_handler_fine, dof_handler_coarse);
 
index aaa1e21d70be94e8a0608441d68aac5392ebb0a4..c533f890f4a18fde43db62a0061741b69e51cbe0 100644 (file)
@@ -79,10 +79,10 @@ do_test(const FiniteElement<dim> &fe_fine, const FiniteElement<dim> &fe_coarse)
 
   // setup transfer operator
   MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>> transfer;
-  transfer.reinit_polynomial_transfer(dof_handler_fine,
-                                      dof_handler_coarse,
-                                      constraint_fine,
-                                      constraint_coarse);
+  transfer.reinit(dof_handler_fine,
+                  dof_handler_coarse,
+                  constraint_fine,
+                  constraint_coarse);
 
   test_transfer_operator(transfer, dof_handler_fine, dof_handler_coarse);
 }
index 004ee665a9ce94e0e06ca6045048d2c46e84c0d8..1dda8054fd0b36c7889e5c72ab30ed574306655d 100644 (file)
@@ -120,10 +120,10 @@ do_test(const FiniteElement<dim> &fe_fine, const FiniteElement<dim> &fe_coarse)
 
   // setup transfer operator
   MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>> transfer;
-  transfer.reinit_polynomial_transfer(dof_handler_fine,
-                                      dof_handler_coarse,
-                                      constraint_fine,
-                                      constraint_coarse);
+  transfer.reinit(dof_handler_fine,
+                  dof_handler_coarse,
+                  constraint_fine,
+                  constraint_coarse);
 
   test_transfer_operator(transfer, dof_handler_fine, dof_handler_coarse);
 }
index f3a7bf8cfaa49c052cc71a6fa98142ba97acb4cb..60acdc1a7c42a58691d0f593a1d8c88aac0877aa 100644 (file)
@@ -92,12 +92,12 @@ do_test(const FiniteElement<dim> &fe_fine, const FiniteElement<dim> &fe_coarse)
     deallog.push("coarse");
     MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>
       transfer;
-    transfer.reinit_polynomial_transfer(dof_handler_fine,
-                                        dof_handler_coarse,
-                                        constraint_fine,
-                                        constraint_coarse,
-                                        0,
-                                        0);
+    transfer.reinit(dof_handler_fine,
+                    dof_handler_coarse,
+                    constraint_fine,
+                    constraint_coarse,
+                    0,
+                    0);
 
     test_transfer_operator(transfer, dof_handler_fine, dof_handler_coarse);
     deallog.pop();
@@ -106,10 +106,10 @@ do_test(const FiniteElement<dim> &fe_fine, const FiniteElement<dim> &fe_coarse)
     deallog.push("active");
     MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>
       transfer;
-    transfer.reinit_polynomial_transfer(dof_handler_fine,
-                                        dof_handler_coarse,
-                                        constraint_fine,
-                                        constraint_coarse);
+    transfer.reinit(dof_handler_fine,
+                    dof_handler_coarse,
+                    constraint_fine,
+                    constraint_coarse);
 
     test_transfer_operator(transfer, dof_handler_fine, dof_handler_coarse);
     deallog.pop();
index fbdc709a904b89c2d2d04266944cf4468d7d33de..62c04637a3ae5de592b23d70b01ae4f3692532a1 100644 (file)
@@ -95,10 +95,10 @@ test(const unsigned int n_refinements,
 
   // set up transfer operator
   for (unsigned int l = min_level; l < max_level; ++l)
-    transfers[l + 1].reinit_geometric_transfer(dof_handlers[l + 1],
-                                               dof_handlers[l],
-                                               constraints[l + 1],
-                                               constraints[l]);
+    transfers[l + 1].reinit(dof_handlers[l + 1],
+                            dof_handlers[l],
+                            constraints[l + 1],
+                            constraints[l]);
 
   MGTransferGlobalCoarsening<dim, VectorType> transfer(
     transfers,
index fc107374fd9ab1234b7fd7214f1d405f7511e022..044374f6ebb3dbe5d2358b3d3ceb7503362d558d 100644 (file)
@@ -123,10 +123,10 @@ test(const unsigned int n_refinements,
 
   // set up transfer operator
   for (unsigned int l = min_level; l < max_level; ++l)
-    transfers[l + 1].reinit_polynomial_transfer(dof_handlers[l + 1],
-                                                dof_handlers[l],
-                                                constraints[l + 1],
-                                                constraints[l]);
+    transfers[l + 1].reinit(dof_handlers[l + 1],
+                            dof_handlers[l],
+                            constraints[l + 1],
+                            constraints[l]);
 
   MGTransferGlobalCoarsening<dim, VectorType> transfer(
     transfers,
index 2b7a764838b704845f95af56a5a59d0b2f335802..00aa4430edfd0baa5cd24182601d8d3c7755e3b1 100644 (file)
@@ -98,12 +98,12 @@ test(const unsigned int n_refinements, const unsigned int fe_degree_fine)
 
   // set up transfer operator
   for (unsigned int l = min_level; l < max_level; ++l)
-    transfers[l + 1].reinit_polynomial_transfer(dof_handlers[l + 1],
-                                                dof_handlers[l],
-                                                constraints[l + 1],
-                                                constraints[l],
-                                                0 /*level*/,
-                                                0 /*level*/);
+    transfers[l + 1].reinit(dof_handlers[l + 1],
+                            dof_handlers[l],
+                            constraints[l + 1],
+                            constraints[l],
+                            0 /*level*/,
+                            0 /*level*/);
 
   MGTransferGlobalCoarsening<dim, VectorType> transfer(
     transfers,

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.