]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix TriaAccessor::is_ghost_on_level for case without MPI 13945/head
authorMartin Kronbichler <martin.kronbichler@uni-a.de>
Thu, 9 Jun 2022 09:26:36 +0000 (11:26 +0200)
committerMartin Kronbichler <martin.kronbichler@uni-a.de>
Thu, 9 Jun 2022 12:50:39 +0000 (14:50 +0200)
include/deal.II/grid/tria_accessor.templates.h
source/grid/tria_accessor.cc
tests/multigrid/is_locally_owned_on_level.cc [new file with mode: 0644]
tests/multigrid/is_locally_owned_on_level.output [new file with mode: 0644]

index 56c4409f82b079afa1b5d220ba7f745b0087cb5c..3f51d0902a2f3822b5035727fca47e094ae68b5b 100644 (file)
@@ -3728,9 +3728,10 @@ CellAccessor<dim, spacedim>::is_locally_owned() const
   // subdomain, so the first condition checks whether we have a serial
   // triangulation, in which case all cells are locally owned. The second
   // condition compares the subdomain id in the parallel case.
-  return (this->tria->locally_owned_subdomain() ==
-            numbers::invalid_subdomain_id ||
-          this->subdomain_id() == this->tria->locally_owned_subdomain());
+  const types::subdomain_id locally_owned_subdomain =
+    this->tria->locally_owned_subdomain();
+  return (locally_owned_subdomain == numbers::invalid_subdomain_id ||
+          this->subdomain_id() == locally_owned_subdomain);
 
 #endif
 }
@@ -3748,9 +3749,10 @@ CellAccessor<dim, spacedim>::is_locally_owned_on_level() const
   // subdomain, so the first condition checks whether we have a serial
   // triangulation, in which case all cells are locally owned. The second
   // condition compares the subdomain id in the parallel case.
-  return (this->tria->locally_owned_subdomain() ==
-            numbers::invalid_subdomain_id ||
-          this->level_subdomain_id() == this->tria->locally_owned_subdomain());
+  const types::subdomain_id locally_owned_subdomain =
+    this->tria->locally_owned_subdomain();
+  return (locally_owned_subdomain == numbers::invalid_subdomain_id ||
+          this->level_subdomain_id() == locally_owned_subdomain);
 
 #endif
 }
@@ -3774,10 +3776,12 @@ CellAccessor<dim, spacedim>::is_ghost() const
   // serial triangulation are locally owned and none is ghosted. The second
   // and third conditions check whether the cell's subdomain is not the
   // locally owned one and not artificial.
-  return (this->tria->locally_owned_subdomain() !=
-            numbers::invalid_subdomain_id &&
-          this->subdomain_id() != this->tria->locally_owned_subdomain() &&
-          this->subdomain_id() != numbers::artificial_subdomain_id);
+  const types::subdomain_id locally_owned_subdomain =
+    this->tria->locally_owned_subdomain();
+  const types::subdomain_id subdomain_id = this->subdomain_id();
+  return (locally_owned_subdomain != numbers::invalid_subdomain_id &&
+          subdomain_id != locally_owned_subdomain &&
+          subdomain_id != numbers::artificial_subdomain_id);
 
 #endif
 }
@@ -3788,17 +3792,19 @@ inline bool
 CellAccessor<dim, spacedim>::is_ghost_on_level() const
 {
 #ifndef DEAL_II_WITH_MPI
-  return true;
+  return false;
 #else
 
   // Serial triangulations report invalid_subdomain_id as their locally owned
   // subdomain, so the first condition checks whether we have a serial
   // triangulation, in which case all cells are locally owned. The second
   // condition compares the subdomain id in the parallel case.
-  return (
-    this->tria->locally_owned_subdomain() == numbers::invalid_subdomain_id ||
-    (this->level_subdomain_id() != this->tria->locally_owned_subdomain() &&
-     this->level_subdomain_id() != numbers::artificial_subdomain_id));
+  const types::subdomain_id locally_owned_subdomain =
+    this->tria->locally_owned_subdomain();
+  const types::subdomain_id subdomain_id = this->level_subdomain_id();
+  return (locally_owned_subdomain != numbers::invalid_subdomain_id &&
+          subdomain_id != locally_owned_subdomain &&
+          subdomain_id != numbers::artificial_subdomain_id);
 
 #endif
 }
