From: Denis Davydov Date: Tue, 27 Mar 2018 11:03:05 +0000 (+0200) Subject: la::p::Vector: split compress() and update_ghosts() calls into chunks X-Git-Tag: v9.0.0-rc1~273^2~1 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=19a74297a3fe845dec58b0e2891857fbe46ca136;p=dealii.git la::p::Vector: split compress() and update_ghosts() calls into chunks --- diff --git a/include/deal.II/lac/la_parallel_block_vector.h b/include/deal.II/lac/la_parallel_block_vector.h index 999b28a14f..4c6aa0ff13 100644 --- a/include/deal.II/lac/la_parallel_block_vector.h +++ b/include/deal.II/lac/la_parallel_block_vector.h @@ -81,7 +81,13 @@ namespace LinearAlgebra class BlockVector : public BlockVectorBase >, public VectorSpaceVector { + /** + * The chunks size to split communication in update_ghost_values() + * and compress() calls. + */ + static constexpr unsigned int communication_block_size = 20; public: + /** * Typedef the base class for simpler access to its own typedefs. */ diff --git a/include/deal.II/lac/la_parallel_block_vector.templates.h b/include/deal.II/lac/la_parallel_block_vector.templates.h index df6ad24706..990c0a75cd 100644 --- a/include/deal.II/lac/la_parallel_block_vector.templates.h +++ b/include/deal.II/lac/la_parallel_block_vector.templates.h @@ -241,15 +241,22 @@ namespace LinearAlgebra void BlockVector::compress (::dealii::VectorOperation::values operation) { - // start all requests for all blocks before finishing the transfers as - // this saves repeated synchronizations. In order to avoid conflict with - // possible other ongoing communication requests (from - // LA::distributed::Vector that supports unfinished requests), add an - // arbitrary number 8273 to the communication tag - for (unsigned int block=0; blockn_blocks(); ++block) - this->block(block).compress_start(block + 8273, operation); - for (unsigned int block=0; blockn_blocks(); ++block) - this->block(block).compress_finish(operation); + const unsigned int n_chunks = (this->n_blocks()+communication_block_size-1)/communication_block_size; + for (unsigned int c = 0; c < n_chunks; ++c) + { + const unsigned int start = c*communication_block_size; + const unsigned int end = std::min(start+communication_block_size, this->n_blocks()); + + // start all requests for all blocks before finishing the transfers as + // this saves repeated synchronizations. In order to avoid conflict with + // possible other ongoing communication requests (from + // LA::distributed::Vector that supports unfinished requests), add an + // arbitrary number 8273 to the communication tag + for (unsigned int block=start; blockblock(block).compress_start(block + 8273 - start, operation); + for (unsigned int block=start; blockblock(block).compress_finish(operation); + } } @@ -258,13 +265,20 @@ namespace LinearAlgebra void BlockVector::update_ghost_values () const { - // In order to avoid conflict with possible other ongoing communication - // requests (from LA::distributed::Vector that supports unfinished - // requests), add an arbitrary number 9923 to the communication tag - for (unsigned int block=0; blockn_blocks(); ++block) - this->block(block).update_ghost_values_start(block + 9923); - for (unsigned int block=0; blockn_blocks(); ++block) - this->block(block).update_ghost_values_finish(); + const unsigned int n_chunks = (this->n_blocks()+communication_block_size-1)/communication_block_size; + for (unsigned int c = 0; c < n_chunks; ++c) + { + const unsigned int start = c*communication_block_size; + const unsigned int end = std::min(start+communication_block_size, this->n_blocks()); + + // In order to avoid conflict with possible other ongoing communication + // requests (from LA::distributed::Vector that supports unfinished + // requests), add an arbitrary number 9923 to the communication tag + for (unsigned int block=start; blockblock(block).update_ghost_values_start(block - start + 9923); + for (unsigned int block=start; blockblock(block).update_ghost_values_finish(); + } } diff --git a/tests/mpi/parallel_block_vector_13.cc b/tests/mpi/parallel_block_vector_13.cc new file mode 100644 index 0000000000..7be0ded44f --- /dev/null +++ b/tests/mpi/parallel_block_vector_13.cc @@ -0,0 +1,176 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2011 - 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. +// +// --------------------------------------------------------------------- + + +// check ghost handling on parallel block vectors for large +// number of blocks with split compres_start()/update_ghosts_start(). +// almost copy-paste of parallel_block_vector_02.cc + +#include "../tests.h" +#include +#include +#include +#include +#include +#include + + +void test () +{ + unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD); + unsigned int numproc = Utilities::MPI::n_mpi_processes (MPI_COMM_WORLD); + + if (myid==0) deallog << "numproc=" << numproc << std::endl; + + + // each processor from processor 1 to 8 owns 2 indices (the other processors + // do not own any dof), and all processors are ghosting element 1 + IndexSet local_owned(std::min(16U,numproc*2)); + if (myid < 8) + local_owned.add_range(myid*2,myid*2+2); + IndexSet local_relevant(numproc*2); + local_relevant = local_owned; + local_relevant.add_range(1,2); + + LinearAlgebra::distributed::Vector v(local_owned, local_relevant, + MPI_COMM_WORLD); + + // set local values + if (myid < 8) + { + v(myid*2)=myid*2.0; + v(myid*2+1)=myid*2.0+1.0; + } + v.compress(VectorOperation::insert); + + const unsigned int n_blocks = 107; + + LinearAlgebra::distributed::BlockVector w(n_blocks); + for (unsigned int i=0; i x(w); + Assert(x.has_ghost_elements() == false, ExcInternalError()); + for (unsigned int i=0; i should not have ghosts enabled + x = w; + Assert(x.has_ghost_elements() == false, ExcInternalError()); + for (unsigned int i=0; i