From 46b8df752960172b9c4f0339b0afed23c2648144 Mon Sep 17 00:00:00 2001 From: Timo Heister Date: Thu, 14 Nov 2024 06:32:51 -0500 Subject: [PATCH] enable MemorySpace in GMG transfer --- cmake/config/template-arguments.in | 9 +- include/deal.II/multigrid/mg_transfer.h | 47 ++--- .../deal.II/multigrid/mg_transfer.templates.h | 57 +++--- .../multigrid/mg_transfer_global_coarsening.h | 105 ++++++----- .../mg_transfer_global_coarsening.templates.h | 168 ++++++++++-------- source/multigrid/mg_level_global_transfer.cc | 24 ++- .../mg_level_global_transfer.inst.in | 74 +++++--- .../mg_transfer_global_coarsening.inst.in | 10 +- source/multigrid/mg_transfer_prebuilt.inst.in | 4 +- 9 files changed, 284 insertions(+), 214 deletions(-) diff --git a/cmake/config/template-arguments.in b/cmake/config/template-arguments.in index 96c0175391..19e7b3a571 100644 --- a/cmake/config/template-arguments.in +++ b/cmake/config/template-arguments.in @@ -169,18 +169,13 @@ EXTERNAL_PARALLEL_VECTORS := { @DEAL_II_EXPAND_TRILINOS_MPI_VECTOR@; @DEAL_II_EXPAND_PETSC_MPI_BLOCKVECTOR@ } -// TODO: I don't understand this one. LA::distributed::Vector does not have a -// native matrix type, especially if we don't compile with Trilinos. Currently -// only used: mg_transfer_prebuilt.inst.in and -// mg_level_global_transfer.inst.in -VECTORS_WITH_MATRIX := { Vector; +// Vectors we can do GMG with excluding the LinearAlgebra::distributed::Vector<> (which is handled separately): +VECTORS_WITHOUT_LAVEC := { Vector; Vector ; BlockVector; BlockVector; - LinearAlgebra::distributed::Vector; - @DEAL_II_EXPAND_TRILINOS_MPI_VECTOR@; @DEAL_II_EXPAND_EPETRA_VECTOR@; @DEAL_II_EXPAND_TPETRA_VECTOR_DOUBLE@; diff --git a/include/deal.II/multigrid/mg_transfer.h b/include/deal.II/multigrid/mg_transfer.h index 909300fdbe..9b038bc510 100644 --- a/include/deal.II/multigrid/mg_transfer.h +++ b/include/deal.II/multigrid/mg_transfer.h @@ -407,9 +407,11 @@ private: * routines as compared to the %parallel vectors in the PETScWrappers and * TrilinosWrappers namespaces. */ -template -class MGLevelGlobalTransfer> - : public MGTransferBase> +template +class MGLevelGlobalTransfer< + LinearAlgebra::distributed::Vector> + : public MGTransferBase< + LinearAlgebra::distributed::Vector> { public: /** @@ -426,9 +428,10 @@ public: */ template void - copy_to_mg(const DoFHandler &dof_handler, - MGLevelObject> &dst, - const LinearAlgebra::distributed::Vector &src) const; + copy_to_mg( + const DoFHandler &dof_handler, + MGLevelObject> &dst, + const LinearAlgebra::distributed::Vector &src) const; /** * Transfer from multi-level vector to normal vector. @@ -440,9 +443,10 @@ public: template void copy_from_mg( - const DoFHandler &dof_handler, - LinearAlgebra::distributed::Vector &dst, - const MGLevelObject> &src) const; + const DoFHandler &dof_handler, + LinearAlgebra::distributed::Vector &dst, + const MGLevelObject> + &src) const; /** * Add a multi-level vector to a normal vector. @@ -452,9 +456,10 @@ public: template void copy_from_mg_add( - const DoFHandler &dof_handler, - LinearAlgebra::distributed::Vector &dst, - const MGLevelObject> &src) const; + const DoFHandler &dof_handler, + LinearAlgebra::distributed::Vector &dst, + const MGLevelObject> + &src) const; /** * If this object operates on BlockVector objects, we need to describe how @@ -493,10 +498,11 @@ protected: */ template void - copy_to_mg(const DoFHandler &dof_handler, - MGLevelObject> &dst, - const LinearAlgebra::distributed::Vector &src, - const bool solution_transfer) const; + copy_to_mg( + const DoFHandler &dof_handler, + MGLevelObject> &dst, + const LinearAlgebra::distributed::Vector &src, + const bool solution_transfer) const; /** * Internal function to @p fill copy_indices*. Called by derived classes. @@ -581,26 +587,27 @@ protected: * global vector for inserting into the level vectors. This vector is * populated with those entries. */ - mutable LinearAlgebra::distributed::Vector ghosted_global_vector; + mutable LinearAlgebra::distributed::Vector + ghosted_global_vector; /** * Same as above but used when working with solution vectors. */ - mutable LinearAlgebra::distributed::Vector + mutable LinearAlgebra::distributed::Vector solution_ghosted_global_vector; /** * In the function copy_from_mg, we access all level vectors with certain * ghost entries for inserting the result into a global vector. */ - mutable MGLevelObject> + mutable MGLevelObject> ghosted_level_vector; /** * Function to initialize internal level vectors. */ std::function &)> + LinearAlgebra::distributed::Vector &)> initialize_dof_vector; private: diff --git a/include/deal.II/multigrid/mg_transfer.templates.h b/include/deal.II/multigrid/mg_transfer.templates.h index 4b89ad6ee9..effad1c943 100644 --- a/include/deal.II/multigrid/mg_transfer.templates.h +++ b/include/deal.II/multigrid/mg_transfer.templates.h @@ -394,27 +394,29 @@ MGLevelGlobalTransfer::assert_built( /* --------- MGLevelGlobalTransfer ------- */ -template +template template void -MGLevelGlobalTransfer>::copy_to_mg( - const DoFHandler &dof_handler, - MGLevelObject> &dst, - const LinearAlgebra::distributed::Vector &src) const +MGLevelGlobalTransfer>:: + copy_to_mg( + const DoFHandler &dof_handler, + MGLevelObject> &dst, + const LinearAlgebra::distributed::Vector &src) const { assert_built(dof_handler); copy_to_mg(dof_handler, dst, src, false); } -template +template template void -MGLevelGlobalTransfer>::copy_to_mg( - const DoFHandler &dof_handler, - MGLevelObject> &dst, - const LinearAlgebra::distributed::Vector &src, - const bool solution_transfer) const +MGLevelGlobalTransfer>:: + copy_to_mg( + const DoFHandler &dof_handler, + MGLevelObject> &dst, + const LinearAlgebra::distributed::Vector &src, + const bool solution_transfer) const { assert_built(dof_handler); LinearAlgebra::distributed::Vector &this_ghosted_global_vector = @@ -513,13 +515,15 @@ MGLevelGlobalTransfer>::copy_to_mg( -template +template template void -MGLevelGlobalTransfer>::copy_from_mg( - const DoFHandler &dof_handler, - LinearAlgebra::distributed::Vector &dst, - const MGLevelObject> &src) const +MGLevelGlobalTransfer>:: + copy_from_mg( + const DoFHandler &dof_handler, + LinearAlgebra::distributed::Vector &dst, + const MGLevelObject> + &src) const { assert_built(dof_handler); AssertIndexRange(src.max_level(), @@ -598,14 +602,15 @@ MGLevelGlobalTransfer>::copy_from_mg( -template +template template void -MGLevelGlobalTransfer>:: +MGLevelGlobalTransfer>:: copy_from_mg_add( - const DoFHandler &dof_handler, - LinearAlgebra::distributed::Vector &dst, - const MGLevelObject> &src) const + const DoFHandler &dof_handler, + LinearAlgebra::distributed::Vector &dst, + const MGLevelObject> + &src) const { assert_built(dof_handler); // For non-DG: degrees of freedom in the refinement face may need special @@ -645,19 +650,19 @@ MGLevelGlobalTransfer>:: -template +template void -MGLevelGlobalTransfer>:: +MGLevelGlobalTransfer>:: set_component_to_block_map(const std::vector &map) { component_to_block_map = map; } -template +template template void -MGLevelGlobalTransfer>::assert_built( - const DoFHandler &dof_handler) const +MGLevelGlobalTransfer>:: + assert_built(const DoFHandler &dof_handler) const { (void)dof_handler; Assert(copy_indices.size() == diff --git a/include/deal.II/multigrid/mg_transfer_global_coarsening.h b/include/deal.II/multigrid/mg_transfer_global_coarsening.h index b982d60e66..76d469e23b 100644 --- a/include/deal.II/multigrid/mg_transfer_global_coarsening.h +++ b/include/deal.II/multigrid/mg_transfer_global_coarsening.h @@ -15,6 +15,7 @@ #ifndef dealii_mg_transfer_global_coarsening_h #define dealii_mg_transfer_global_coarsening_h +#include #include #include #include @@ -50,7 +51,7 @@ namespace RepartitioningPolicyTools class Base; } -template +template class MGTransferMF; #endif @@ -227,7 +228,12 @@ public: static_assert( std::is_same_v< VectorType, - LinearAlgebra::distributed::Vector>, + LinearAlgebra::distributed::Vector> || + std::is_same_v< + VectorType, + LinearAlgebra::distributed::Vector>, "This class is currently only implemented for vectors of " "type LinearAlgebra::distributed::Vector."); @@ -286,17 +292,13 @@ protected: * Perform prolongation on vectors with correct ghosting. */ virtual void - prolongate_and_add_internal( - LinearAlgebra::distributed::Vector &dst, - const LinearAlgebra::distributed::Vector &src) const = 0; + prolongate_and_add_internal(VectorType &dst, const VectorType &src) const = 0; /** * Perform restriction on vectors with correct ghosting. */ virtual void - restrict_and_add_internal( - LinearAlgebra::distributed::Vector &dst, - const LinearAlgebra::distributed::Vector &src) const = 0; + restrict_and_add_internal(VectorType &dst, const VectorType &src) const = 0; /** * A wrapper around update_ghost_values() optimized in case the @@ -360,7 +362,7 @@ protected: /** * Internal vector on which the actual prolongation/restriction is performed. */ - mutable LinearAlgebra::distributed::Vector vec_coarse; + mutable VectorType vec_coarse; /** * Internal vector needed for collecting all degrees of freedom of the fine @@ -369,7 +371,7 @@ protected: * accessed by the given vectors in the prolongate/restrict functions), * otherwise it is left at size zero. */ - mutable LinearAlgebra::distributed::Vector vec_fine; + mutable VectorType vec_fine; /** * Bool indicating whether fine vector has relevant ghost values. @@ -436,7 +438,12 @@ public: static_assert( std::is_same_v< VectorType, - LinearAlgebra::distributed::Vector>, + LinearAlgebra::distributed::Vector> || + std::is_same_v< + VectorType, + LinearAlgebra::distributed::Vector>, "This class is currently only implemented for vectors of " "type LinearAlgebra::distributed::Vector."); @@ -725,7 +732,7 @@ private: friend class internal::MGTwoLevelTransferImplementation; - friend class MGTransferMF; + friend class MGTransferMF; }; @@ -956,7 +963,7 @@ private: */ std::vector level_dof_indices_fine_ptrs; - friend class MGTransferMF; + friend class MGTransferMF; }; @@ -989,15 +996,15 @@ private: * * The implementation of this class is explained in detail in @cite munch2022gc. */ -template +template class MGTransferMF : public dealii::MGLevelGlobalTransfer< - LinearAlgebra::distributed::Vector> + LinearAlgebra::distributed::Vector> { public: /** * Value type. */ - using VectorType = LinearAlgebra::distributed::Vector; + using VectorType = LinearAlgebra::distributed::Vector; /** * Default constructor. @@ -1325,13 +1332,17 @@ private: */ template class MGTransferBlockMF - : public MGTransferBlockMatrixFreeBase> + : public MGTransferBlockMatrixFreeBase< + dim, + Number, + MGTransferMF> { public: /** * Constructor. */ - MGTransferBlockMF(const MGTransferMF &transfer_operator); + MGTransferBlockMF(const MGTransferMF + &transfer_operator); /** * Constructor. @@ -1388,27 +1399,30 @@ public: build(const std::vector *> &dof_handler); protected: - const MGTransferMF & + const MGTransferMF & get_matrix_free_transfer(const unsigned int b) const override; private: /** * Internal non-block version of transfer operation. */ - std::vector> transfer_operators_internal; + std::vector> + transfer_operators_internal; /** * Non-block version of transfer operation. */ - std::vector>> + std::vector>> transfer_operators; }; template -using MGTransferGlobalCoarsening = - MGTransferMF; +using MGTransferGlobalCoarsening = MGTransferMF; template using MGTransferBlockGlobalCoarsening = @@ -1423,9 +1437,9 @@ using MGTransferBlockGlobalCoarsening = -template +template template -MGTransferMF::MGTransferMF( +MGTransferMF::MGTransferMF( const MGLevelObject &transfer, const std::function &initialize_dof_vector) @@ -1439,10 +1453,10 @@ MGTransferMF::MGTransferMF( -template +template template void -MGTransferMF::initialize_two_level_transfers( +MGTransferMF::initialize_two_level_transfers( const MGLevelObject &transfer) { this->initialize_transfer_references(transfer); @@ -1450,10 +1464,10 @@ MGTransferMF::initialize_two_level_transfers( -template +template template void -MGTransferMF::initialize_transfer_references( +MGTransferMF::initialize_transfer_references( const MGLevelObject &transfer) { const unsigned int min_level = transfer.min_level(); @@ -1470,10 +1484,10 @@ MGTransferMF::initialize_transfer_references( -template +template template void -MGTransferMF::initialize_dof_vector( +MGTransferMF::initialize_dof_vector( const unsigned int level, VectorType &vec, const InVector &vec_reference, @@ -1518,12 +1532,13 @@ MGTransferMF::initialize_dof_vector( -template +template template void -MGTransferMF::copy_to_mg(const DoFHandler &dof_handler, - MGLevelObject &dst, - const InVector &src) const +MGTransferMF::copy_to_mg( + const DoFHandler &dof_handler, + MGLevelObject &dst, + const InVector &src) const { assert_dof_handler(dof_handler); @@ -1576,10 +1591,10 @@ MGTransferMF::copy_to_mg(const DoFHandler &dof_handler, -template +template template void -MGTransferMF::copy_from_mg( +MGTransferMF::copy_from_mg( const DoFHandler &dof_handler, OutVector &dst, const MGLevelObject &src) const @@ -1630,11 +1645,12 @@ MGTransferMF::copy_from_mg( -template +template template void -MGTransferMF::interpolate_to_mg(MGLevelObject &dst, - const InVector &src) const +MGTransferMF::interpolate_to_mg( + MGLevelObject &dst, + const InVector &src) const { DoFHandler dof_handler_dummy; @@ -1643,12 +1659,13 @@ MGTransferMF::interpolate_to_mg(MGLevelObject &dst, -template +template template void -MGTransferMF::interpolate_to_mg(const DoFHandler &dof_handler, - MGLevelObject &dst, - const InVector &src) const +MGTransferMF::interpolate_to_mg( + const DoFHandler &dof_handler, + MGLevelObject &dst, + const InVector &src) const { assert_dof_handler(dof_handler); diff --git a/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h b/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h index f34bc2bb4e..09ed15f80d 100644 --- a/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h +++ b/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h @@ -18,6 +18,7 @@ #include +#include "deal.II/base/memory_space.h" #include #include @@ -1471,13 +1472,14 @@ namespace internal /** * Compute weights. */ - template + template static void setup_weights( const dealii::AffineConstraints &constraints_fine, - MGTwoLevelTransfer> - &transfer, - const bool is_feq) + MGTwoLevelTransfer< + dim, + LinearAlgebra::distributed::Vector> &transfer, + const bool is_feq) { if (transfer.fine_element_is_continuous == false) return; // nothing to do @@ -1613,7 +1615,7 @@ namespace internal - template + template static void reinit_geometric_transfer( const DoFHandler &dof_handler_fine, @@ -1622,8 +1624,9 @@ namespace internal const dealii::AffineConstraints &constraints_coarse, const unsigned int mg_level_fine, const unsigned int mg_level_coarse, - MGTwoLevelTransfer> - &transfer) + MGTwoLevelTransfer< + dim, + LinearAlgebra::distributed::Vector> &transfer) { Assert((mg_level_fine == numbers::invalid_unsigned_int && mg_level_coarse == numbers::invalid_unsigned_int) || @@ -2729,10 +2732,9 @@ namespace internal { const auto &vector_partitioner = vec.get_partitioner(); - ArrayView ghost_array( - const_cast &>(vec).begin() + - vector_partitioner->locally_owned_size(), - vector_partitioner->n_ghost_indices()); + ArrayView ghost_array(const_cast(vec).begin() + + vector_partitioner->locally_owned_size(), + vector_partitioner->n_ghost_indices()); for (const auto &my_ghosts : embedded_partitioner->ghost_indices_within_larger_ghost_set()) @@ -4255,8 +4257,8 @@ MGTwoLevelTransferBase::zero_out_ghost_values( -template -MGTransferMF::MGTransferMF() +template +MGTransferMF::MGTransferMF() { this->transfer.clear(); this->internal_transfer.clear(); @@ -4264,8 +4266,8 @@ MGTransferMF::MGTransferMF() -template -MGTransferMF::MGTransferMF( +template +MGTransferMF::MGTransferMF( const MGConstrainedDoFs &mg_constrained_dofs) { this->transfer.clear(); @@ -4275,9 +4277,9 @@ MGTransferMF::MGTransferMF( -template +template void -MGTransferMF::initialize_constraints( +MGTransferMF::initialize_constraints( const MGConstrainedDoFs &mg_constrained_dofs) { this->mg_constrained_dofs = &mg_constrained_dofs; @@ -4285,9 +4287,9 @@ MGTransferMF::initialize_constraints( -template +template void -MGTransferMF::initialize_internal_transfer( +MGTransferMF::initialize_internal_transfer( const DoFHandler &dof_handler, const ObserverPointer &mg_constrained_dofs) { @@ -4317,41 +4319,48 @@ MGTransferMF::initialize_internal_transfer( -template +template std::pair *, unsigned int> -MGTransferMF::get_dof_handler_fine() const +MGTransferMF::get_dof_handler_fine() const { if (this->transfer.n_levels() <= 1) // single level: the information cannot be retrieved return {nullptr, numbers::invalid_unsigned_int}; - if (const auto t = dynamic_cast< - const MGTwoLevelTransfer> *>( + if (const auto t = dynamic_cast> *>( this->transfer[this->transfer.max_level()].get())) { return {t->dof_handler_fine, t->mg_level_fine}; } - else if (const auto t = dynamic_cast> *>( - this->transfer[this->transfer.max_level()].get())) - { - return {t->dof_handler_fine, t->mg_level_fine}; - } - else + + if constexpr (std::is_same_v) { - DEAL_II_NOT_IMPLEMENTED(); - return {nullptr, numbers::invalid_unsigned_int}; + // MGTwoLevelTransferNonNested transfer is only instantiated for Host + // memory: + if (const auto t = dynamic_cast> *>( + this->transfer[this->transfer.max_level()].get())) + { + return {t->dof_handler_fine, t->mg_level_fine}; + } } + + { + DEAL_II_NOT_IMPLEMENTED(); + return {nullptr, numbers::invalid_unsigned_int}; + } } -template +template void -MGTransferMF::fill_and_communicate_copy_indices_global_coarsening( - const DoFHandler &dof_handler_out) +MGTransferMF:: + fill_and_communicate_copy_indices_global_coarsening( + const DoFHandler &dof_handler_out) { const auto dof_handler_and_level_in = get_dof_handler_fine(); const auto dof_handler_in = dof_handler_and_level_in.first; @@ -4448,9 +4457,9 @@ MGTransferMF::fill_and_communicate_copy_indices_global_coarsening( -template +template void -MGTransferMF::build( +MGTransferMF::build( const std::vector> &external_partitioners) { @@ -4489,9 +4498,9 @@ MGTransferMF::build( -template +template void -MGTransferMF::build( +MGTransferMF::build( const std::function &initialize_dof_vector) { @@ -4506,7 +4515,8 @@ MGTransferMF::build( for (unsigned int l = min_level; l <= max_level; ++l) { - LinearAlgebra::distributed::Vector + LinearAlgebra::distributed::Vector vector; initialize_dof_vector(l, vector); external_partitioners[l - min_level] = vector.get_partitioner(); @@ -4522,9 +4532,9 @@ MGTransferMF::build( -template +template void -MGTransferMF::build( +MGTransferMF::build( const DoFHandler &dof_handler, const std::vector> &external_partitioners) @@ -4549,9 +4559,9 @@ MGTransferMF::build( -template +template void -MGTransferMF::build( +MGTransferMF::build( const DoFHandler &dof_handler, const std::function &initialize_dof_vector) @@ -4576,11 +4586,11 @@ MGTransferMF::build( -template +template void -MGTransferMF::prolongate(const unsigned int to_level, - VectorType &dst, - const VectorType &src) const +MGTransferMF::prolongate(const unsigned int to_level, + VectorType &dst, + const VectorType &src) const { dst = Number(0.0); prolongate_and_add(to_level, dst, src); @@ -4588,31 +4598,33 @@ MGTransferMF::prolongate(const unsigned int to_level, -template +template void -MGTransferMF::prolongate_and_add(const unsigned int to_level, - VectorType &dst, - const VectorType &src) const +MGTransferMF::prolongate_and_add( + const unsigned int to_level, + VectorType &dst, + const VectorType &src) const { this->transfer[to_level]->prolongate_and_add(dst, src); } -template +template void -MGTransferMF::restrict_and_add(const unsigned int from_level, - VectorType &dst, - const VectorType &src) const +MGTransferMF::restrict_and_add( + const unsigned int from_level, + VectorType &dst, + const VectorType &src) const { this->transfer[from_level]->restrict_and_add(dst, src); } -template +template void -MGTransferMF::assert_dof_handler( +MGTransferMF::assert_dof_handler( const DoFHandler &dof_handler_out) const { #ifndef DEBUG @@ -4699,9 +4711,9 @@ MGTransferMF::assert_dof_handler( -template +template std::size_t -MGTransferMF::memory_consumption() const +MGTransferMF::memory_consumption() const { std::size_t size = 0; @@ -4716,26 +4728,26 @@ MGTransferMF::memory_consumption() const -template +template inline unsigned int -MGTransferMF::min_level() const +MGTransferMF::min_level() const { return transfer.min_level(); } -template +template inline unsigned int -MGTransferMF::max_level() const +MGTransferMF::max_level() const { return transfer.max_level(); } -template +template inline void -MGTransferMF::clear() +MGTransferMF::clear() { MGLevelGlobalTransfer::clear(); @@ -4748,8 +4760,12 @@ MGTransferMF::clear() template MGTransferBlockMF::MGTransferBlockMF( - const MGTransferMF &transfer_operator) - : MGTransferBlockMatrixFreeBase>(true) + const MGTransferMF + &transfer_operator) + : MGTransferBlockMatrixFreeBase< + dim, + Number, + MGTransferMF>(true) { this->transfer_operators = {&transfer_operator}; } @@ -4759,7 +4775,10 @@ MGTransferBlockMF::MGTransferBlockMF( template MGTransferBlockMF::MGTransferBlockMF( const MGConstrainedDoFs &mg_constrained_dofs) - : MGTransferBlockMatrixFreeBase>(true) + : MGTransferBlockMatrixFreeBase< + dim, + Number, + MGTransferMF>(true) { initialize_constraints(mg_constrained_dofs); } @@ -4769,7 +4788,10 @@ MGTransferBlockMF::MGTransferBlockMF( template MGTransferBlockMF::MGTransferBlockMF( const std::vector &mg_constrained_dofs) - : MGTransferBlockMatrixFreeBase>(false) + : MGTransferBlockMatrixFreeBase< + dim, + Number, + MGTransferMF>(false) { initialize_constraints(mg_constrained_dofs); } @@ -4843,7 +4865,7 @@ MGTransferBlockMF::build( template -const MGTransferMF & +const MGTransferMF & MGTransferBlockMF::get_matrix_free_transfer( const unsigned int b) const { diff --git a/source/multigrid/mg_level_global_transfer.cc b/source/multigrid/mg_level_global_transfer.cc index 34a7b49a1e..85943e942f 100644 --- a/source/multigrid/mg_level_global_transfer.cc +++ b/source/multigrid/mg_level_global_transfer.cc @@ -275,10 +275,10 @@ namespace } } // namespace -template +template template void -MGLevelGlobalTransfer>:: +MGLevelGlobalTransfer>:: fill_and_communicate_copy_indices(const DoFHandler &mg_dof) { const MPI_Comm mpi_communicator = mg_dof.get_mpi_communicator(); @@ -382,9 +382,10 @@ MGLevelGlobalTransfer>:: -template +template void -MGLevelGlobalTransfer>::clear() +MGLevelGlobalTransfer< + LinearAlgebra::distributed::Vector>::clear() { sizes.resize(0); copy_indices.clear(); @@ -404,9 +405,9 @@ MGLevelGlobalTransfer>::clear() -template +template void -MGLevelGlobalTransfer>:: +MGLevelGlobalTransfer>:: print_indices(std::ostream &os) const { for (unsigned int level = 0; level < copy_indices.size(); ++level) @@ -435,10 +436,10 @@ MGLevelGlobalTransfer>:: -template +template std::size_t -MGLevelGlobalTransfer< - LinearAlgebra::distributed::Vector>::memory_consumption() const +MGLevelGlobalTransfer>:: + memory_consumption() const { std::size_t result = sizeof(*this); result += MemoryConsumption::memory_consumption(sizes); @@ -459,9 +460,4 @@ MGLevelGlobalTransfer< // explicit instantiation #include "mg_level_global_transfer.inst" -// create an additional instantiation currently not supported by the automatic -// template instantiation scheme -template class MGLevelGlobalTransfer>; - - DEAL_II_NAMESPACE_CLOSE diff --git a/source/multigrid/mg_level_global_transfer.inst.in b/source/multigrid/mg_level_global_transfer.inst.in index 8f6ee1a1e7..3cbb47b81f 100644 --- a/source/multigrid/mg_level_global_transfer.inst.in +++ b/source/multigrid/mg_level_global_transfer.inst.in @@ -14,12 +14,18 @@ -for (V1 : VECTORS_WITH_MATRIX) +for (V1 : VECTORS_WITHOUT_LAVEC) { template class MGLevelGlobalTransfer; } -for (deal_II_dimension : DIMENSIONS; V1 : VECTORS_WITH_MATRIX) +for (S1 : REAL_SCALARS) + { + template class MGLevelGlobalTransfer< + LinearAlgebra::distributed::Vector>; + } + +for (deal_II_dimension : DIMENSIONS; V1 : VECTORS_WITHOUT_LAVEC) { template void MGLevelGlobalTransfer:: fill_and_communicate_copy_indices( @@ -43,42 +49,60 @@ for (deal_II_dimension : DIMENSIONS; V1, V2 : DEAL_II_VEC_TEMPLATES; const MGLevelObject> &) const; } -for (deal_II_dimension : DIMENSIONS) +for (deal_II_dimension : DIMENSIONS; S1 : REAL_SCALARS) { - template void - MGLevelGlobalTransfer>:: + template void MGLevelGlobalTransfer< + LinearAlgebra::distributed::Vector>:: fill_and_communicate_copy_indices( const DoFHandler &mg_dof); } for (deal_II_dimension : DIMENSIONS; S1, S2 : REAL_SCALARS) { - template void - MGLevelGlobalTransfer>::copy_to_mg( - const DoFHandler &, - MGLevelObject> &, - const LinearAlgebra::distributed::Vector &, - const bool) const; + template void MGLevelGlobalTransfer< + LinearAlgebra::distributed::Vector>:: + copy_to_mg( + const DoFHandler &, + MGLevelObject> + &, + const LinearAlgebra::distributed::Vector &, + const bool) const; + + // template void + // MGLevelGlobalTransfer>::copy_to_mg( const + // DoFHandler &, + // MGLevelObject> + // &, + // const LinearAlgebra::distributed::Vector &, + // const bool) const; } for (deal_II_dimension : DIMENSIONS; S1, S2 : REAL_SCALARS) { - template void - MGLevelGlobalTransfer>::copy_to_mg( - const DoFHandler &, - MGLevelObject> &, - const LinearAlgebra::distributed::Vector &) const; - template void - MGLevelGlobalTransfer>::copy_from_mg( - const DoFHandler &, - LinearAlgebra::distributed::Vector &, - const MGLevelObject> &) const; - template void - MGLevelGlobalTransfer>:: + template void MGLevelGlobalTransfer< + LinearAlgebra::distributed::Vector>:: + copy_to_mg( + const DoFHandler &, + MGLevelObject> + &, + const LinearAlgebra::distributed::Vector &) + const; + template void MGLevelGlobalTransfer< + LinearAlgebra::distributed::Vector>:: + copy_from_mg( + const DoFHandler &, + LinearAlgebra::distributed::Vector &, + const MGLevelObject< + LinearAlgebra::distributed::Vector> &) const; + template void MGLevelGlobalTransfer< + LinearAlgebra::distributed::Vector>:: copy_from_mg_add( const DoFHandler &, - LinearAlgebra::distributed::Vector &, - const MGLevelObject> &) const; + LinearAlgebra::distributed::Vector &, + const MGLevelObject< + LinearAlgebra::distributed::Vector> &) const; } for (deal_II_dimension : DIMENSIONS) diff --git a/source/multigrid/mg_transfer_global_coarsening.inst.in b/source/multigrid/mg_transfer_global_coarsening.inst.in index 058ac9fc0d..f4e40c0fdf 100644 --- a/source/multigrid/mg_transfer_global_coarsening.inst.in +++ b/source/multigrid/mg_transfer_global_coarsening.inst.in @@ -16,7 +16,7 @@ for (S1 : REAL_SCALARS) { template class MGTwoLevelTransferBase< - LinearAlgebra::distributed::Vector>; + LinearAlgebra::distributed::Vector>; } for (deal_II_dimension : DIMENSIONS; S1 : REAL_SCALARS) @@ -26,12 +26,16 @@ for (deal_II_dimension : DIMENSIONS; S1 : REAL_SCALARS) template class MGTwoLevelTransferNonNested< deal_II_dimension, LinearAlgebra::distributed::Vector>; + template class MGTransferBlockMatrixFreeBase< deal_II_dimension, S1, - MGTransferMF>; - template class MGTransferMF; + MGTransferMF>; template class MGTransferBlockMF; + + template class MGTransferMF; } for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS) diff --git a/source/multigrid/mg_transfer_prebuilt.inst.in b/source/multigrid/mg_transfer_prebuilt.inst.in index 9881e1bc9c..330426dff3 100644 --- a/source/multigrid/mg_transfer_prebuilt.inst.in +++ b/source/multigrid/mg_transfer_prebuilt.inst.in @@ -14,12 +14,12 @@ -for (V1 : VECTORS_WITH_MATRIX) +for (V1 : VECTORS_WITHOUT_LAVEC) { template class MGTransferPrebuilt; } -for (deal_II_dimension : DIMENSIONS; V1 : VECTORS_WITH_MATRIX) +for (deal_II_dimension : DIMENSIONS; V1 : VECTORS_WITHOUT_LAVEC) { template void MGTransferPrebuilt::build( const DoFHandler &mg_dof); -- 2.39.5