From: Marc Fehling Date: Thu, 17 Mar 2022 01:39:44 +0000 (-0600) Subject: cell_weight: added particle weighting test. X-Git-Tag: v9.4.0-rc1~352^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=08bfbc7057d55309c2bebbc44a20205ecafb540c;p=dealii.git cell_weight: added particle weighting test. --- diff --git a/tests/particles/particle_weights_01.cc b/tests/particles/particle_weights_01.cc new file mode 100644 index 0000000000..8d4825d48e --- /dev/null +++ b/tests/particles/particle_weights_01.cc @@ -0,0 +1,128 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2022 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. +// +// --------------------------------------------------------------------- + + + +// verify particle weighting mechanism + + +#include +#include + +#include + +#include + +#include +#include +#include + +#include +#include + +#include "../tests.h" + + +template +void +test() +{ + // initialize grid + parallel::distributed::Triangulation triangulation(MPI_COMM_WORLD); + { + GridGenerator::hyper_cube(triangulation); + triangulation.refine_global(2); + } + + // initialize particle handler + const MappingQ1 mapping; + Particles::ParticleHandler particle_handler(triangulation, mapping); + { + const auto local_bounding_box = + GridTools::compute_mesh_predicate_bounding_box( + triangulation, IteratorFilters::LocallyOwnedCell()); + const auto global_bounding_boxes = + Utilities::MPI::all_gather(triangulation.get_communicator(), + local_bounding_box); + + Particles::Generators::quadrature_points(triangulation, + QMidpoint(), + global_bounding_boxes, + particle_handler); + + Assert(particle_handler.n_locally_owned_particles() == + triangulation.n_locally_owned_active_cells(), + ExcInternalError()); + Assert(particle_handler.n_global_particles() == + triangulation.n_global_active_cells(), + ExcInternalError()); + } + + // register weighting function that returns the number of particles per cell + triangulation.signals.weight.connect( + [&particle_handler]( + const typename parallel::distributed::Triangulation::cell_iterator + &cell, + const typename parallel::distributed::Triangulation::CellStatus + status) -> unsigned int { + unsigned int n_particles_in_cell = 0; + switch (status) + { + case parallel::distributed::Triangulation::CELL_PERSIST: + case parallel::distributed::Triangulation::CELL_REFINE: + n_particles_in_cell = particle_handler.n_particles_in_cell(cell); + break; + + case parallel::distributed::Triangulation::CELL_INVALID: + break; + + case parallel::distributed::Triangulation::CELL_COARSEN: + for (const auto &child : cell->child_iterators()) + n_particles_in_cell += + particle_handler.n_particles_in_cell(child); + break; + + default: + Assert(false, ExcInternalError()); + break; + } + deallog << ' ' << n_particles_in_cell << std::endl; + return n_particles_in_cell; + }); + + // refine one cell, and coarsen one group of siblings + for (const auto &cell : triangulation.active_cell_iterators() | + IteratorFilters::LocallyOwnedCell()) + { + if (cell->id().to_string() == "0_2:00") + cell->set_refine_flag(); + else if (cell->parent()->id().to_string() == "0_1:3") + cell->set_coarsen_flag(); + } + + deallog << "weights before repartitioning:" << std::endl; + triangulation.execute_coarsening_and_refinement(); +} + + +int +main(int argc, char *argv[]) +{ + Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1); + + MPILogInitAll all; + + test<2>(); +} diff --git a/tests/particles/particle_weights_01.with_p4est=true.mpirun=1.output b/tests/particles/particle_weights_01.with_p4est=true.mpirun=1.output new file mode 100644 index 0000000000..2117194459 --- /dev/null +++ b/tests/particles/particle_weights_01.with_p4est=true.mpirun=1.output @@ -0,0 +1,18 @@ + +DEAL:0::weights before repartitioning: +DEAL:0:: 1 +DEAL:0:: 0 +DEAL:0:: 0 +DEAL:0:: 0 +DEAL:0:: 1 +DEAL:0:: 1 +DEAL:0:: 1 +DEAL:0:: 1 +DEAL:0:: 1 +DEAL:0:: 1 +DEAL:0:: 1 +DEAL:0:: 1 +DEAL:0:: 1 +DEAL:0:: 1 +DEAL:0:: 1 +DEAL:0:: 4 diff --git a/tests/particles/particle_weights_01.with_p4est=true.mpirun=4.output b/tests/particles/particle_weights_01.with_p4est=true.mpirun=4.output new file mode 100644 index 0000000000..b9467f0ea4 --- /dev/null +++ b/tests/particles/particle_weights_01.with_p4est=true.mpirun=4.output @@ -0,0 +1,27 @@ + +DEAL:0::weights before repartitioning: +DEAL:0:: 1 +DEAL:0:: 0 +DEAL:0:: 0 +DEAL:0:: 0 +DEAL:0:: 1 +DEAL:0:: 1 +DEAL:0:: 1 + +DEAL:1::weights before repartitioning: +DEAL:1:: 1 +DEAL:1:: 1 +DEAL:1:: 1 +DEAL:1:: 1 + + +DEAL:2::weights before repartitioning: +DEAL:2:: 1 +DEAL:2:: 1 +DEAL:2:: 1 +DEAL:2:: 1 + + +DEAL:3::weights before repartitioning: +DEAL:3:: 4 +