From 3e8d4251eeb3a6d0d6508775757942cae730b549 Mon Sep 17 00:00:00 2001 From: heister Date: Mon, 29 Aug 2011 22:18:26 +0000 Subject: [PATCH] document bug in parallel::distributed::refine_and_coarsen*() with more CPUs than cells git-svn-id: https://svn.dealii.org/trunk@24210 0785d39b-7218-0410-832d-ea1e28bc413d --- .../mpi/refine_and_coarsen_fixed_number_06.cc | 91 +++++++++++++++++++ .../ncpu_1/cmp/generic | 4 + .../ncpu_10/cmp/generic | 4 + .../ncpu_4/cmp/generic | 4 + 4 files changed, 103 insertions(+) create mode 100644 tests/mpi/refine_and_coarsen_fixed_number_06.cc create mode 100644 tests/mpi/refine_and_coarsen_fixed_number_06/ncpu_1/cmp/generic create mode 100644 tests/mpi/refine_and_coarsen_fixed_number_06/ncpu_10/cmp/generic create mode 100644 tests/mpi/refine_and_coarsen_fixed_number_06/ncpu_4/cmp/generic diff --git a/tests/mpi/refine_and_coarsen_fixed_number_06.cc b/tests/mpi/refine_and_coarsen_fixed_number_06.cc new file mode 100644 index 0000000000..0b0326654d --- /dev/null +++ b/tests/mpi/refine_and_coarsen_fixed_number_06.cc @@ -0,0 +1,91 @@ +//--------------------------------------------------------------------------- +// $Id: simple_mpi_01.cc 23327 2011-02-11 03:19:07Z bangerth $ +// Version: $Name$ +// +// Copyright (C) 2009 by the deal.II authors +// +// This file is subject to QPL and may not be distributed +// without copyright and license information. Please refer +// to the file deal.II/doc/license.html for the text and +// further information on this license. +// +//--------------------------------------------------------------------------- + + +// document bug in parallel::distributed::GridRefinement +// ::refine_and_coarsen_fixed_number() and fixed_fraction() with one CPU with +// 0 cells: +//#0 0x00007ffff67297e2 in dealii::(anonymous namespace)::min_element ( +// criteria=...) +// at /w/heister/deal-trunk/deal.II/source/distributed/grid_refinement.cc:57 +//#1 0x00007ffff6727672 in dealii::(anonymous namespace)::compute_global_min_and_max_at_root (criteria=..., mpi_communicator=0x6fc860) +// at /w/heister/deal-trunk/deal.II/source/distributed/grid_refinement.cc:82 +//#2 0x00007ffff672b4ef in dealii::parallel::distributed::GridRefinement::refine_and_coarsen_fixed_number<2, dealii::Vector, 2> + +#include "../tests.h" + +#include +#include +#include +#include +#include +#include +#include +#include +//#include + +template +void test() +{ + unsigned int myid = Utilities::System::get_this_mpi_process (MPI_COMM_WORLD); + unsigned int numprocs = Utilities::System::get_n_mpi_processes (MPI_COMM_WORLD); + + parallel::distributed::Triangulation triangulation(MPI_COMM_WORLD); + GridGenerator::hyper_cube (triangulation); + + Vector estimated_error_per_cell (triangulation.n_active_cells()); + parallel::distributed::GridRefinement:: + refine_and_coarsen_fixed_number (triangulation, + estimated_error_per_cell, + 0.3, 0.03); + parallel::distributed::GridRefinement:: + refine_and_coarsen_fixed_fraction (triangulation, + estimated_error_per_cell, + 0.3, 0.03); + triangulation.execute_coarsening_and_refinement (); + + if (myid==0) + deallog << "n_global_active_cells=" + << triangulation.n_global_active_cells() + << std::endl; + + if (myid==0) + deallog << "OK" << std::endl; +} + + +int main(int argc, char *argv[]) +{ +#ifdef DEAL_II_COMPILER_SUPPORTS_MPI + MPI_Init (&argc,&argv); +#else + (void)argc; + (void)argv; +#endif + + if (Utilities::System::get_this_mpi_process (MPI_COMM_WORLD) == 0) + { + std::ofstream logfile(output_file_for_mpi("refine_and_coarsen_fixed_number_06").c_str()); + deallog.attach(logfile); + deallog.depth_console(0); + deallog.threshold_double(1.e-10); + + test<2>(); + } + else + test<2>(); + +#ifdef DEAL_II_COMPILER_SUPPORTS_MPI + MPI_Finalize(); +#endif +} diff --git a/tests/mpi/refine_and_coarsen_fixed_number_06/ncpu_1/cmp/generic b/tests/mpi/refine_and_coarsen_fixed_number_06/ncpu_1/cmp/generic new file mode 100644 index 0000000000..76a3add384 --- /dev/null +++ b/tests/mpi/refine_and_coarsen_fixed_number_06/ncpu_1/cmp/generic @@ -0,0 +1,4 @@ + +DEAL::Running on 1 CPU(s). +DEAL::2nd try +DEAL::done diff --git a/tests/mpi/refine_and_coarsen_fixed_number_06/ncpu_10/cmp/generic b/tests/mpi/refine_and_coarsen_fixed_number_06/ncpu_10/cmp/generic new file mode 100644 index 0000000000..729abc0fba --- /dev/null +++ b/tests/mpi/refine_and_coarsen_fixed_number_06/ncpu_10/cmp/generic @@ -0,0 +1,4 @@ + +DEAL::Running on 10 CPU(s). +DEAL::2nd try +DEAL::done diff --git a/tests/mpi/refine_and_coarsen_fixed_number_06/ncpu_4/cmp/generic b/tests/mpi/refine_and_coarsen_fixed_number_06/ncpu_4/cmp/generic new file mode 100644 index 0000000000..25e315071a --- /dev/null +++ b/tests/mpi/refine_and_coarsen_fixed_number_06/ncpu_4/cmp/generic @@ -0,0 +1,4 @@ + +DEAL::Running on 4 CPU(s). +DEAL::2nd try +DEAL::done -- 2.39.5