Skip to content

Commit

Permalink
tests for mg set_cell_selection
Browse files Browse the repository at this point in the history
- implement DataOut::get_cell_selection
- add tests to check using set_cell_selection with non-active tests

Basis for of #9100
  • Loading branch information
tjhei committed Jan 31, 2020
1 parent 7762a06 commit 82fb578
Show file tree
Hide file tree
Showing 7 changed files with 843 additions and 0 deletions.
22 changes: 22 additions & 0 deletions include/deal.II/numerics/data_out.h
Expand Up @@ -164,6 +164,21 @@ class DataOut : public DataOut_DoFData<DoFHandlerType,
DoFHandlerType::dimension,
DoFHandlerType::space_dimension>::cell_iterator;

/**
* The type of the function object returning the first cell as used in
* set_cell_selection().
*/
using FirstCellFunctionType =
typename std::function<cell_iterator(const Triangulation<dim, spacedim> &)>;

/**
* The type of the function object returning the next cell as used in
* set_cell_selection().
*/
using NextCellFunctionType =
typename std::function<cell_iterator(const Triangulation<dim, spacedim> &,
const cell_iterator &)>;

/**
* Enumeration describing the part of the domain in which cells
* should be written with curved boundaries. In reality, no file
Expand Down Expand Up @@ -400,6 +415,13 @@ class DataOut : public DataOut_DoFData<DoFHandlerType,
void
set_cell_selection(const FilteredIterator<cell_iterator> &filtered_iterator);

/**
* Return the two function objects that are in use for determining the first
* and the next cell as set by set_cell_selection().
*/
const std::pair<FirstCellFunctionType, NextCellFunctionType>
get_cell_selection() const;

