]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add some more tests for iterator filters
authorJean-Paul Pelteret <jppelteret@gmail.com>
Wed, 5 Jan 2022 18:15:21 +0000 (19:15 +0100)
committerJean-Paul Pelteret <jppelteret@gmail.com>
Thu, 6 Jan 2022 08:13:13 +0000 (09:13 +0100)
tests/grid/filtered_iterator_07.cc [new file with mode: 0644]
tests/grid/filtered_iterator_07.output [new file with mode: 0644]
tests/grid/filtered_iterator_07_operator.cc [new file with mode: 0644]
tests/grid/filtered_iterator_07_operator.output [new file with mode: 0644]

diff --git a/tests/grid/filtered_iterator_07.cc b/tests/grid/filtered_iterator_07.cc
new file mode 100644 (file)
index 0000000..e4d300d
--- /dev/null
@@ -0,0 +1,129 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2021 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 that we can create iterators filtered using ActiveFEIndexEqualTo.
+
+
+#include <deal.II/dofs/dof_handler.h>
+
+#include <deal.II/fe/fe_q.h>
+
+#include <deal.II/grid/filtered_iterator.h>
+#include <deal.II/grid/grid_generator.h>
+
+#include <deal.II/hp/fe_collection.h>
+
+#include "../tests.h"
+
+
+template <int dim, int spacedim>
+void
+test()
+{
+  Triangulation<dim, spacedim> tria;
+  GridGenerator::hyper_cube(tria, 0, 1, true);
+  tria.refine_global(1);
+
+  const FE_Q<dim, spacedim> fe(1);
+
+  DoFHandler<dim, spacedim> dof_handler(tria);
+  dof_handler.distribute_dofs(fe);
+
+  int n_cells_visited = 0;
+  for (const auto &cell :
+       filter_iterators(dof_handler.active_cell_iterators(),
+                        IteratorFilters::ActiveFEIndexEqualTo(0)))
+    {
+      ++n_cells_visited;
+      (void)cell;
+    }
+
+  Assert(n_cells_visited == Utilities::fixed_power<dim>(2),
+         ExcMessage("Wrong number of cells visited."));
+
+  n_cells_visited = 0;
+  for (const auto &cell :
+       filter_iterators(dof_handler.active_cell_iterators(),
+                        IteratorFilters::ActiveFEIndexEqualTo(1)))
+    {
+      ++n_cells_visited;
+      (void)cell;
+    }
+
+  Assert(n_cells_visited == 0, ExcMessage("Wrong number of cells visited."));
+}
+
+
+template <int dim, int spacedim>
+void
+test_hp()
+{
+  Triangulation<dim, spacedim> tria;
+  GridGenerator::hyper_cube(tria, 0, 1, true);
+  tria.refine_global(1);
+
+  const hp::FECollection<dim, spacedim> fe_collection{FE_Q<dim, spacedim>(1),
+                                                      FE_Q<dim, spacedim>(1)};
+
+  DoFHandler<dim, spacedim> dof_handler(tria);
+
+  for (auto &cell : dof_handler.active_cell_iterators())
+    if (cell == dof_handler.begin_active())
+      cell->set_active_fe_index(1);
+    else
+      cell->set_active_fe_index(0);
+
+  dof_handler.distribute_dofs(fe_collection);
+
+  int n_cells_visited = 0;
+  for (const auto &cell :
+       filter_iterators(dof_handler.active_cell_iterators(),
+                        IteratorFilters::ActiveFEIndexEqualTo(0)))
+    {
+      ++n_cells_visited;
+      (void)cell;
+    }
+
+  Assert(n_cells_visited == Utilities::fixed_power<dim>(2) - 1,
+         ExcMessage("Wrong number of cells visited."));
+
+  n_cells_visited = 0;
+  for (const auto &cell :
+       filter_iterators(dof_handler.active_cell_iterators(),
+                        IteratorFilters::ActiveFEIndexEqualTo(1)))
+    {
+      ++n_cells_visited;
+      (void)cell;
+    }
+
+  Assert(n_cells_visited == 1, ExcMessage("Wrong number of cells visited."));
+}
+
+
+int
+main()
+{
+  initlog();
+
+  test<2, 2>();
+  test<2, 3>();
+  test<3, 3>();
+
+  test_hp<2, 2>();
+  test_hp<2, 3>();
+  test_hp<3, 3>();
+
+  deallog << "OK" << std::endl;
+}
diff --git a/tests/grid/filtered_iterator_07.output b/tests/grid/filtered_iterator_07.output
new file mode 100644 (file)
index 0000000..0fd8fc1
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL::OK
diff --git a/tests/grid/filtered_iterator_07_operator.cc b/tests/grid/filtered_iterator_07_operator.cc
new file mode 100644 (file)
index 0000000..6af9dc7
--- /dev/null
@@ -0,0 +1,128 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2021 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 that we can create iterators filtered using ActiveFEIndexEqualTo.
+//
+// Compared to filtered_iterator_07, this test checks the availability
+// of operator| to create the filter.
+
+
+#include <deal.II/dofs/dof_handler.h>
+
+#include <deal.II/fe/fe_q.h>
+
+#include <deal.II/grid/filtered_iterator.h>
+#include <deal.II/grid/grid_generator.h>
+
+#include <deal.II/hp/fe_collection.h>
+
+#include "../tests.h"
+
+
+template <int dim, int spacedim>
+void
+test()
+{
+  Triangulation<dim, spacedim> tria;
+  GridGenerator::hyper_cube(tria, 0, 1, true);
+  tria.refine_global(1);
+
+  const FE_Q<dim, spacedim> fe(1);
+
+  DoFHandler<dim, spacedim> dof_handler(tria);
+  dof_handler.distribute_dofs(fe);
+
+  int n_cells_visited = 0;
+  for (const auto &cell : dof_handler.active_cell_iterators() |
+                            IteratorFilters::ActiveFEIndexEqualTo(0))
+    {
+      ++n_cells_visited;
+      (void)cell;
+    }
+
+  Assert(n_cells_visited == Utilities::fixed_power<dim>(2),
+         ExcMessage("Wrong number of cells visited."));
+
+  n_cells_visited = 0;
+  for (const auto &cell : dof_handler.active_cell_iterators() |
+                            IteratorFilters::ActiveFEIndexEqualTo(1))
+    {
+      ++n_cells_visited;
+      (void)cell;
+    }
+
+  Assert(n_cells_visited == 0, ExcMessage("Wrong number of cells visited."));
+}
+
+
+template <int dim, int spacedim>
+void
+test_hp()
+{
+  Triangulation<dim, spacedim> tria;
+  GridGenerator::hyper_cube(tria, 0, 1, true);
+  tria.refine_global(1);
+
+  const hp::FECollection<dim, spacedim> fe_collection{FE_Q<dim, spacedim>(1),
+                                                      FE_Q<dim, spacedim>(1)};
+
+  DoFHandler<dim, spacedim> dof_handler(tria);
+
+  for (auto &cell : dof_handler.active_cell_iterators())
+    if (cell == dof_handler.begin_active())
+      cell->set_active_fe_index(1);
+    else
+      cell->set_active_fe_index(0);
+
+  dof_handler.distribute_dofs(fe_collection);
+
+  int n_cells_visited = 0;
+  for (const auto &cell : dof_handler.active_cell_iterators() |
+                            IteratorFilters::ActiveFEIndexEqualTo(0))
+    {
+      ++n_cells_visited;
+      (void)cell;
+    }
+
+  Assert(n_cells_visited == Utilities::fixed_power<dim>(2) - 1,
+         ExcMessage("Wrong number of cells visited."));
+
+  n_cells_visited = 0;
+  for (const auto &cell : dof_handler.active_cell_iterators() |
+                            IteratorFilters::ActiveFEIndexEqualTo(1))
+    {
+      ++n_cells_visited;
+      (void)cell;
+    }
+
+  Assert(n_cells_visited == 1, ExcMessage("Wrong number of cells visited."));
+}
+
+
+int
+main()
+{
+  initlog();
+
+  test<2, 2>();
+  test<2, 3>();
+  test<3, 3>();
+
+  test_hp<2, 2>();
+  test_hp<2, 3>();
+  test_hp<3, 3>();
+
+  deallog << "OK" << std::endl;
+}
diff --git a/tests/grid/filtered_iterator_07_operator.output b/tests/grid/filtered_iterator_07_operator.output
new file mode 100644 (file)
index 0000000..0fd8fc1
--- /dev/null
@@ -0,0 +1,2 @@
+
+DEAL::OK

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.