From: Denis Davydov Date: Sat, 24 Mar 2018 19:26:25 +0000 (+0100) Subject: improve MPI::Partitioner to support empty ArrayView for owned DoFs X-Git-Tag: v9.0.0-rc1~283^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=b469f2bc05ad9be2abb740d4783aaf0ce4ea02f6;p=dealii.git improve MPI::Partitioner to support empty ArrayView for owned DoFs --- diff --git a/include/deal.II/base/partitioner.templates.h b/include/deal.II/base/partitioner.templates.h index 5ba72ce87d..66d24d5f1f 100644 --- a/include/deal.II/base/partitioner.templates.h +++ b/include/deal.II/base/partitioner.templates.h @@ -42,7 +42,6 @@ namespace Utilities const ArrayView &ghost_array, std::vector &requests) const { - AssertDimension(locally_owned_array.size(), local_size()); AssertDimension(temporary_storage.size(), n_import_indices()); Assert(ghost_array.size() == n_ghost_indices() || ghost_array.size() == n_ghost_indices_in_larger_set, @@ -52,6 +51,9 @@ namespace Utilities const unsigned int n_import_targets = import_targets_data.size(); const unsigned int n_ghost_targets = ghost_targets_data.size(); + if (n_import_targets>0) + AssertDimension(locally_owned_array.size(), local_size()); + Assert(requests.size() == 0, ExcMessage("Another operation seems to still be running. " "Call update_ghost_values_finish() first.")); @@ -316,7 +318,6 @@ namespace Utilities const ArrayView &ghost_array, std::vector &requests) const { - AssertDimension(locally_owned_array.size(), local_size()); AssertDimension(temporary_storage.size(), n_import_indices()); Assert(ghost_array.size() == n_ghost_indices() || ghost_array.size() == n_ghost_indices_in_larger_set, @@ -353,6 +354,7 @@ namespace Utilities // first wait for the receive to complete if (requests.size() > 0 && n_import_targets > 0) { + AssertDimension(locally_owned_array.size(), local_size()); const int ierr = MPI_Waitall (n_import_targets, requests.data(), MPI_STATUSES_IGNORE); AssertThrowMPI(ierr); diff --git a/tests/mpi/parallel_partitioner_08.cc b/tests/mpi/parallel_partitioner_08.cc new file mode 100644 index 0000000000..cd0a8b3105 --- /dev/null +++ b/tests/mpi/parallel_partitioner_08.cc @@ -0,0 +1,172 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2017 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 MPI::Partitioner update_ghosts() and compress() in case we have +// empty owned DoFs + +#include "../tests.h" +#include +#include +#include +#include + +template +void test() +{ + const unsigned int rank = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD); + + // setup index sets + // subset: is2 + // ghost: 8 9 10 11 is3 + // rank 0 : 00 01 02 03 04 05 06 07 00 00 00 00 + // rank 1 : 00 00 00 00 00 00 00 00 00 00 00 00 + // ghost: 0 1 2 3 is3 + // subset: 1 2 is2 + // + // expected result udpate ghosts() + // + // rank 0 : 00 01 02 03 04 05 06 07 00 00 00 00 + // rank 1 : 00 00 00 00 00 00 00 00 00 01 02 00 + // + // compress(insert) -- does not change anything but zero ghosts + // + // set rank1 ghosts to: 00 10 20 00 + // compress(add) + // + // rank 0 : 00 11 22 03 04 05 06 07 00 00 00 00 + // rank 1 : 00 00 00 00 00 00 00 00 00 10 20 00 + + + IndexSet is1(16), is2(16), is3(16); + + if (rank==0) + { + is1.add_range(0,8); + // note: empty is2 + is3.add_range(8,12); + } + else if (rank==1) + { + is1.add_range(8,16); + is2.add_index(1); + is2.add_index(2); + is3.add_range(0,4); + } + + // create partitioner + std::shared_ptr partitioner + (new Utilities::MPI::Partitioner(is1, MPI_COMM_WORLD)); + partitioner->set_ghost_indices(is3); + std::shared_ptr tight_partitioner + (new Utilities::MPI::Partitioner(is1, MPI_COMM_WORLD)); + tight_partitioner->set_ghost_indices(is2,is3); + + // create vector + AlignedVector owned(rank==0?8:0); + AlignedVector ghost(4); + + for (unsigned int i = 0; i < 4; ++i) + ghost[i] = 0.; + + if (rank==0) + for (int i = 0; i < 8; i++) + owned[i] = i; + + // update ghost values + // vector of requests + std::vector requests; + std::vector compress_requests; + + // allocate temporal array + AlignedVector tmp_data; + tmp_data.resize_fast(tight_partitioner->n_import_indices()); + + // begin exchange, and ... + tight_partitioner->export_to_ghosted_array_start( + 0, + ArrayView(owned.begin(), owned.size()), + ArrayView(tmp_data.begin(), tight_partitioner->n_import_indices()), + ArrayView(ghost.begin(), ghost.size()), + requests); + + // ... finish exchange + tight_partitioner->export_to_ghosted_array_finish( + ArrayView(ghost.begin(), ghost.size()), + requests); + + auto print = [&] () + { + deallog << "owned:" << std::endl; + for (auto el : owned) + deallog << el << " "; + deallog << std::endl + << "ghost:" << std::endl; + for (auto el : ghost) + deallog << el << " "; + deallog << std::endl; + }; + + deallog << "update ghosts()"< import_data; + import_data.resize_fast(tight_partitioner->n_import_indices()); + + // now do insert: + auto compress = [&] (VectorOperation::values operation) + { + const unsigned int counter = 0; + tight_partitioner->import_from_ghosted_array_start + (operation, counter, + ArrayView(ghost.begin(), ghost.size()), + ArrayView(import_data.begin(), tight_partitioner->n_import_indices()), + compress_requests); + + tight_partitioner->import_from_ghosted_array_finish + (operation, + ArrayView(import_data.begin(), tight_partitioner->n_import_indices()), + ArrayView(owned.begin(), owned.size()), + ArrayView(ghost.begin(), ghost.size()), + compress_requests); + }; + + deallog << "compress(insert)"<