/**
* Return the first cell which we want output for. The default
* implementation returns the first active cell, but you might want to
Expand Down
10 changes: 10 additions & 0 deletions source/numerics/data_out.cc
Expand Up @@ -1107,6 +1107,16 @@ DataOut<dim, DoFHandlerType>::set_cell_selection(



template <int dim, typename DoFHandlerType>
const std::pair<typename DataOut<dim, DoFHandlerType>::FirstCellFunctionType,
typename DataOut<dim, DoFHandlerType>::NextCellFunctionType>
DataOut<dim, DoFHandlerType>::get_cell_selection() const
{
return std::make_pair(first_cell_function, next_cell_function);
}



template <int dim, typename DoFHandlerType>
typename DataOut<dim, DoFHandlerType>::cell_iterator
DataOut<dim, DoFHandlerType>::first_cell()
Expand Down
99 changes: 99 additions & 0 deletions tests/multigrid/mg_data_out_01.cc
@@ -0,0 +1,99 @@
// ---------------------------------------------------------------------
//
// Copyright (C) 2006 - 2016 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 at
// the top level of the deal.II distribution.
//
// ---------------------------------------------------------------------

// Test sequential DataOut::set_cell_selection for multilevel cells

#include <deal.II/dofs/dof_handler.h>

#include <deal.II/grid/grid_generator.h>
#include <deal.II/grid/tria.h>
#include <deal.II/grid/tria_accessor.h>
#include <deal.II/grid/tria_iterator.h>

#include <deal.II/numerics/data_out.h>

#include "../tests.h"

using namespace dealii;

template <int dim>
void
print(DataOut<dim> &data_out, const Triangulation<dim> &tria)
{
auto &p = data_out.get_cell_selection();

for (auto cell = p.first(tria); cell.state() == IteratorState::valid;
cell = p.second(tria, cell))
{
deallog << " cell lvl=" << cell->level() << " active=" << cell->active()
<< " id=" << cell->id().to_string() << std::endl;
}
}

template <int dim>
void
do_test()
{
Triangulation<dim> triangulation(
Triangulation<dim>::limit_level_difference_at_vertices);
GridGenerator::hyper_cube(triangulation);
triangulation.refine_global(1);
triangulation.begin_active()->set_refine_flag();
triangulation.execute_coarsening_and_refinement();

deallog << "dim= " << dim << " cells=" << triangulation.n_cells()
<< std::endl;

DataOut<dim> data_out;
data_out.attach_triangulation(triangulation);
deallog << "* default:" << std::endl;
print(data_out, triangulation);

deallog << "* all cells:" << std::endl;
data_out.set_cell_selection(
[](const typename Triangulation<dim>::cell_iterator &cell) {
return true;
});
print(data_out, triangulation);

deallog << "* all cells with level <=1:" << std::endl;
data_out.set_cell_selection(
[](const typename Triangulation<dim>::cell_iterator &cell) {
return (cell->level() <= 1);
});
print(data_out, triangulation);

for (unsigned int level = 0; level < triangulation.n_levels(); ++level)
{
deallog << "* LevelEqualTo " << level << std::endl;
DataOut<dim> data_out;
data_out.attach_triangulation(triangulation);

data_out.set_cell_selection(IteratorFilters::LevelEqualTo(level));
print(data_out, triangulation);
}
}


int
main(int argc, char **argv)
{
Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1);
mpi_initlog();

do_test<2>();
do_test<3>();
return 0;
}
103 changes: 103 additions & 0 deletions tests/multigrid/mg_data_out_01.output
@@ -0,0 +1,103 @@

DEAL::dim= 2 cells=9
DEAL::* default:
DEAL:: cell lvl=1 active=1 id=0_1:1
DEAL:: cell lvl=1 active=1 id=0_1:2
DEAL:: cell lvl=1 active=1 id=0_1:3
DEAL:: cell lvl=2 active=1 id=0_2:00
DEAL:: cell lvl=2 active=1 id=0_2:01
DEAL:: cell lvl=2 active=1 id=0_2:02
DEAL:: cell lvl=2 active=1 id=0_2:03
DEAL::* all cells:
DEAL:: cell lvl=0 active=0 id=0_0:
DEAL:: cell lvl=1 active=0 id=0_1:0
DEAL:: cell lvl=1 active=1 id=0_1:1
DEAL:: cell lvl=1 active=1 id=0_1:2
DEAL:: cell lvl=1 active=1 id=0_1:3
DEAL:: cell lvl=2 active=1 id=0_2:00
DEAL:: cell lvl=2 active=1 id=0_2:01
DEAL:: cell lvl=2 active=1 id=0_2:02
DEAL:: cell lvl=2 active=1 id=0_2:03
DEAL::* all cells with level <=1:
DEAL:: cell lvl=0 active=0 id=0_0:
DEAL:: cell lvl=1 active=0 id=0_1:0
DEAL:: cell lvl=1 active=1 id=0_1:1
DEAL:: cell lvl=1 active=1 id=0_1:2
DEAL:: cell lvl=1 active=1 id=0_1:3
DEAL::* LevelEqualTo 0
DEAL:: cell lvl=0 active=0 id=0_0:
DEAL::* LevelEqualTo 1
DEAL:: cell lvl=1 active=0 id=0_1:0
DEAL:: cell lvl=1 active=1 id=0_1:1
DEAL:: cell lvl=1 active=1 id=0_1:2
DEAL:: cell lvl=1 active=1 id=0_1:3
DEAL::* LevelEqualTo 2
DEAL:: cell lvl=2 active=1 id=0_2:00
DEAL:: cell lvl=2 active=1 id=0_2:01
DEAL:: cell lvl=2 active=1 id=0_2:02
DEAL:: cell lvl=2 active=1 id=0_2:03
DEAL::dim= 3 cells=17
DEAL::* default:
DEAL:: cell lvl=1 active=1 id=0_1:1
DEAL:: cell lvl=1 active=1 id=0_1:2
DEAL:: cell lvl=1 active=1 id=0_1:3
DEAL:: cell lvl=1 active=1 id=0_1:4
DEAL:: cell lvl=1 active=1 id=0_1:5
DEAL:: cell lvl=1 active=1 id=0_1:6
DEAL:: cell lvl=1 active=1 id=0_1:7
DEAL:: cell lvl=2 active=1 id=0_2:00
DEAL:: cell lvl=2 active=1 id=0_2:01
DEAL:: cell lvl=2 active=1 id=0_2:02
DEAL:: cell lvl=2 active=1 id=0_2:03
DEAL:: cell lvl=2 active=1 id=0_2:04
DEAL:: cell lvl=2 active=1 id=0_2:05
DEAL:: cell lvl=2 active=1 id=0_2:06
DEAL:: cell lvl=2 active=1 id=0_2:07
DEAL::* all cells:
DEAL:: cell lvl=0 active=0 id=0_0:
DEAL:: cell lvl=1 active=0 id=0_1:0
DEAL:: cell lvl=1 active=1 id=0_1:1
DEAL:: cell lvl=1 active=1 id=0_1:2
DEAL:: cell lvl=1 active=1 id=0_1:3
DEAL:: cell lvl=1 active=1 id=0_1:4
DEAL:: cell lvl=1 active=1 id=0_1:5
DEAL:: cell lvl=1 active=1 id=0_1:6
DEAL:: cell lvl=1 active=1 id=0_1:7
DEAL:: cell lvl=2 active=1 id=0_2:00
DEAL:: cell lvl=2 active=1 id=0_2:01
DEAL:: cell lvl=2 active=1 id=0_2:02
DEAL:: cell lvl=2 active=1 id=0_2:03
DEAL:: cell lvl=2 active=1 id=0_2:04
DEAL:: cell lvl=2 active=1 id=0_2:05
DEAL:: cell lvl=2 active=1 id=0_2:06
DEAL:: cell lvl=2 active=1 id=0_2:07
DEAL::* all cells with level <=1:
DEAL:: cell lvl=0 active=0 id=0_0:
DEAL:: cell lvl=1 active=0 id=0_1:0
DEAL:: cell lvl=1 active=1 id=0_1:1
DEAL:: cell lvl=1 active=1 id=0_1:2
DEAL:: cell lvl=1 active=1 id=0_1:3
DEAL:: cell lvl=1 active=1 id=0_1:4
DEAL:: cell lvl=1 active=1 id=0_1:5
DEAL:: cell lvl=1 active=1 id=0_1:6
DEAL:: cell lvl=1 active=1 id=0_1:7
DEAL::* LevelEqualTo 0
DEAL:: cell lvl=0 active=0 id=0_0:
DEAL::* LevelEqualTo 1
DEAL:: cell lvl=1 active=0 id=0_1:0
DEAL:: cell lvl=1 active=1 id=0_1:1
DEAL:: cell lvl=1 active=1 id=0_1:2
DEAL:: cell lvl=1 active=1 id=0_1:3
DEAL:: cell lvl=1 active=1 id=0_1:4
DEAL:: cell lvl=1 active=1 id=0_1:5
DEAL:: cell lvl=1 active=1 id=0_1:6
DEAL:: cell lvl=1 active=1 id=0_1:7
DEAL::* LevelEqualTo 2
DEAL:: cell lvl=2 active=1 id=0_2:00
DEAL:: cell lvl=2 active=1 id=0_2:01
DEAL:: cell lvl=2 active=1 id=0_2:02
DEAL:: cell lvl=2 active=1 id=0_2:03
DEAL:: cell lvl=2 active=1 id=0_2:04
DEAL:: cell lvl=2 active=1 id=0_2:05
DEAL:: cell lvl=2 active=1 id=0_2:06
DEAL:: cell lvl=2 active=1 id=0_2:07
109 changes: 109 additions & 0 deletions tests/multigrid/mg_data_out_02.cc
@@ -0,0 +1,109 @@
// ---------------------------------------------------------------------
//
// Copyright (C) 2006 - 2016 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 at
// the top level of the deal.II distribution.
//
// ---------------------------------------------------------------------

// Test parallel DataOut::set_cell_selection for multilevel cells

#include <deal.II/distributed/tria.h>

#include <deal.II/dofs/dof_handler.h>

#include <deal.II/grid/grid_generator.h>
#include <deal.II/grid/tria_accessor.h>
#include <deal.II/grid/tria_iterator.h>

#include <deal.II/numerics/data_out.h>

#include "../tests.h"

using namespace dealii;

template <int dim>
void
print(DataOut<dim> &data_out, const Triangulation<dim> &tria)
{
auto &p = data_out.get_cell_selection();

for (auto cell = p.first(tria); cell.state() == IteratorState::valid;
cell = p.second(tria, cell))
{
deallog << " cell lvl=" << cell->level() << " active=" << cell->active()
<< " owner=" << static_cast<int>(cell->level_subdomain_id())
<< " id=" << cell->id().to_string() << std::endl;
}
}

template <int dim>
void
do_test()
{
parallel::distributed::Triangulation<dim> triangulation(
MPI_COMM_WORLD,
dealii::Triangulation<dim>::none,
parallel::distributed::Triangulation<dim>::construct_multigrid_hierarchy);

GridGenerator::hyper_cube(triangulation);
triangulation.refine_global(1);
triangulation.begin_active()->set_refine_flag();
triangulation.execute_coarsening_and_refinement();

deallog << "dim= " << dim << " cells=" << triangulation.n_cells()
<< std::endl;

DataOut<dim> data_out;
data_out.attach_triangulation(triangulation);

deallog << "* default:" << std::endl;
print(data_out, triangulation);

deallog << "* all cells:" << std::endl;
data_out.set_cell_selection(
[](const typename Triangulation<dim>::cell_iterator &cell) {
return true;
});
print(data_out, triangulation);

deallog << "* LocallyOwnedLevelCell:" << std::endl;
data_out.set_cell_selection(IteratorFilters::LocallyOwnedLevelCell());
print(data_out, triangulation);

for (unsigned int level = 0; level < triangulation.n_global_levels(); ++level)
{
DataOut<dim> data_out;
data_out.attach_triangulation(triangulation);

deallog << "* LevelEqualTo " << level << std::endl;
data_out.set_cell_selection(IteratorFilters::LevelEqualTo(level));
print(data_out, triangulation);

deallog << "* owned on level " << level << std::endl;
data_out.set_cell_selection(
[level](const typename Triangulation<dim>::cell_iterator &cell) {
return (cell->level() == level && cell->is_locally_owned_on_level());
});
print(data_out, triangulation);
}
}


int
main(int argc, char **argv)
{
Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1);
MPILogInitAll log;

do_test<2>();
do_test<3>();
return 0;
}

0 comments on commit 82fb578

Please sign in to comment.