Skip to content

Commit

Permalink
Merge pull request #13170 from jppelteret/filtered_iterator_test_01
Browse files Browse the repository at this point in the history
Add some more tests for iterator filters
  • Loading branch information
bangerth committed Jan 6, 2022
2 parents bcfa3f6 + 72203a5 commit d6c45b5
Show file tree
Hide file tree
Showing 8 changed files with 389 additions and 14 deletions.
6 changes: 6 additions & 0 deletions doc/news/changes/minor/2022010531SimonStickoJean-PaulPelteret
Original file line number Diff line number Diff line change
@@ -0,0 +1,6 @@
Fixed: The FEInterfaceValues class is now made compatible with filtered
iterators. Specifically, it is now possible to call FEInterfaceValues::reinit()
with the active `cell` being a filtered iterator and the type returned by
`cell->neighbor()` being some other non-filtered iterator type.
<br>
(Simon Sticko, Jean-Paul Pelteret, 2022/01/05)
28 changes: 14 additions & 14 deletions include/deal.II/fe/fe_interface_values.h
Original file line number Diff line number Diff line change
Expand Up @@ -1405,14 +1405,14 @@ class FEInterfaceValues
* @param[in] sub_face_no_neighbor Like `sub_face_no`, just for the
* neighboring cell.
*/
template <class CellIteratorType>
template <class CellIteratorType, class CellNeighborIteratorType>
void
reinit(const CellIteratorType & cell,
const unsigned int face_no,
const unsigned int sub_face_no,
const typename identity<CellIteratorType>::type &cell_neighbor,
const unsigned int face_no_neighbor,
const unsigned int sub_face_no_neighbor);
reinit(const CellIteratorType & cell,
const unsigned int face_no,
const unsigned int sub_face_no,
const CellNeighborIteratorType &cell_neighbor,
const unsigned int face_no_neighbor,
const unsigned int sub_face_no_neighbor);

/**
* Re-initialize this object to be used on an interface given by a single face
Expand Down Expand Up @@ -2028,15 +2028,15 @@ FEInterfaceValues<dim, spacedim>::FEInterfaceValues(


template <int dim, int spacedim>
template <class CellIteratorType>
template <class CellIteratorType, class CellNeighborIteratorType>
void
FEInterfaceValues<dim, spacedim>::reinit(
const CellIteratorType & cell,
const unsigned int face_no,
const unsigned int sub_face_no,
const typename identity<CellIteratorType>::type &cell_neighbor,
const unsigned int face_no_neighbor,
const unsigned int sub_face_no_neighbor)
const CellIteratorType & cell,
const unsigned int face_no,
const unsigned int sub_face_no,
const CellNeighborIteratorType &cell_neighbor,
const unsigned int face_no_neighbor,
const unsigned int sub_face_no_neighbor)
{
if (sub_face_no == numbers::invalid_unsigned_int)
{
Expand Down
129 changes: 129 additions & 0 deletions tests/grid/filtered_iterator_07.cc
Original file line number Diff line number Diff line change
@@ -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;
}
2 changes: 2 additions & 0 deletions tests/grid/filtered_iterator_07.output
Original file line number Diff line number Diff line change
@@ -0,0 +1,2 @@

DEAL::OK
128 changes: 128 additions & 0 deletions tests/grid/filtered_iterator_07_operator.cc
Original file line number Diff line number Diff line change
@@ -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;
}
2 changes: 2 additions & 0 deletions tests/grid/filtered_iterator_07_operator.output
Original file line number Diff line number Diff line change
@@ -0,0 +1,2 @@

DEAL::OK

0 comments on commit d6c45b5

Please sign in to comment.