From: Denis Davydov Date: Fri, 11 Mar 2016 22:27:46 +0000 (+0100) Subject: fix a bug where MPI_LOGICAL fortran type was used X-Git-Tag: v8.5.0-rc1~1223^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F2343%2Fhead;p=dealii.git fix a bug where MPI_LOGICAL fortran type was used add a simple trilinos test with vec.reinit() which previously triggered MPI_LOGICAL-related error in IndexSet class. --- diff --git a/source/base/index_set.cc b/source/base/index_set.cc index 5850d56eb0..f0045e83b8 100644 --- a/source/base/index_set.cc +++ b/source/base/index_set.cc @@ -513,14 +513,8 @@ IndexSet::make_trilinos_map (const MPI_Comm &communicator, // Check that all the processors have a contiguous range of values. Otherwise, // we risk to call different Epetra_Map on different processors and the code // hangs. -#ifdef DEAL_II_WITH_MPI - bool all_contiguous = is_contiguous(); - MPI_Allreduce(MPI_IN_PLACE, &all_contiguous, 1, MPI_LOGICAL, MPI_LAND, - communicator); - if ((all_contiguous == true) && (!overlapping)) -#else - if ((is_contiguous() == true) && (!overlapping)) -#endif + bool all_contiguous = (Utilities::MPI::min (is_contiguous() ? 1 : 0, communicator) == 1); + if ((all_contiguous) && (!overlapping)) return Epetra_Map (TrilinosWrappers::types::int_type(size()), TrilinosWrappers::types::int_type(n_elements()), 0, diff --git a/tests/trilinos/vector_reinit.cc b/tests/trilinos/vector_reinit.cc new file mode 100644 index 0000000000..178c58eeb1 --- /dev/null +++ b/tests/trilinos/vector_reinit.cc @@ -0,0 +1,77 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2004 - 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. +// +// --------------------------------------------------------------------- + +// run TrilinosWrappers::MPI::Vector.renit() in MPI environment +// to test IndexSet::make_trilinos_map(), which contained an MPI-related bug. +// Namely, a Fortran type MPI_LOGICAL was used instead of MPI_INT. + + +#include "../tests.h" + +#include +#include +#include +#include +#include +#include + +#include + +#include +#include + +using namespace dealii; + +static const unsigned int dim = 2; + +void test () +{ + MPI_Comm mpi_communicator (MPI_COMM_WORLD); + + parallel::distributed::Triangulation tria(mpi_communicator, + typename Triangulation::MeshSmoothing + (Triangulation::smoothing_on_refinement | + Triangulation::smoothing_on_coarsening)); + + GridGenerator::hyper_cube (tria, -1,0); + tria.refine_global (2); + + const unsigned int poly_degree = 1; + FE_Q fe(poly_degree); + + DoFHandler dof_handler(tria); + dof_handler.distribute_dofs(fe); + + IndexSet locally_owned_dofs = dof_handler.locally_owned_dofs (); + + TrilinosWrappers::MPI::Vector vector_Re; + vector_Re.reinit(locally_owned_dofs, mpi_communicator); + + deallog << "OK" << std::endl; +} + + +int main (int argc, char *argv[]) +{ + std::ofstream logfile("output"); + deallog.attach(logfile); + deallog.depth_console(0); + deallog.threshold_double(1.e-10); + + Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1); + { + test (); + } +} diff --git a/tests/trilinos/vector_reinit.mpirun=3.output b/tests/trilinos/vector_reinit.mpirun=3.output new file mode 100644 index 0000000000..0fd8fc12f0 --- /dev/null +++ b/tests/trilinos/vector_reinit.mpirun=3.output @@ -0,0 +1,2 @@ + +DEAL::OK