@@ -3834,7 +3840,9 @@ CellAccessor<dim, spacedim>::is_artificial_on_level() const
 #ifndef DEAL_II_WITH_MPI
   return false;
 #else
-  return (is_locally_owned_on_level() || is_ghost_on_level()) == false;
+  return (this->tria->locally_owned_subdomain() !=
+            numbers::invalid_subdomain_id &&
+          this->level_subdomain_id() == numbers::artificial_subdomain_id);
 #endif
 }
 
@@ -3853,6 +3861,17 @@ CellAccessor<dim, spacedim>::subdomain_id() const
 
 
 
+template <int dim, int spacedim>
+inline types::subdomain_id
+CellAccessor<dim, spacedim>::level_subdomain_id() const
+{
+  Assert(this->used(), TriaAccessorExceptions::ExcCellNotUsed());
+  return this->tria->levels[this->present_level]
+    ->level_subdomain_ids[this->present_index];
+}
+
+
+
 template <int dim, int spacedim>
 inline unsigned int
 CellAccessor<dim, spacedim>::neighbor_face_no(const unsigned int neighbor) const
index b9d192e15431dc8b0b3289530aaeb0bef07a040e..2b0b436f6f606746cdbd52132f69ebcaf79b181e 100644 (file)
@@ -2227,17 +2227,6 @@ CellAccessor<dim, spacedim>::set_subdomain_id(
 
 
 
-template <int dim, int spacedim>
-types::subdomain_id
-CellAccessor<dim, spacedim>::level_subdomain_id() const
-{
-  Assert(this->used(), TriaAccessorExceptions::ExcCellNotUsed());
-  return this->tria->levels[this->present_level]
-    ->level_subdomain_ids[this->present_index];
-}
-
-
-
 template <int dim, int spacedim>
 void
 CellAccessor<dim, spacedim>::set_level_subdomain_id(
diff --git a/tests/multigrid/is_locally_owned_on_level.cc b/tests/multigrid/is_locally_owned_on_level.cc
new file mode 100644 (file)
index 0000000..c1af646
--- /dev/null
@@ -0,0 +1,74 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2022 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+
+// check CellAccessor::is_locally_owned_on_level and is_ghost_on_level for a
+// serial triangulation (where trivially all cells are locally owned and none
+// is ghosted/artificial), going through a code path that is not tested
+// otherwise
+
+#include <deal.II/base/utilities.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_out.h>
+#include <deal.II/grid/grid_tools.h>
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/tria_accessor.h>
+
+#include "../tests.h"
+
+
+
+template <int dim>
+void
+test()
+{
+  Triangulation<dim> tr;
+
+  GridGenerator::subdivided_hyper_cube(tr, 2);
+
+  for (const auto &cell : tr.cell_iterators())
+    {
+      if (cell->is_locally_owned_on_level())
+        {
+          deallog << cell << ": locally owned" << std::endl;
+          Assert(!cell->is_ghost_on_level() && !cell->is_artificial_on_level(),
+                 ExcInternalError());
+        }
+      if (cell->is_ghost_on_level())
+        {
+          deallog << cell << ": ghost" << std::endl;
+          Assert(!cell->is_locally_owned_on_level() &&
+                   !cell->is_artificial_on_level(),
+                 ExcInternalError());
+        }
+      if (cell->is_artificial_on_level())
+        {
+          deallog << cell << ": artificial" << std::endl;
+          Assert(!cell->is_locally_owned_on_level() &&
+                   !cell->is_ghost_on_level(),
+                 ExcInternalError());
+        }
+    }
+}
+
+
+int
+main()
+{
+  initlog();
+  test<2>();
+}
diff --git a/tests/multigrid/is_locally_owned_on_level.output b/tests/multigrid/is_locally_owned_on_level.output
new file mode 100644 (file)
index 0000000..9b6e4d4
--- /dev/null
@@ -0,0 +1,5 @@
+
+DEAL::0.0: locally owned
+DEAL::0.1: locally owned
+DEAL::0.2: locally owned
+DEAL::0.3: locally owned

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.