@@ -2746,15 +2746,15 @@ namespace internal
27462746
27472747
27482748template <typename VectorType>
2749- MGTwoLevelTransferBase <VectorType>::MGTwoLevelTransferBase ()
2749+ MGTwoLevelTransferCore <VectorType>::MGTwoLevelTransferCore ()
27502750 : vec_fine_needs_ghost_update(true )
27512751{}
27522752
27532753
27542754
27552755template <typename VectorType>
27562756void
2757- MGTwoLevelTransferBase <VectorType>::prolongate_and_add(
2757+ MGTwoLevelTransferCore <VectorType>::prolongate_and_add(
27582758 VectorType &dst,
27592759 const VectorType &src) const
27602760{
@@ -2793,7 +2793,6 @@ MGTwoLevelTransferBase<VectorType>::prolongate_and_add(
27932793}
27942794
27952795
2796-
27972796namespace internal
27982797{
27992798 // Helper class to compute correct weights, which works by simply using
@@ -3038,7 +3037,7 @@ MGTwoLevelTransfer<dim, VectorType>::prolongate_and_add_internal(
30383037
30393038template <typename VectorType>
30403039void
3041- MGTwoLevelTransferBase <VectorType>::restrict_and_add(
3040+ MGTwoLevelTransferCore <VectorType>::restrict_and_add(
30423041 VectorType &dst,
30433042 const VectorType &src) const
30443043{
@@ -3524,7 +3523,7 @@ namespace internal
35243523template <typename VectorType>
35253524template <int dim, std::size_t width, typename IndexType>
35263525std::pair<bool , bool >
3527- MGTwoLevelTransferBase <VectorType>::
3526+ MGTwoLevelTransferCore <VectorType>::
35283527 internal_enable_inplace_operations_if_possible (
35293528 const std::shared_ptr<const Utilities::MPI::Partitioner>
35303529 &external_partitioner_coarse,
@@ -3949,7 +3948,7 @@ MGTwoLevelTransfer<dim, VectorType>::memory_consumption() const
39493948
39503949template <typename VectorType>
39513950void
3952- MGTwoLevelTransferBase <VectorType>::update_ghost_values(
3951+ MGTwoLevelTransferCore <VectorType>::update_ghost_values(
39533952 const VectorType &vec) const
39543953{
39553954 if ((vec.get_partitioner ().get () == this ->partitioner_coarse .get ()) &&
@@ -3970,7 +3969,7 @@ MGTwoLevelTransferBase<VectorType>::update_ghost_values(
39703969
39713970template <typename VectorType>
39723971void
3973- MGTwoLevelTransferBase <VectorType>::compress(
3972+ MGTwoLevelTransferCore <VectorType>::compress(
39743973 VectorType &vec,
39753974 const VectorOperation::values op) const
39763975{
@@ -3994,7 +3993,7 @@ MGTwoLevelTransferBase<VectorType>::compress(
39943993
39953994template <typename VectorType>
39963995void
3997- MGTwoLevelTransferBase <VectorType>::zero_out_ghost_values(
3996+ MGTwoLevelTransferCore <VectorType>::zero_out_ghost_values(
39983997 const VectorType &vec) const
39993998{
40003999 if ((vec.get_partitioner ().get () == this ->partitioner_coarse .get ()) &&
@@ -4091,6 +4090,7 @@ MGTransferMatrixFree<dim, Number, MemorySpace>::get_dof_handler_fine() const
40914090 return {t->dof_handler_fine , t->mg_level_fine };
40924091 }
40934092
4093+
40944094 if constexpr (std::is_same_v<MemorySpace, ::dealii::MemorySpace::Host>)
40954095 {
40964096 // MGTwoLevelTransferNonNested transfer is only instantiated for Host
@@ -4104,6 +4104,19 @@ MGTransferMatrixFree<dim, Number, MemorySpace>::get_dof_handler_fine() const
41044104 }
41054105 }
41064106
4107+ if constexpr (std::is_same_v<MemorySpace, ::dealii::MemorySpace::Default>)
4108+ {
4109+ // MGTwoLevelTransferCopyToHost should only be instantiated on device
4110+ // memory:
4111+ if (const auto t = dynamic_cast <const MGTwoLevelTransferCopyToHost<
4112+ dim,
4113+ LinearAlgebra::distributed::Vector<Number, MemorySpace>> *>(
4114+ this ->transfer [this ->transfer .max_level ()].get ()))
4115+ {
4116+ return {(t->host_transfer ).dof_handler_fine ,
4117+ (t->host_transfer ).mg_level_fine };
4118+ }
4119+ }
41074120 {
41084121 DEAL_II_NOT_IMPLEMENTED ();
41094122 return {nullptr , numbers::invalid_unsigned_int};
@@ -4297,22 +4310,33 @@ MGTransferMatrixFree<dim, Number, MemorySpace>::build(
42974310 const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
42984311 &external_partitioners)
42994312{
4300- const bool use_local_smoothing =
4301- this ->transfer .n_levels () == 0 || this ->internal_transfer .n_levels () > 0 ;
4302-
4303- if (use_local_smoothing)
4313+ // Local smoothing is not supported on the device
4314+ if constexpr (std::is_same_v<MemorySpace, dealii::MemorySpace::Host>)
43044315 {
4305- this ->initialize_internal_transfer (dof_handler,
4306- this ->mg_constrained_dofs );
4307- this ->initialize_transfer_references (internal_transfer);
4316+ bool use_local_smoothing = this ->transfer .n_levels () == 0 ||
4317+ this ->internal_transfer .n_levels () > 0 ;
4318+
4319+ if (use_local_smoothing)
4320+ {
4321+ this ->initialize_internal_transfer (dof_handler,
4322+ this ->mg_constrained_dofs );
4323+ this ->initialize_transfer_references (internal_transfer);
4324+ }
4325+
4326+ this ->build (external_partitioners);
4327+
4328+ {
4329+ if (use_local_smoothing)
4330+ this ->fill_and_communicate_copy_indices (dof_handler);
4331+ else
4332+ this ->fill_and_communicate_copy_indices_global_coarsening (
4333+ dof_handler);
4334+ }
43084335 }
43094336
43104337 this ->build (external_partitioners);
43114338
4312- if (use_local_smoothing)
4313- this ->fill_and_communicate_copy_indices (dof_handler);
4314- else
4315- this ->fill_and_communicate_copy_indices_global_coarsening (dof_handler);
4339+ this ->fill_and_communicate_copy_indices_global_coarsening (dof_handler);
43164340}
43174341
43184342
@@ -4324,22 +4348,33 @@ MGTransferMatrixFree<dim, Number, MemorySpace>::build(
43244348 const std::function<void (const unsigned int , VectorType &)>
43254349 &initialize_dof_vector)
43264350{
4327- const bool use_local_smoothing =
4328- this ->transfer .n_levels () == 0 || this ->internal_transfer .n_levels () > 0 ;
4329-
4330- if (use_local_smoothing)
4351+ // Local smoothing is not supported on the device
4352+ if constexpr (std::is_same_v<MemorySpace, dealii::MemorySpace::Host>)
43314353 {
4332- this ->initialize_internal_transfer (dof_handler,
4333- this ->mg_constrained_dofs );
4334- this ->initialize_transfer_references (internal_transfer);
4354+ bool use_local_smoothing = this ->transfer .n_levels () == 0 ||
4355+ this ->internal_transfer .n_levels () > 0 ;
4356+
4357+ if (use_local_smoothing)
4358+ {
4359+ this ->initialize_internal_transfer (dof_handler,
4360+ this ->mg_constrained_dofs );
4361+ this ->initialize_transfer_references (internal_transfer);
4362+ }
4363+
4364+ this ->build (initialize_dof_vector);
4365+
4366+ {
4367+ if (use_local_smoothing)
4368+ this ->fill_and_communicate_copy_indices (dof_handler);
4369+ else
4370+ this ->fill_and_communicate_copy_indices_global_coarsening (
4371+ dof_handler);
4372+ }
43354373 }
43364374
43374375 this ->build (initialize_dof_vector);
43384376
4339- if (use_local_smoothing)
4340- this ->fill_and_communicate_copy_indices (dof_handler);
4341- else
4342- this ->fill_and_communicate_copy_indices_global_coarsening (dof_handler);
4377+ this ->fill_and_communicate_copy_indices_global_coarsening (dof_handler);
43434378}
43444379
43454380
0 commit comments