From: Peter Munch Date: Mon, 14 Mar 2022 20:11:11 +0000 (+0100) Subject: MGTransferGlobalCoarsening::initialize_dof_vector() skip zeroing X-Git-Tag: v9.4.0-rc1~368^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F13539%2Fhead;p=dealii.git MGTransferGlobalCoarsening::initialize_dof_vector() skip zeroing --- diff --git a/include/deal.II/multigrid/mg_transfer_global_coarsening.h b/include/deal.II/multigrid/mg_transfer_global_coarsening.h index cc6b5552a9..b72b71e879 100644 --- a/include/deal.II/multigrid/mg_transfer_global_coarsening.h +++ b/include/deal.II/multigrid/mg_transfer_global_coarsening.h @@ -597,15 +597,22 @@ public: memory_consumption() const; private: + /** + * Function to initialize internal level vectors. + */ + void + initialize_dof_vector(const unsigned int level, VectorType &vector) const; + /** * Collection of the two-level transfer operators. */ MGLevelObject> transfer; /** - * Function to initialize internal level vectors. + * External partitioners used during initialize_dof_vector(). */ - std::function initialize_dof_vector; + std::vector> + external_partitioners; }; @@ -623,7 +630,7 @@ MGTransferGlobalCoarsening::MGTransferGlobalCoarsening( &initialize_dof_vector) : transfer(transfer) { - build(initialize_dof_vector); + this->build(initialize_dof_vector); } @@ -634,28 +641,19 @@ MGTransferGlobalCoarsening::build( const std::vector> &external_partitioners) { - if (external_partitioners.size() > 0) - { - Assert( - this->initialize_dof_vector == nullptr, - ExcMessage( - "A initialize_dof_vector function has already been registered in the constructor!")); + this->external_partitioners = external_partitioners; + if (this->external_partitioners.size() > 0) + { const unsigned int min_level = transfer.min_level(); const unsigned int max_level = transfer.max_level(); - AssertDimension(external_partitioners.size(), transfer.n_levels()); - - this->initialize_dof_vector = - [external_partitioners, min_level]( - const unsigned int level, - LinearAlgebra::distributed::Vector - &vec) { vec.reinit(external_partitioners[level - min_level]); }; + AssertDimension(this->external_partitioners.size(), transfer.n_levels()); for (unsigned int l = min_level + 1; l <= max_level; ++l) transfer[l].enable_inplace_operations_if_possible( - external_partitioners[l - 1 - min_level], - external_partitioners[l - min_level]); + this->external_partitioners[l - 1 - min_level], + this->external_partitioners[l - min_level]); } } @@ -684,7 +682,7 @@ MGTransferGlobalCoarsening::build( external_partitioners[l - min_level] = vector.get_partitioner(); } - build(external_partitioners); + this->build(external_partitioners); } else { @@ -694,6 +692,24 @@ MGTransferGlobalCoarsening::build( +template +void +MGTransferGlobalCoarsening::initialize_dof_vector( + const unsigned int level, + VectorType & vec) const +{ + AssertDimension(transfer.n_levels(), external_partitioners.size()); + AssertIndexRange(level, transfer.max_level() + 1); + + const auto &external_partitioner = + external_partitioners[level - transfer.min_level()]; + + if (vec.get_partitioner().get() != external_partitioner.get()) + vec.reinit(external_partitioner); +} + + + template void MGTransferGlobalCoarsening::prolongate( @@ -741,13 +757,6 @@ MGTransferGlobalCoarsening::copy_to_mg( { (void)dof_handler; - Assert( - initialize_dof_vector, - ExcMessage( - "To be able to use this function, a function to initialize an internal " - "DoF vector has to be provided in the constructor of " - "MGTransferGlobalCoarsening.")); - for (unsigned int level = dst.min_level(); level <= dst.max_level(); ++level) { initialize_dof_vector(level, dst[level]); @@ -783,13 +792,6 @@ MGTransferGlobalCoarsening::interpolate_to_mg( MGLevelObject &dst, const InVector & src) const { - Assert( - initialize_dof_vector, - ExcMessage( - "To be able to use this function, a function to initialize an internal " - "DoF vector has to be provided in the constructor of " - "MGTransferGlobalCoarsening.")); - const unsigned int min_level = transfer.min_level(); const unsigned int max_level = transfer.max_level();