From: Peter Munch Date: Wed, 21 Oct 2020 20:38:48 +0000 (+0200) Subject: Wrap partitioners within MatrixFree X-Git-Tag: v9.3.0-rc1~976^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=5e3eec386538d6792bdd6cf7a232c2564cdb2e6c;p=dealii.git Wrap partitioners within MatrixFree --- diff --git a/include/deal.II/base/memory_space_data.h b/include/deal.II/base/memory_space_data.h index da957755d6..e2feb375db 100644 --- a/include/deal.II/base/memory_space_data.h +++ b/include/deal.II/base/memory_space_data.h @@ -75,6 +75,11 @@ namespace MemorySpace * Pointer to data on the device. */ std::unique_ptr values_dev; + + /** + * Pointers to the data of the processes sharing the same memory. + */ + std::vector> values_sm; }; @@ -118,6 +123,8 @@ namespace MemorySpace // This is not used but it allows to simplify the code until we start using // CUDA-aware MPI. std::unique_ptr values_dev; + + std::vector> values_sm; }; @@ -165,6 +172,11 @@ namespace MemorySpace std::unique_ptr> values; std::unique_ptr values_dev; + + /** + * This is currently not used. + */ + std::vector> values_sm; }; diff --git a/include/deal.II/lac/la_parallel_vector.h b/include/deal.II/lac/la_parallel_vector.h index 16e1251a2d..f21263188b 100644 --- a/include/deal.II/lac/la_parallel_vector.h +++ b/include/deal.II/lac/la_parallel_vector.h @@ -1120,6 +1120,13 @@ namespace LinearAlgebra void set_ghost_state(const bool ghosted) const; + /** + * Get pointers to the beginning of the values of the other + * processes of the same shared-memory domain. + */ + const std::vector> & + shared_vector_data() const; + //@} /** @@ -1497,6 +1504,15 @@ namespace LinearAlgebra + template + const std::vector> & + Vector::shared_vector_data() const + { + return data.values_sm; + } + + + template inline Number Vector::operator()(const size_type global_index) const diff --git a/include/deal.II/matrix_free/dof_info.h b/include/deal.II/matrix_free/dof_info.h index 36002f9a90..551d6fe4e6 100644 --- a/include/deal.II/matrix_free/dof_info.h +++ b/include/deal.II/matrix_free/dof_info.h @@ -33,6 +33,7 @@ #include #include #include +#include #include #include @@ -517,7 +518,15 @@ namespace internal std::shared_ptr vector_partitioner; /** - * This partitioning selects a subset of ghost indices to the full + * Vector exchanger compatible with vector_partitioner. + */ + std::shared_ptr< + const internal::MatrixFreeFunctions::VectorDataExchange::Base> + vector_exchanger; + + /** + * Vector exchanger compatible with partitioners that select a subset of + * ghost indices to the full * vector partitioner stored in @p vector_partitioner. These * partitioners are used in specialized loops that only import parts of * the ghosted region for reducing the amount of communication. There @@ -534,8 +543,11 @@ namespace internal * values and the gradients on all faces adjacent to the locally owned * cells. */ - std::array, 5> - vector_partitioner_face_variants; + std::array< + std::shared_ptr< + const internal::MatrixFreeFunctions::VectorDataExchange::Base>, + 5> + vector_exchanger_face_variants; /** * This stores a (sorted) list of all locally owned degrees of freedom diff --git a/include/deal.II/matrix_free/dof_info.templates.h b/include/deal.II/matrix_free/dof_info.templates.h index c94c42819c..7d7c189904 100644 --- a/include/deal.II/matrix_free/dof_info.templates.h +++ b/include/deal.II/matrix_free/dof_info.templates.h @@ -496,6 +496,10 @@ namespace internal Utilities::MPI::Partitioner *vec_part = const_cast(vector_partitioner.get()); vec_part->set_ghost_indices(ghost_indices); + + vector_exchanger = std::make_shared< + internal::MatrixFreeFunctions::VectorDataExchange::PartitionerWrapper>( + vector_partitioner); } @@ -1209,17 +1213,22 @@ namespace internal Utilities::MPI::min(compressed_set.n_elements() == part.ghost_indices().n_elements(), part.get_mpi_communicator()) != 0; + + std::shared_ptr temp_0; + if (all_ghosts_equal) - vector_partitioner_face_variants[0] = vector_partitioner; + temp_0 = vector_partitioner; else { - vector_partitioner_face_variants[0] = - std::make_shared( - part.locally_owned_range(), part.get_mpi_communicator()); - const_cast( - vector_partitioner_face_variants[0].get()) + temp_0 = std::make_shared( + part.locally_owned_range(), part.get_mpi_communicator()); + const_cast(temp_0.get()) ->set_ghost_indices(compressed_set, part.ghost_indices()); } + + vector_exchanger_face_variants[0] = + std::make_shared(temp_0); } // construct a numbering of faces @@ -1433,34 +1442,43 @@ namespace internal } }; + std::shared_ptr temp_1, temp_2, temp_3, + temp_4; + // partitioner 1: values on faces - process_values(vector_partitioner_face_variants[1], loop_over_faces); + process_values(temp_1, loop_over_faces); // partitioner 2: values and gradients on faces - process_gradients(vector_partitioner_face_variants[1], - vector_partitioner_face_variants[2], - loop_over_faces); + process_gradients(temp_1, temp_2, loop_over_faces); if (fill_cell_centric) { ghost_indices.clear(); // partitioner 3: values on all faces - process_values(vector_partitioner_face_variants[3], - loop_over_all_faces); + process_values(temp_3, loop_over_all_faces); // partitioner 4: values and gradients on faces - process_gradients(vector_partitioner_face_variants[3], - vector_partitioner_face_variants[4], - loop_over_all_faces); + process_gradients(temp_3, temp_4, loop_over_all_faces); } else { - vector_partitioner_face_variants[3] = - std::make_shared( - part.locally_owned_range(), part.get_mpi_communicator()); - vector_partitioner_face_variants[4] = - std::make_shared( - part.locally_owned_range(), part.get_mpi_communicator()); + temp_3 = std::make_shared( + part.locally_owned_range(), part.get_mpi_communicator()); + temp_4 = std::make_shared( + part.locally_owned_range(), part.get_mpi_communicator()); } + + vector_exchanger_face_variants[1] = std::make_shared< + internal::MatrixFreeFunctions::VectorDataExchange::PartitionerWrapper>( + temp_1); + vector_exchanger_face_variants[2] = std::make_shared< + internal::MatrixFreeFunctions::VectorDataExchange::PartitionerWrapper>( + temp_2); + vector_exchanger_face_variants[3] = std::make_shared< + internal::MatrixFreeFunctions::VectorDataExchange::PartitionerWrapper>( + temp_3); + vector_exchanger_face_variants[4] = std::make_shared< + internal::MatrixFreeFunctions::VectorDataExchange::PartitionerWrapper>( + temp_4); } diff --git a/include/deal.II/matrix_free/matrix_free.h b/include/deal.II/matrix_free/matrix_free.h index ccee8dbd8a..9240fa77b6 100644 --- a/include/deal.II/matrix_free/matrix_free.h +++ b/include/deal.II/matrix_free/matrix_free.h @@ -3108,7 +3108,7 @@ namespace internal DataAccessOnFaces::unspecified) for (unsigned int c = 0; c < matrix_free.n_components(); ++c) AssertDimension( - matrix_free.get_dof_info(c).vector_partitioner_face_variants.size(), + matrix_free.get_dof_info(c).vector_exchanger_face_variants.size(), 5); } @@ -3161,39 +3161,39 @@ namespace internal * Get partitioner for the given @p mf_component taking into * account vector_face_access set in constructor. */ - const Utilities::MPI::Partitioner & + const internal::MatrixFreeFunctions::VectorDataExchange::Base & get_partitioner(const unsigned int mf_component) const { AssertDimension(matrix_free.get_dof_info(mf_component) - .vector_partitioner_face_variants.size(), + .vector_exchanger_face_variants.size(), 5); if (vector_face_access == dealii::MatrixFree:: DataAccessOnFaces::none) return *matrix_free.get_dof_info(mf_component) - .vector_partitioner_face_variants[0]; + .vector_exchanger_face_variants[0]; else if (vector_face_access == dealii::MatrixFree:: DataAccessOnFaces::values) return *matrix_free.get_dof_info(mf_component) - .vector_partitioner_face_variants[1]; + .vector_exchanger_face_variants[1]; else if (vector_face_access == dealii::MatrixFree:: DataAccessOnFaces::gradients) return *matrix_free.get_dof_info(mf_component) - .vector_partitioner_face_variants[2]; + .vector_exchanger_face_variants[2]; else if (vector_face_access == dealii::MatrixFree:: DataAccessOnFaces::values_all_faces) return *matrix_free.get_dof_info(mf_component) - .vector_partitioner_face_variants[3]; + .vector_exchanger_face_variants[3]; else if (vector_face_access == dealii::MatrixFree:: DataAccessOnFaces::gradients_all_faces) return *matrix_free.get_dof_info(mf_component) - .vector_partitioner_face_variants[4]; + .vector_exchanger_face_variants[4]; else - return *matrix_free.get_dof_info(mf_component).vector_partitioner.get(); + return *matrix_free.get_dof_info(mf_component).vector_exchanger.get(); } @@ -3299,12 +3299,13 @@ namespace internal part.export_to_ghosted_array_start( component_in_block_vector + channel_shift, ArrayView(vec.begin(), part.local_size()), - ArrayView(tmp_data[component_in_block_vector]->begin(), - part.n_import_indices()), + vec.shared_vector_data(), ArrayView(const_cast(vec.begin()) + part.local_size(), matrix_free.get_dof_info(mf_component) .vector_partitioner->n_ghost_indices()), + ArrayView(tmp_data[component_in_block_vector]->begin(), + part.n_import_indices()), this->requests[component_in_block_vector]); # endif } @@ -3377,9 +3378,12 @@ namespace internal const auto &part = get_partitioner(mf_component); - if (part.n_ghost_indices() != 0 || part.n_import_indices() != 0) + if (part.n_ghost_indices() != 0 || part.n_import_indices() != 0 || + part.n_import_sm_procs() != 0) { part.export_to_ghosted_array_finish( + ArrayView(vec.begin(), part.local_size()), + vec.shared_vector_data(), ArrayView(const_cast(vec.begin()) + part.local_size(), matrix_free.get_dof_info(mf_component) @@ -3480,7 +3484,8 @@ namespace internal const auto &part = get_partitioner(mf_component); - if (part.n_ghost_indices() == 0 && part.n_import_indices() == 0) + if (part.n_ghost_indices() == 0 && part.n_import_indices() == 0 && + part.n_import_sm_procs() == 0) return; tmp_data[component_in_block_vector] = @@ -3492,6 +3497,8 @@ namespace internal part.import_from_ghosted_array_start( dealii::VectorOperation::add, component_in_block_vector + channel_shift, + ArrayView(vec.begin(), part.local_size()), + vec.shared_vector_data(), ArrayView(vec.begin() + part.local_size(), matrix_free.get_dof_info(mf_component) .vector_partitioner->n_ghost_indices()), @@ -3557,7 +3564,6 @@ namespace internal std::is_same::value, "Type mismatch between VectorType and VectorDataExchange"); (void)component_in_block_vector; - if (vec.size() != 0) { # ifdef DEAL_II_WITH_MPI @@ -3568,23 +3574,25 @@ namespace internal const auto &part = get_partitioner(mf_component); - if (part.n_ghost_indices() == 0 && part.n_import_indices() == 0) - return; - - part.import_from_ghosted_array_finish( - VectorOperation::add, - ArrayView( - tmp_data[component_in_block_vector]->begin(), - part.n_import_indices()), - ArrayView(vec.begin(), part.local_size()), - ArrayView(vec.begin() + part.local_size(), - matrix_free.get_dof_info(mf_component) - .vector_partitioner->n_ghost_indices()), - this->requests[component_in_block_vector]); + if (part.n_ghost_indices() != 0 || part.n_import_indices() != 0 || + part.n_import_sm_procs() != 0) + { + part.import_from_ghosted_array_finish( + VectorOperation::add, + ArrayView(vec.begin(), part.local_size()), + vec.shared_vector_data(), + ArrayView(vec.begin() + part.local_size(), + matrix_free.get_dof_info(mf_component) + .vector_partitioner->n_ghost_indices()), + ArrayView( + tmp_data[component_in_block_vector]->begin(), + part.n_import_indices()), + this->requests[component_in_block_vector]); - matrix_free.release_scratch_data_non_threadsafe( - tmp_data[component_in_block_vector]); - tmp_data[component_in_block_vector] = nullptr; + matrix_free.release_scratch_data_non_threadsafe( + tmp_data[component_in_block_vector]); + tmp_data[component_in_block_vector] = nullptr; + } # endif } } @@ -3640,43 +3648,29 @@ namespace internal if (ghosts_were_set == true) return; - if (vector_face_access == - dealii::MatrixFree:: - DataAccessOnFaces::unspecified || - vec.size() == 0) - vec.zero_out_ghosts(); - else + if (vec.size() != 0) { # ifdef DEAL_II_WITH_MPI AssertDimension(requests.size(), tmp_data.size()); const unsigned int mf_component = find_vector_in_mf(vec); - const Utilities::MPI::Partitioner &part = - get_partitioner(mf_component); - if (&part == - matrix_free.get_dof_info(mf_component).vector_partitioner.get()) - vec.zero_out_ghosts(); - else if (part.n_ghost_indices() > 0) + + const auto &part = get_partitioner(mf_component); + + if (part.n_ghost_indices() > 0) { - for (std::vector>:: - const_iterator my_ghosts = - part.ghost_indices_within_larger_ghost_set().begin(); - my_ghosts != - part.ghost_indices_within_larger_ghost_set().end(); - ++my_ghosts) - for (unsigned int j = my_ghosts->first; j < my_ghosts->second; - j++) - { - const_cast &>( - vec) - .local_element(j + part.local_size()) = 0.; - } + part.reset_ghost_values(ArrayView( + const_cast &>(vec) + .begin() + + part.local_size(), + matrix_free.get_dof_info(mf_component) + .vector_partitioner->n_ghost_indices())); } - // let vector know that it's not ghosted anymore - vec.set_ghost_state(false); # endif } + // let vector know that it's not ghosted anymore + vec.set_ghost_state(false); } diff --git a/include/deal.II/matrix_free/matrix_free.templates.h b/include/deal.II/matrix_free/matrix_free.templates.h index b6152e0076..39651508a7 100644 --- a/include/deal.II/matrix_free/matrix_free.templates.h +++ b/include/deal.II/matrix_free/matrix_free.templates.h @@ -812,6 +812,11 @@ namespace internal std::make_shared(locally_owned_dofs[no], task_info.communicator); + dof_info[no].vector_exchanger = + std::make_shared( + dof_info[no].vector_partitioner); + // initialize the arrays for indices const unsigned int n_components_total = dof_info[no].start_components.back(); diff --git a/include/deal.II/matrix_free/vector_data_exchange.h b/include/deal.II/matrix_free/vector_data_exchange.h new file mode 100644 index 0000000000..38d11aafa9 --- /dev/null +++ b/include/deal.II/matrix_free/vector_data_exchange.h @@ -0,0 +1,256 @@ +// --------------------------------------------------------------------- +// +// 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_matrix_free_vector_data_exchange_h +#define dealii_matrix_free_vector_data_exchange_h + + +#include + +#include + + +DEAL_II_NAMESPACE_OPEN + +namespace internal +{ + namespace MatrixFreeFunctions + { + /** + * Namespace containing classes for inter-process data exchange (i.e., + * for update_ghost_values and compress) in MatrixFree. + */ + namespace VectorDataExchange + { + /** + * Interface needed by MatrixFree. + */ + class Base + { + public: + virtual ~Base() = default; + + virtual unsigned int + local_size() const = 0; + + virtual unsigned int + n_ghost_indices() const = 0; + + virtual unsigned int + n_import_indices() const = 0; + + virtual unsigned int + n_import_sm_procs() const = 0; + + virtual types::global_dof_index + size() const = 0; + + virtual void + export_to_ghosted_array_start( + const unsigned int communication_channel, + const ArrayView & locally_owned_array, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + const ArrayView & temporary_storage, + std::vector & requests) const = 0; + + virtual void + export_to_ghosted_array_finish( + const ArrayView & locally_owned_array, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + std::vector & requests) const = 0; + + virtual void + import_from_ghosted_array_start( + const VectorOperation::values vector_operation, + const unsigned int communication_channel, + const ArrayView & locally_owned_array, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + const ArrayView & temporary_storage, + std::vector & requests) const = 0; + + virtual void + import_from_ghosted_array_finish( + const VectorOperation::values vector_operation, + const ArrayView & locally_owned_storage, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + const ArrayView & temporary_storage, + std::vector & requests) const = 0; + + virtual void + reset_ghost_values(const ArrayView &ghost_array) const = 0; + + virtual void + export_to_ghosted_array_start( + const unsigned int communication_channel, + const ArrayView & locally_owned_array, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + const ArrayView & temporary_storage, + std::vector & requests) const = 0; + + virtual void + export_to_ghosted_array_finish( + const ArrayView & locally_owned_array, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + std::vector & requests) const = 0; + + virtual void + import_from_ghosted_array_start( + const VectorOperation::values vector_operation, + const unsigned int communication_channel, + const ArrayView & locally_owned_array, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + const ArrayView & temporary_storage, + std::vector & requests) const = 0; + + virtual void + import_from_ghosted_array_finish( + const VectorOperation::values vector_operation, + const ArrayView & locally_owned_storage, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + const ArrayView & temporary_storage, + std::vector & requests) const = 0; + + virtual void + reset_ghost_values(const ArrayView &ghost_array) const = 0; + }; + + + /** + * Class that simply delegates the task to a Utilities::MPI::Partitioner. + */ + class PartitionerWrapper : public Base + { + public: + PartitionerWrapper( + const std::shared_ptr + &partitioner); + + virtual ~PartitionerWrapper() = default; + + unsigned int + local_size() const override; + + unsigned int + n_ghost_indices() const override; + + unsigned int + n_import_indices() const override; + + unsigned int + n_import_sm_procs() const override; + + types::global_dof_index + size() const override; + + void + export_to_ghosted_array_start( + const unsigned int communication_channel, + const ArrayView & locally_owned_array, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + const ArrayView & temporary_storage, + std::vector & requests) const override; + + void + export_to_ghosted_array_finish( + const ArrayView & locally_owned_array, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + std::vector & requests) const override; + + void + import_from_ghosted_array_start( + const VectorOperation::values vector_operation, + const unsigned int communication_channel, + const ArrayView & locally_owned_array, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + const ArrayView & temporary_storage, + std::vector & requests) const override; + + void + import_from_ghosted_array_finish( + const VectorOperation::values vector_operation, + const ArrayView & locally_owned_storage, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + const ArrayView & temporary_storage, + std::vector & requests) const override; + + void + reset_ghost_values(const ArrayView &ghost_array) const override; + + void + export_to_ghosted_array_start( + const unsigned int communication_channel, + const ArrayView & locally_owned_array, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + const ArrayView & temporary_storage, + std::vector & requests) const override; + + void + export_to_ghosted_array_finish( + const ArrayView & locally_owned_array, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + std::vector & requests) const override; + + void + import_from_ghosted_array_start( + const VectorOperation::values vector_operation, + const unsigned int communication_channel, + const ArrayView & locally_owned_array, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + const ArrayView & temporary_storage, + std::vector & requests) const override; + + void + import_from_ghosted_array_finish( + const VectorOperation::values vector_operation, + const ArrayView & locally_owned_storage, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + const ArrayView & temporary_storage, + std::vector & requests) const override; + + void + reset_ghost_values(const ArrayView &ghost_array) const override; + + private: + template + void + reset_ghost_values_impl(const ArrayView &ghost_array) const; + + const std::shared_ptr partitioner; + }; + + } // namespace VectorDataExchange + } // end of namespace MatrixFreeFunctions +} // end of namespace internal + +DEAL_II_NAMESPACE_CLOSE + +#endif diff --git a/source/matrix_free/CMakeLists.txt b/source/matrix_free/CMakeLists.txt index 1fe2e83174..a993b8391c 100644 --- a/source/matrix_free/CMakeLists.txt +++ b/source/matrix_free/CMakeLists.txt @@ -29,6 +29,7 @@ SET(_src matrix_free.cc shape_info.cc task_info.cc + vector_data_exchange.cc ) SET(_inst diff --git a/source/matrix_free/vector_data_exchange.cc b/source/matrix_free/vector_data_exchange.cc new file mode 100644 index 0000000000..2193d30417 --- /dev/null +++ b/source/matrix_free/vector_data_exchange.cc @@ -0,0 +1,324 @@ +// --------------------------------------------------------------------- +// +// 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. +// +// --------------------------------------------------------------------- + +#include +#include +#include +#include +#include +#include + +#include + +#ifdef DEAL_II_WITH_64BIT_INDICES +# include +#endif + +#include +#include + + +DEAL_II_NAMESPACE_OPEN + +namespace internal +{ + namespace MatrixFreeFunctions + { + namespace VectorDataExchange + { + PartitionerWrapper::PartitionerWrapper( + const std::shared_ptr &partitioner) + : partitioner(partitioner) + {} + + + + unsigned int + PartitionerWrapper::local_size() const + { + return partitioner->local_size(); + } + + + + unsigned int + PartitionerWrapper::n_ghost_indices() const + { + return partitioner->n_ghost_indices(); + } + + + + unsigned int + PartitionerWrapper::n_import_indices() const + { + return partitioner->n_import_indices(); + } + + + + unsigned int + PartitionerWrapper::n_import_sm_procs() const + { + return 0; + } + + + + types::global_dof_index + PartitionerWrapper::size() const + { + return partitioner->size(); + } + + + + void + PartitionerWrapper::export_to_ghosted_array_start( + const unsigned int communication_channel, + const ArrayView & locally_owned_array, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + const ArrayView & temporary_storage, + std::vector & requests) const + { + (void)shared_arrays; +#ifndef DEAL_II_WITH_MPI + (void)communication_channel; + (void)locally_owned_array; + (void)ghost_array; + (void)temporary_storage; + (void)requests; +#else + partitioner->export_to_ghosted_array_start(communication_channel, + locally_owned_array, + temporary_storage, + ghost_array, + requests); +#endif + } + + + + void + PartitionerWrapper::export_to_ghosted_array_finish( + const ArrayView & locally_owned_array, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + std::vector & requests) const + { + (void)locally_owned_array; + (void)shared_arrays; +#ifndef DEAL_II_WITH_MPI + (void)ghost_array; + (void)requests; +#else + partitioner->export_to_ghosted_array_finish(ghost_array, requests); +#endif + } + + + + void + PartitionerWrapper::import_from_ghosted_array_start( + const VectorOperation::values vector_operation, + const unsigned int communication_channel, + const ArrayView & locally_owned_array, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + const ArrayView & temporary_storage, + std::vector & requests) const + { + (void)locally_owned_array; + (void)shared_arrays; +#ifndef DEAL_II_WITH_MPI + (void)vector_operation; + (void)communication_channel; + (void)ghost_array; + (void)temporary_storage; + (void)requests; +#else + partitioner->import_from_ghosted_array_start(vector_operation, + communication_channel, + ghost_array, + temporary_storage, + requests); +#endif + } + + + + void + PartitionerWrapper::import_from_ghosted_array_finish( + const VectorOperation::values vector_operation, + const ArrayView & locally_owned_storage, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + const ArrayView & temporary_storage, + std::vector & requests) const + { + (void)shared_arrays; +#ifndef DEAL_II_WITH_MPI + (void)vector_operation; + (void)locally_owned_storage; + (void)ghost_array; + (void)temporary_storage; + (void)requests; +#else + partitioner->import_from_ghosted_array_finish(vector_operation, + temporary_storage, + locally_owned_storage, + ghost_array, + requests); +#endif + } + + + + void + PartitionerWrapper::reset_ghost_values( + const ArrayView &ghost_array) const + { + reset_ghost_values_impl(ghost_array); + } + + + + void + PartitionerWrapper::export_to_ghosted_array_start( + const unsigned int communication_channel, + const ArrayView & locally_owned_array, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + const ArrayView & temporary_storage, + std::vector & requests) const + { + (void)shared_arrays; +#ifndef DEAL_II_WITH_MPI + (void)communication_channel; + (void)locally_owned_array; + (void)ghost_array; + (void)temporary_storage; + (void)requests; +#else + partitioner->export_to_ghosted_array_start(communication_channel, + locally_owned_array, + temporary_storage, + ghost_array, + requests); +#endif + } + + + + void + PartitionerWrapper::export_to_ghosted_array_finish( + const ArrayView & locally_owned_array, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + std::vector & requests) const + { + (void)locally_owned_array; + (void)shared_arrays; +#ifndef DEAL_II_WITH_MPI + (void)ghost_array; + (void)requests; +#else + partitioner->export_to_ghosted_array_finish(ghost_array, requests); +#endif + } + + + + void + PartitionerWrapper::import_from_ghosted_array_start( + const VectorOperation::values vector_operation, + const unsigned int communication_channel, + const ArrayView & locally_owned_array, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + const ArrayView & temporary_storage, + std::vector & requests) const + { + (void)locally_owned_array; + (void)shared_arrays; +#ifndef DEAL_II_WITH_MPI + (void)vector_operation; + (void)communication_channel; + (void)ghost_array; + (void)temporary_storage; + (void)requests; +#else + partitioner->import_from_ghosted_array_start(vector_operation, + communication_channel, + ghost_array, + temporary_storage, + requests); +#endif + } + + + + void + PartitionerWrapper::import_from_ghosted_array_finish( + const VectorOperation::values vector_operation, + const ArrayView & locally_owned_storage, + const std::vector> &shared_arrays, + const ArrayView & ghost_array, + const ArrayView & temporary_storage, + std::vector & requests) const + { + (void)shared_arrays; +#ifndef DEAL_II_WITH_MPI + (void)vector_operation; + (void)locally_owned_storage; + (void)ghost_array; + (void)temporary_storage; + (void)requests; +#else + partitioner->import_from_ghosted_array_finish(vector_operation, + temporary_storage, + locally_owned_storage, + ghost_array, + requests); +#endif + } + + + + void + PartitionerWrapper::reset_ghost_values( + const ArrayView &ghost_array) const + { + reset_ghost_values_impl(ghost_array); + } + + template + void + PartitionerWrapper::reset_ghost_values_impl( + const ArrayView &ghost_array) const + { + for (const auto &my_ghosts : + partitioner->ghost_indices_within_larger_ghost_set()) + for (unsigned int j = my_ghosts.first; j < my_ghosts.second; ++j) + ghost_array[j] = 0.; + } + + } // namespace VectorDataExchange + } // namespace MatrixFreeFunctions +} // namespace internal + + +DEAL_II_NAMESPACE_CLOSE