]> https://gitweb.dealii.org/ - dealii.git/commitdiff
add a test for multilevel access of FEInterfaceValues
authorTimo Heister <timo.heister@gmail.com>
Wed, 10 Jun 2020 14:29:33 +0000 (10:29 -0400)
committerTimo Heister <timo.heister@gmail.com>
Wed, 10 Jun 2020 14:50:08 +0000 (10:50 -0400)
part of #8884

adds a test for change in #9978

tests/feinterface/multigrid_01.cc [new file with mode: 0644]
tests/feinterface/multigrid_01.output [new file with mode: 0644]

diff --git a/tests/feinterface/multigrid_01.cc b/tests/feinterface/multigrid_01.cc
new file mode 100644 (file)
index 0000000..24702a7
--- /dev/null
@@ -0,0 +1,219 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2018 - 2020 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 basic properties of FEInterfaceValues for multilevel
+
+#include <deal.II/base/quadrature_lib.h>
+
+#include <deal.II/fe/fe_dgq.h>
+#include <deal.II/fe/fe_interface_values.h>
+#include <deal.II/fe/mapping_q.h>
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_refinement.h>
+#include <deal.II/grid/tria.h>
+
+#include <fstream>
+#include <iostream>
+
+#include "../tests.h"
+
+
+template <int dim>
+void
+inspect_fiv(FEInterfaceValues<dim> &fiv)
+{
+  deallog << "at_boundary(): " << fiv.at_boundary() << "\n"
+          << "n_current_interface_dofs(): " << fiv.n_current_interface_dofs()
+          << "\n";
+
+  std::vector<types::global_dof_index> indices =
+    fiv.get_interface_dof_indices();
+  Assert(indices.size() == fiv.n_current_interface_dofs(), ExcInternalError());
+
+  deallog << "interface_dof_indices: ";
+  for (auto i : indices)
+    deallog << i << " ";
+  deallog << "\n";
+
+
+  unsigned int idx = 0;
+  for (auto v : indices)
+    {
+      deallog << "  index " << idx << " global_dof_index:" << v << ":\n";
+
+      const auto pair = fiv.interface_dof_to_dof_indices(idx);
+      deallog << "    dof indices: " << pair[0] << " | " << pair[1] << "\n";
+
+      ++idx;
+    }
+
+  deallog << "\n";
+  deallog << std::endl;
+}
+
+
+template <int dim>
+void
+make_2_cells(Triangulation<dim> &tria);
+
+template <>
+void make_2_cells<2>(Triangulation<2> &tria)
+{
+  const unsigned int        dim         = 2;
+  std::vector<unsigned int> repetitions = {2, 1};
+  Point<dim>                p1;
+  Point<dim>                p2(2.0, 1.0);
+
+  GridGenerator::subdivided_hyper_rectangle(tria, repetitions, p1, p2);
+}
+
+template <>
+void make_2_cells<3>(Triangulation<3> &tria)
+{
+  const unsigned int        dim         = 3;
+  std::vector<unsigned int> repetitions = {2, 1, 1};
+  Point<dim>                p1;
+  Point<dim>                p2(2.0, 1.0, 1.0);
+
+  GridGenerator::subdivided_hyper_rectangle(tria, repetitions, p1, p2);
+}
+
+
+template <int dim>
+void
+test()
+{
+  Triangulation<dim> tria(
+    Triangulation<dim>::limit_level_difference_at_vertices);
+  make_2_cells(tria);
+  tria.refine_global();
+
+
+  DoFHandler<dim> dofh(tria);
+  FE_DGQ<dim>     fe(1);
+  dofh.distribute_dofs(fe);
+  dofh.distribute_mg_dofs();
+
+  deallog << "n_dofs = " << dofh.n_dofs() << " n_dofs(0)= " << dofh.n_dofs(0)
+          << " n_dofs(1)= " << dofh.n_dofs(1) << std::endl;
+
+  MappingQ<dim> mapping(1);
+  UpdateFlags   update_flags = update_values | update_gradients |
+                             update_quadrature_points | update_JxW_values;
+
+  FEInterfaceValues<dim> fiv(mapping,
+                             fe,
+                             QGauss<dim - 1>(fe.degree + 1),
+                             update_flags);
+
+
+  auto cell = dofh.begin_mg(0);
+
+  deallog << "** interface between cell 0 and 1 on level 0: **\n";
+
+  for (const unsigned int f : GeometryInfo<dim>::face_indices())
+    if (!cell->at_boundary(f))
+      {
+        fiv.reinit(cell,
+                   f,
+                   numbers::invalid_unsigned_int,
+                   cell->neighbor(f),
+                   cell->neighbor_of_neighbor(f),
+                   numbers::invalid_unsigned_int);
+
+        Assert(fiv.get_fe_face_values(0).get_cell() == cell,
+               ExcInternalError());
+        Assert(fiv.get_fe_face_values(1).get_cell() == cell->neighbor(f),
+               ExcInternalError());
+        Assert(fiv.n_current_interface_dofs() == 2 * fe.n_dofs_per_cell(),
+               ExcInternalError());
+        Assert(!fiv.at_boundary(), ExcInternalError());
+
+        auto mycell = cell;
+        for (unsigned int c = 0; c < 2; ++c)
+          {
+            std::vector<types::global_dof_index> indices(fe.n_dofs_per_cell());
+            mycell->get_mg_dof_indices(indices);
+            deallog << "cell " << c << ": ";
+            for (auto i : indices)
+              deallog << i << " ";
+            deallog << "\n";
+            ++mycell;
+          }
+
+        inspect_fiv(fiv);
+      }
+
+
+  deallog << "** interface on level 1: **\n";
+  cell = dofh.begin_mg(1);
+
+  for (const unsigned int f : GeometryInfo<dim>::face_indices())
+    if (!cell->at_boundary(f))
+      {
+        fiv.reinit(cell,
+                   f,
+                   numbers::invalid_unsigned_int,
+                   cell->neighbor(f),
+                   cell->neighbor_of_neighbor(f),
+                   numbers::invalid_unsigned_int);
+
+        Assert(fiv.get_fe_face_values(0).get_cell() == cell,
+               ExcInternalError());
+        Assert(fiv.get_fe_face_values(1).get_cell() == cell->neighbor(f),
+               ExcInternalError());
+        Assert(fiv.n_current_interface_dofs() == 2 * fe.n_dofs_per_cell(),
+               ExcInternalError());
+        Assert(!fiv.at_boundary(), ExcInternalError());
+
+        auto mycell = cell;
+        for (unsigned int c = 0; c < 2; ++c)
+          {
+            std::vector<types::global_dof_index> indices(fe.n_dofs_per_cell());
+            mycell->get_mg_dof_indices(indices);
+            deallog << "cell " << c << ": ";
+            for (auto i : indices)
+              deallog << i << " ";
+            deallog << "\n";
+            mycell = cell->neighbor(f);
+          }
+
+        inspect_fiv(fiv);
+        break; // only look at the first internal face
+      }
+
+  deallog << "** boundary interface on cell 1 **\n";
+  cell = dofh.begin_mg(0);
+  {
+    ++cell;
+    fiv.reinit(cell, 1);
+    Assert(fiv.get_fe_face_values(0).get_cell() == cell, ExcInternalError());
+    Assert(fiv.n_current_interface_dofs() == fe.n_dofs_per_cell(),
+           ExcInternalError());
+    Assert(fiv.at_boundary(), ExcInternalError());
+    inspect_fiv(fiv);
+  }
+}
+
+
+
+int
+main()
+{
+  initlog();
+  test<2>();
+}
diff --git a/tests/feinterface/multigrid_01.output b/tests/feinterface/multigrid_01.output
new file mode 100644 (file)
index 0000000..96fdf64
--- /dev/null
@@ -0,0 +1,64 @@
+
+DEAL::n_dofs = 32 n_dofs(0)= 8 n_dofs(1)= 32
+DEAL::** interface between cell 0 and 1 on level 0: **
+cell 0: 0 1 2 3 
+cell 1: 4 5 6 7 
+at_boundary(): 0
+n_current_interface_dofs(): 8
+interface_dof_indices: 0 1 2 3 4 5 6 7 
+  index 0 global_dof_index:0:
+    dof indices: 0 | 4294967295
+  index 1 global_dof_index:1:
+    dof indices: 1 | 4294967295
+  index 2 global_dof_index:2:
+    dof indices: 2 | 4294967295
+  index 3 global_dof_index:3:
+    dof indices: 3 | 4294967295
+  index 4 global_dof_index:4:
+    dof indices: 4294967295 | 0
+  index 5 global_dof_index:5:
+    dof indices: 4294967295 | 1
+  index 6 global_dof_index:6:
+    dof indices: 4294967295 | 2
+  index 7 global_dof_index:7:
+    dof indices: 4294967295 | 3
+
+
+DEAL::** interface on level 1: **
+cell 0: 0 1 2 3 
+cell 1: 4 5 6 7 
+at_boundary(): 0
+n_current_interface_dofs(): 8
+interface_dof_indices: 0 1 2 3 4 5 6 7 
+  index 0 global_dof_index:0:
+    dof indices: 0 | 4294967295
+  index 1 global_dof_index:1:
+    dof indices: 1 | 4294967295
+  index 2 global_dof_index:2:
+    dof indices: 2 | 4294967295
+  index 3 global_dof_index:3:
+    dof indices: 3 | 4294967295
+  index 4 global_dof_index:4:
+    dof indices: 4294967295 | 0
+  index 5 global_dof_index:5:
+    dof indices: 4294967295 | 1
+  index 6 global_dof_index:6:
+    dof indices: 4294967295 | 2
+  index 7 global_dof_index:7:
+    dof indices: 4294967295 | 3
+
+
+DEAL::** boundary interface on cell 1 **
+at_boundary(): 1
+n_current_interface_dofs(): 4
+interface_dof_indices: 4 5 6 7 
+  index 0 global_dof_index:4:
+    dof indices: 0 | 4294967295
+  index 1 global_dof_index:5:
+    dof indices: 1 | 4294967295
+  index 2 global_dof_index:6:
+    dof indices: 2 | 4294967295
+  index 3 global_dof_index:7:
+    dof indices: 3 | 4294967295
+
+

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.