From 127d9a4a4e66ae7dc2e4aba0769d8e7e79a6e59a Mon Sep 17 00:00:00 2001 From: Daniel Arndt Date: Thu, 9 Aug 2018 11:21:02 +0200 Subject: [PATCH] Add test for level difference across periodic faces on serial Triangulations --- tests/grid/periodicity_01.cc | 91 ++++++++++++++++++++++++++++++++ tests/grid/periodicity_01.output | 3 ++ 2 files changed, 94 insertions(+) create mode 100644 tests/grid/periodicity_01.cc create mode 100644 tests/grid/periodicity_01.output diff --git a/tests/grid/periodicity_01.cc b/tests/grid/periodicity_01.cc new file mode 100644 index 0000000000..588e79f458 --- /dev/null +++ b/tests/grid/periodicity_01.cc @@ -0,0 +1,91 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 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. +// +// --------------------------------------------------------------------- + +// Check that there is at most a one level difference between cells +// across periodic boundary faces for serial Triangulations. + +#include + +#include +#include +#include +#include +#include +#include + +#include "../tests.h" + +template +void +test() +{ + double left = -1.; + double right = 1.; + + Triangulation tria; + + GridGenerator::hyper_cube(tria, left, right, true); + + // setup periodic boundary conditions + std::vector< + GridTools::PeriodicFacePair::cell_iterator>> + periodicity_vector; + + // for the hypercube faces + unsigned int id1 = 0; // id of the left boundary + unsigned int id2 = 1; // id of the right boundary boundary + GridTools::collect_periodic_faces( + tria, id1, id2, /*x-direction*/ 0, periodicity_vector); + + tria.add_periodicity(periodicity_vector); + + const unsigned int n_refinements = 5; + const double tolerance = 1.e-3; + + for (unsigned int n = 0; n < n_refinements; ++n) + { + for (const auto &cell : tria.active_cell_iterators()) + if (cell->center()[0] > -tolerance) + cell->set_refine_flag(); + + tria.execute_coarsening_and_refinement(); + } + + // Let's check how many cells have a periodic neighbor. + unsigned int neighbor_count = 0; + + for (const auto &cell : tria.active_cell_iterators()) + for (unsigned int f = 0; f < GeometryInfo<2>::faces_per_cell; ++f) + if (cell->has_periodic_neighbor(f)) + ++neighbor_count; + + deallog << "Found " << neighbor_count << " faces with periodic neighbor" + << std::endl; + + const unsigned int expected_neighbor_count = + Utilities::pow(2, (dim - 1) * (n_refinements - 1)) + + Utilities::pow(2, (dim - 1) * n_refinements); + AssertThrow(neighbor_count == expected_neighbor_count, ExcInternalError()); +} + + +int +main() +{ + initlog(); + test<2>(); + test<3>(); + return 0; +} diff --git a/tests/grid/periodicity_01.output b/tests/grid/periodicity_01.output new file mode 100644 index 0000000000..cf04d7bf87 --- /dev/null +++ b/tests/grid/periodicity_01.output @@ -0,0 +1,3 @@ + +DEAL::Found 48 faces with periodic neighbor +DEAL::Found 1280 faces with periodic neighbor -- 2.39.5