From: Peter Munch Date: Sat, 21 Mar 2020 13:32:41 +0000 (+0100) Subject: Add empty partitioner X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=be3b63e12ac350231650ecf6f9bb535df7622137;p=dealii.git Add empty partitioner --- diff --git a/include/deal.II/lac/la_sm_partitioner.h b/include/deal.II/lac/la_sm_partitioner.h new file mode 100644 index 0000000000..5245b3198f --- /dev/null +++ b/include/deal.II/lac/la_sm_partitioner.h @@ -0,0 +1,67 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2020 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. +// +// --------------------------------------------------------------------- + +#ifndef dealii_la_sm_partitioner_h +#define dealii_la_sm_partitioner_h + + +#include + +#include + +DEAL_II_NAMESPACE_OPEN + + +namespace LinearAlgebra +{ + namespace SharedMPI + { + template + class Partitioner : public LinearAlgebra::CommunicationPatternBase + { + public: + Partitioner(const MPI_Comm &comm, const MPI_Comm &comm_sm) + : comm(comm) + , comm_sm(comm_sm) + {} + + const MPI_Comm & + get_mpi_communicator() const override + { + return comm; + } + + void + reinit(const IndexSet &indexset_has, + const IndexSet &indexset_want, + const MPI_Comm &communicator) override + { + (void)indexset_has; + (void)indexset_want; + (void)communicator; + } + + private: + const MPI_Comm &comm; + const MPI_Comm &comm_sm; + }; + + } // end of namespace SharedMPI +} // end of namespace LinearAlgebra + + +DEAL_II_NAMESPACE_CLOSE + +#endif \ No newline at end of file diff --git a/include/deal.II/lac/la_sm_vector.h b/include/deal.II/lac/la_sm_vector.h index 748b73db31..5a60b1dfab 100644 --- a/include/deal.II/lac/la_sm_vector.h +++ b/include/deal.II/lac/la_sm_vector.h @@ -24,6 +24,7 @@ #include #include +#include #include #include #include @@ -157,7 +158,8 @@ namespace LinearAlgebra void reinit( - const std::shared_ptr &partitioner); + const std::shared_ptr &partitioner, + const std::shared_ptr> &partitioner_sm); void swap(Vector &v); @@ -475,6 +477,8 @@ namespace LinearAlgebra std::shared_ptr partitioner; + std::shared_ptr> partitioner_sm; + size_type allocated_size; mutable MemorySpaceData data; diff --git a/include/deal.II/lac/la_sm_vector.templates.h b/include/deal.II/lac/la_sm_vector.templates.h index fd4ca5de3b..2e46004f13 100644 --- a/include/deal.II/lac/la_sm_vector.templates.h +++ b/include/deal.II/lac/la_sm_vector.templates.h @@ -13,8 +13,8 @@ // // --------------------------------------------------------------------- -#ifndef dealii_la_parallel_vector_templates_h -#define dealii_la_parallel_vector_templates_h +#ifndef dealii_la_sm_vector_templates_h +#define dealii_la_sm_vector_templates_h #include @@ -170,14 +170,16 @@ namespace LinearAlgebra template void Vector::reinit( - const std::shared_ptr &partitioner_in) + const std::shared_ptr &partitioner, + const std::shared_ptr> & partitioner_sm) { clear_mpi_requests(); - partitioner = partitioner_in; + this->partitioner = partitioner; + this->partitioner_sm = partitioner_sm; // set vector size and allocate memory const size_type new_allocated_size = - partitioner->local_size() + partitioner->n_ghost_indices(); + this->partitioner->local_size() + this->partitioner->n_ghost_indices(); resize_val(new_allocated_size); // initialize to zero diff --git a/include/deal.II/matrix_free/matrix_free.h b/include/deal.II/matrix_free/matrix_free.h index 082d549213..a3d9078ea1 100644 --- a/include/deal.II/matrix_free/matrix_free.h +++ b/include/deal.II/matrix_free/matrix_free.h @@ -40,6 +40,7 @@ #include #include #include +#include #include #include @@ -1481,6 +1482,7 @@ public: template void initialize_dof_vector(LinearAlgebra::SharedMPI::Vector &vec, + const MPI_Comm comm_sm, const unsigned int dof_handler_index = 0) const; /** @@ -2090,6 +2092,14 @@ private: */ std::vector dof_info; + // TODO: move it to DoFInfo + mutable std::map< + unsigned int, + std::map< + const MPI_Comm, + std::shared_ptr>>> + partitioner_sm; + /** * Contains the weights for constraints stored in DoFInfo. Filled into a * separate field since several vector components might share similar @@ -2221,10 +2231,17 @@ template inline void MatrixFree::initialize_dof_vector( LinearAlgebra::SharedMPI::Vector &vec, + const MPI_Comm comm_sm, const unsigned int comp) const { AssertIndexRange(comp, n_components()); - vec.reinit(dof_info[comp].vector_partitioner); + + if (partitioner_sm[comp][comm_sm] == nullptr) + partitioner_sm[comp][comm_sm] = + std::make_shared>( + dof_info[comp].vector_partitioner->get_communicator(), comm_sm); + + vec.reinit(dof_info[comp].vector_partitioner, partitioner_sm[comp][comm_sm]); } diff --git a/tests/sm/vector_sm_01.cc b/tests/sm/vector_sm_01.cc index 749c8d033d..17d6097391 100644 --- a/tests/sm/vector_sm_01.cc +++ b/tests/sm/vector_sm_01.cc @@ -60,8 +60,10 @@ test(const int n_refinements, const int degree) dealii::MatrixFree matrix_free; matrix_free.reinit(mapping, dof_handler, constraint, quad, additional_data); + MPI_Comm comm_sm; + LinearAlgebra::SharedMPI::Vector vec; - matrix_free.initialize_dof_vector(vec); + matrix_free.initialize_dof_vector(vec, comm_sm); } int