]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fixed mesh_loop in one dimension.
authorLuca Heltai <luca.heltai@sissa.it>
Wed, 8 Jul 2020 14:46:37 +0000 (16:46 +0200)
committerLuca Heltai <luca.heltai@sissa.it>
Wed, 8 Jul 2020 14:47:47 +0000 (16:47 +0200)
doc/news/changes/minor/20200708LucaHeltai [new file with mode: 0644]
include/deal.II/meshworker/mesh_loop.h
tests/meshworker/mesh_loop_08.cc [new file with mode: 0644]
tests/meshworker/mesh_loop_08.output [new file with mode: 0644]

diff --git a/doc/news/changes/minor/20200708LucaHeltai b/doc/news/changes/minor/20200708LucaHeltai
new file mode 100644 (file)
index 0000000..b2e7239
--- /dev/null
@@ -0,0 +1,3 @@
+Fixed: MeshWorker::mesh_loop() did not work on 1d with refined grids. This is now fixed.
+<br>
+(Luca Heltai, 2020/07/08)
index 79445a8ef9eb692892421ced61de0a64c108498e..2c2267de1cc6d3b6dfe643f5a02a0aabd07a3bdc 100644 (file)
@@ -435,6 +435,42 @@ namespace MeshWorker
                                     copy);
                       }
                   }
+                else if (CellIteratorType::AccessorType::dimension == 1 &&
+                         cell->level() > neighbor->level())
+                  {
+                    // In one dimension, there is no other check to do
+                    const unsigned int neighbor_face_no =
+                      periodic_neighbor ?
+                        cell->periodic_neighbor_face_no(face_no) :
+                        cell->neighbor_face_no(face_no);
+                    Assert(periodic_neighbor ||
+                             neighbor->face(neighbor_face_no) ==
+                               cell->face(face_no),
+                           ExcInternalError());
+
+                    face_worker(cell,
+                                face_no,
+                                numbers::invalid_unsigned_int,
+                                neighbor,
+                                neighbor_face_no,
+                                numbers::invalid_unsigned_int,
+                                scratch,
+                                copy);
+
+                    if (flags & assemble_own_interior_faces_both)
+                      {
+                        // If own faces are to be assembled from both sides,
+                        // call the faceworker again with swapped arguments.
+                        face_worker(neighbor,
+                                    neighbor_face_no,
+                                    numbers::invalid_unsigned_int,
+                                    cell,
+                                    face_no,
+                                    numbers::invalid_unsigned_int,
+                                    scratch,
+                                    copy);
+                      }
+                  }
                 else
                   {
                     // If iterator is active and neighbor is refined, skip
diff --git a/tests/meshworker/mesh_loop_08.cc b/tests/meshworker/mesh_loop_08.cc
new file mode 100644 (file)
index 0000000..035ea56
--- /dev/null
@@ -0,0 +1,109 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2007 - 2018 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.
+//
+// ---------------------------------------------------------------------
+
+
+
+// test mesh_loop on 1d meshes with local refinement
+
+#include <deal.II/base/work_stream.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/tria_iterator.h>
+
+#include <deal.II/meshworker/assemble_flags.h>
+#include <deal.II/meshworker/mesh_loop.h>
+
+#include <fstream>
+
+#include "../tests.h"
+
+struct ScratchData
+{};
+struct CopyData
+{};
+
+using namespace MeshWorker;
+
+template <int dim, int spacedim>
+void
+test()
+{
+  Triangulation<dim, spacedim> tria;
+  GridGenerator::hyper_cube(tria);
+  tria.refine_global(1);
+  tria.begin_active()->set_refine_flag();
+  tria.execute_coarsening_and_refinement();
+
+  ScratchData scratch;
+  CopyData    copy;
+
+  auto cell = tria.begin_active();
+  auto endc = tria.end();
+
+  typedef decltype(cell) Iterator;
+
+
+
+  auto cell_worker = [](const Iterator &cell, ScratchData &s, CopyData &c) {
+    deallog << "Cell worker on : " << cell << std::endl;
+  };
+
+  auto boundary_worker =
+    [](const Iterator &cell, const unsigned int &f, ScratchData &, CopyData &) {
+      deallog << "Boundary worker on : " << cell << ", Face : " << f
+              << std::endl;
+    };
+
+  auto face_worker = [](const Iterator &    cell,
+                        const unsigned int &f,
+                        const unsigned int &sf,
+                        const Iterator &    ncell,
+                        const unsigned int &nf,
+                        const unsigned int &nsf,
+                        ScratchData &       s,
+                        CopyData &          c) {
+    deallog << "Face worker on : " << cell << ", Neighbor cell : " << ncell
+            << ", Face : " << f << ", Neighbor Face : " << nf
+            << ", Subface: " << sf << ", Neighbor Subface: " << nsf
+            << std::endl;
+  };
+
+  auto copier = [](const CopyData &) { deallog << "copier" << std::endl; };
+
+  deallog << "ONE DIMENSION" << std::endl << std::endl;
+
+  mesh_loop(cell,
+            endc,
+            cell_worker,
+            copier,
+            scratch,
+            copy,
+            MeshWorker::assemble_own_cells |
+              MeshWorker::assemble_own_interior_faces_once |
+              MeshWorker::assemble_boundary_faces,
+            boundary_worker,
+            face_worker);
+}
+
+
+int
+main()
+{
+  initlog();
+  MultithreadInfo::set_thread_limit(1); // to make output deterministic
+
+  test<1, 1>();
+}
diff --git a/tests/meshworker/mesh_loop_08.output b/tests/meshworker/mesh_loop_08.output
new file mode 100644 (file)
index 0000000..6ef1634
--- /dev/null
@@ -0,0 +1,13 @@
+
+DEAL::ONE DIMENSION
+DEAL::
+DEAL::Cell worker on : 1.1
+DEAL::Boundary worker on : 1.1, Face : 1
+DEAL::copier
+DEAL::Cell worker on : 2.0
+DEAL::Boundary worker on : 2.0, Face : 0
+DEAL::Face worker on : 2.0, Neighbor cell : 2.1, Face : 1, Neighbor Face : 0, Subface: 4294967295, Neighbor Subface: 4294967295
+DEAL::copier
+DEAL::Cell worker on : 2.1
+DEAL::Face worker on : 2.1, Neighbor cell : 1.1, Face : 1, Neighbor Face : 0, Subface: 4294967295, Neighbor Subface: 4294967295
+DEAL::copier

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.