]> https://gitweb.dealii.org/ - dealii.git/commitdiff
MGTransferMatrixFree: assert that only homogeneous or identity constraints are used 16292/head
authorPeter Munch <peterrmuench@gmail.com>
Fri, 24 Nov 2023 19:52:44 +0000 (20:52 +0100)
committerPeter Munch <peterrmuench@gmail.com>
Fri, 24 Nov 2023 20:07:14 +0000 (21:07 +0100)
include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h
source/multigrid/mg_transfer_internal.cc
tests/multigrid/transfer_matrix_free_13.cc
tests/multigrid/transfer_matrix_free_13.output

index 3011e8ce1ffdf0bc1a9686ab9df87a596f562da4..aa97154a47ed4f3497beb6584c7e40942cb97136 100644 (file)
@@ -1774,6 +1774,9 @@ namespace internal
                (mg_level_coarse + 1 == mg_level_fine),
              ExcNotImplemented());
 
+      AssertDimension(constraints_fine.n_inhomogeneities(), 0);
+      AssertDimension(constraints_coarse.n_inhomogeneities(), 0);
+
       transfer.dof_handler_fine = &dof_handler_fine;
       transfer.mg_level_fine    = mg_level_fine;
 
@@ -2317,6 +2320,9 @@ namespace internal
           "(numbers::invalid_unsigned_int) or on refinement levels without "
           "hanging nodes."));
 
+      AssertDimension(constraints_fine.n_inhomogeneities(), 0);
+      AssertDimension(constraints_coarse.n_inhomogeneities(), 0);
+
       transfer.dof_handler_fine = &dof_handler_fine;
       transfer.mg_level_fine    = mg_level_fine;
 
index 2bbafff46317ad8b1c707b6ffcb3bc95d60e6a0c..584e4846e19c60676fb2e730a046be0d288b74e9 100644 (file)
@@ -637,6 +637,39 @@ namespace internal
       dirichlet_indices.clear();
       weights_on_refined.clear();
 
+#ifdef DEBUG
+      if (mg_constrained_dofs)
+        {
+          const unsigned int n_levels =
+            dof_handler.get_triangulation().n_global_levels();
+
+          for (unsigned int l = 0; l < n_levels; ++l)
+            {
+              const auto &constraints =
+                mg_constrained_dofs->get_user_constraint_matrix(l);
+
+              // no inhomogeneities are supported
+              AssertDimension(constraints.n_inhomogeneities(), 0);
+
+              for (const auto dof : constraints.get_local_lines())
+                {
+                  const auto *entries_ptr =
+                    constraints.get_constraint_entries(dof);
+
+                  if (entries_ptr == nullptr)
+                    continue;
+
+                  // only homogeneous or identity constraints are supported
+                  Assert((entries_ptr->size() == 0) ||
+                           ((entries_ptr->size() == 1) &&
+                            (std::abs((*entries_ptr)[0].second - 1.) <
+                             100 * std::numeric_limits<double>::epsilon())),
+                         ExcNotImplemented());
+                }
+            }
+        }
+#endif
+
       // we collect all child DoFs of a mother cell together. For faster
       // tensorized operations, we align the degrees of freedom
       // lexicographically. We distinguish FE_Q elements and FE_DGQ elements
index 197fa2106b83d7974b8a97193a64f6e3172984fe..ac1955f0264a7fc3d59b3da3b1d66da8cec678c4 100644 (file)
@@ -41,7 +41,7 @@ check()
 
   tr.refine_global(1);
 
-  FE_Q<dim> fe = FE_Q<dim>(1);
+  FE_Q<dim> fe(1);
 
   DoFHandler<dim> mgdof(tr);
   mgdof.distribute_dofs(fe);
@@ -63,7 +63,6 @@ check()
       if (user_constraints.can_store_line(face_dofs[i]))
         {
           user_constraints.constrain_dof_to_zero(face_dofs[i]);
-          user_constraints.set_inhomogeneity(face_dofs[i], 5.0);
         }
     }
   user_constraints.close();
@@ -74,6 +73,7 @@ check()
 
   deallog << "SRC Vector" << std::endl;
   LinearAlgebra::distributed::Vector<double> src_level_0(mgdof.n_dofs(0));
+  src_level_0 = 1.0;
   for (unsigned int i = 0; i < mgdof.n_dofs(0); ++i)
     deallog << src_level_0(i) << ' ';
   deallog << std::endl << std::endl;
index 91b1651826d69dc09aac1788303156bd5078b455..cdb7ec1e21ef55bf7c1604869039c2ac3e54d2ec 100644 (file)
@@ -1,13 +1,13 @@
 
 DEAL::SRC Vector
-DEAL::0.00000 0.00000 0.00000 0.00000 
+DEAL::1.00000 1.00000 1.00000 1.00000 
 DEAL::
 DEAL::DST Vector
-DEAL::5.00000 2.50000 5.00000 2.50000 0.00000 0.00000 5.00000 2.50000 0.00000 
+DEAL::0.00000 0.500000 0.00000 0.500000 1.00000 1.00000 0.00000 0.500000 1.00000 
 DEAL::
 DEAL::SRC Vector
-DEAL::0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 0.00000 
+DEAL::1.00000 1.00000 1.00000 1.00000 1.00000 1.00000 1.00000 1.00000 
 DEAL::
 DEAL::DST Vector
-DEAL::5.00000 2.50000 5.00000 2.50000 5.00000 2.50000 5.00000 2.50000 0.00000 0.00000 0.00000 0.00000 5.00000 2.50000 5.00000 2.50000 0.00000 0.00000 5.00000 2.50000 5.00000 2.50000 0.00000 0.00000 5.00000 2.50000 0.00000 
+DEAL::0.00000 0.500000 0.00000 0.500000 0.00000 0.500000 0.00000 0.500000 1.00000 1.00000 1.00000 1.00000 0.00000 0.500000 0.00000 0.500000 1.00000 1.00000 0.00000 0.500000 0.00000 0.500000 1.00000 1.00000 0.00000 0.500000 1.00000 
 DEAL::

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.