16#ifndef dealii_mg_transfer_global_coarsening_h
17#define dealii_mg_transfer_global_coarsening_h
39 class MGTwoLevelTransferImplementation;
44 template <
int dim,
int spacedim>
88 const unsigned int degree,
96 std::vector<unsigned int>
98 const unsigned int max_degree,
111 template <
int dim,
int spacedim>
112 std::vector<std::shared_ptr<const Triangulation<dim, spacedim>>>
132 template <
int dim,
int spacedim>
133 std::vector<std::shared_ptr<const Triangulation<dim, spacedim>>>
137 const bool preserve_fine_triangulation,
138 const bool repartition_fine_triangulation);
145 template <
int dim,
int spacedim>
146 std::vector<std::shared_ptr<const Triangulation<dim, spacedim>>>
150 const bool repartition_fine_triangulation =
false);
159template <
int dim,
typename VectorType>
190 const std::shared_ptr<const Utilities::MPI::Partitioner>
192 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner_fine);
207template <
int dim,
typename Number>
280 const unsigned int fe_degree_coarse);
313 const std::shared_ptr<const Utilities::MPI::Partitioner>
315 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner_fine);
331 struct MGTransferScheme
459 friend class internal::MGTwoLevelTransferImplementation;
474template <
int dim,
typename VectorType>
481 using Number =
typename VectorType::value_type;
492 const std::function<
void(
const unsigned int, VectorType &)>
504 build(
const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
512 build(
const std::function<
void(
const unsigned int, VectorType &)>
521 const VectorType & src)
const override;
529 const VectorType & src)
const override;
537 const VectorType & src)
const override;
545 template <
class InVector,
int spacedim>
549 const InVector & src)
const;
557 template <
class OutVector,
int spacedim>
575 template <
class InVector>
585 template <
class InVector,
int spacedim>
589 const InVector & src)
const;
615 std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
626template <
int dim,
typename VectorType>
630 typename VectorType::value_type,
631 MGTransferGlobalCoarsening<dim, VectorType>>
659template <
int dim,
typename VectorType>
662 const std::function<
void(
const unsigned int, VectorType &)>
663 &initialize_dof_vector)
666 this->build(initialize_dof_vector);
671template <
int dim,
typename VectorType>
674 const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
675 &external_partitioners)
677 this->external_partitioners = external_partitioners;
679 if (this->external_partitioners.size() > 0)
681 const unsigned int min_level = transfer.min_level();
682 const unsigned int max_level = transfer.max_level();
684 AssertDimension(this->external_partitioners.size(), transfer.n_levels());
686 for (
unsigned int l = min_level + 1;
l <= max_level; ++
l)
687 transfer[l].enable_inplace_operations_if_possible(
688 this->external_partitioners[l - 1 - min_level],
689 this->external_partitioners[l - min_level]);
695template <
int dim,
typename VectorType>
698 const std::function<
void(
const unsigned int, VectorType &)>
699 &initialize_dof_vector)
701 if (initialize_dof_vector)
703 const unsigned int min_level = transfer.min_level();
704 const unsigned int max_level = transfer.max_level();
705 const unsigned int n_levels = transfer.n_levels();
707 std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
708 external_partitioners(n_levels);
710 for (
unsigned int l = min_level;
l <= max_level; ++
l)
714 initialize_dof_vector(l, vector);
718 this->build(external_partitioners);
728template <
int dim,
typename VectorType>
731 const unsigned int level,
732 VectorType & vec)
const
737 const auto &external_partitioner =
738 external_partitioners[
level - transfer.min_level()];
740 if (vec.get_partitioner().get() != external_partitioner.get())
741 vec.reinit(external_partitioner);
746template <
int dim,
typename VectorType>
749 const unsigned int to_level,
751 const VectorType & src)
const
754 prolongate_and_add(to_level, dst, src);
759template <
int dim,
typename VectorType>
762 const unsigned int to_level,
764 const VectorType & src)
const
766 this->transfer[to_level].prolongate_and_add(dst, src);
771template <
int dim,
typename VectorType>
774 const unsigned int from_level,
776 const VectorType & src)
const
778 this->transfer[from_level].restrict_and_add(dst, src);
783template <
int dim,
typename VectorType>
784template <
class InVector,
int spacedim>
789 const InVector & src)
const
798 dst[
level].copy_locally_owned_data_from(src);
806template <
int dim,
typename VectorType>
807template <
class OutVector,
int spacedim>
816 dst.copy_locally_owned_data_from(src[src.
max_level()]);
821template <
int dim,
typename VectorType>
822template <
class InVector>
826 const InVector & src)
const
828 const unsigned int min_level = transfer.min_level();
829 const unsigned int max_level = transfer.max_level();
837 dst[transfer.max_level()].copy_locally_owned_data_from(src);
839 for (
unsigned int l = max_level;
l > min_level; --
l)
845template <
int dim,
typename VectorType>
846template <
class InVector,
int spacedim>
851 const InVector & src)
const
855 this->interpolate_to_mg(dst, src);
860template <
int dim,
typename VectorType>
864 std::size_t size = 0;
866 const unsigned int min_level = transfer.min_level();
867 const unsigned int max_level = transfer.max_level();
869 for (
unsigned int l = min_level + 1;
l <= max_level; ++
l)
unsigned int max_level() const
unsigned int min_level() const
const MGTransferGlobalCoarsening< dim, VectorType > & get_matrix_free_transfer(const unsigned int b) const override
MGTransferBlockGlobalCoarsening(const MGTransferGlobalCoarsening< dim, VectorType > &transfer_operator)
const MGTransferGlobalCoarsening< dim, VectorType > & transfer_operator
MGTransferGlobalCoarsening(const MGLevelObject< MGTwoLevelTransfer< dim, VectorType > > &transfer, const std::function< void(const unsigned int, VectorType &)> &initialize_dof_vector={})
void copy_from_mg(const DoFHandler< dim, spacedim > &dof_handler, OutVector &dst, const MGLevelObject< VectorType > &src) const
void interpolate_to_mg(const DoFHandler< dim, spacedim > &dof_handler, MGLevelObject< VectorType > &dst, const InVector &src) const
void prolongate(const unsigned int to_level, VectorType &dst, const VectorType &src) const override
typename VectorType::value_type Number
std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > external_partitioners
void copy_to_mg(const DoFHandler< dim, spacedim > &dof_handler, MGLevelObject< VectorType > &dst, const InVector &src) const
virtual void restrict_and_add(const unsigned int from_level, VectorType &dst, const VectorType &src) const override
std::size_t memory_consumption() const
void initialize_dof_vector(const unsigned int level, VectorType &vector) const
void build(const std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > &external_partitioners={})
void build(const std::function< void(const unsigned int, VectorType &)> &initialize_dof_vector)
void interpolate_to_mg(MGLevelObject< VectorType > &dst, const InVector &src) const
MGLevelObject< MGTwoLevelTransfer< dim, VectorType > > transfer
void prolongate_and_add(const unsigned int to_level, VectorType &dst, const VectorType &src) const override
void enable_inplace_operations_if_possible(const std::shared_ptr< const Utilities::MPI::Partitioner > &partitioner_coarse, const std::shared_ptr< const Utilities::MPI::Partitioner > &partitioner_fine)
void prolongate_and_add(LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
LinearAlgebra::distributed::Vector< Number > vec_fine
std::vector< Number > weights
LinearAlgebra::distributed::Vector< Number > vec_coarse
std::size_t memory_consumption() const
void reinit_polynomial_transfer(const DoFHandler< dim > &dof_handler_fine, const DoFHandler< dim > &dof_handler_coarse, const AffineConstraints< Number > &constraint_fine=AffineConstraints< Number >(), const AffineConstraints< Number > &constraint_coarse=AffineConstraints< Number >(), const unsigned int mg_level_fine=numbers::invalid_unsigned_int, const unsigned int mg_level_coarse=numbers::invalid_unsigned_int)
static bool fast_polynomial_transfer_supported(const unsigned int fe_degree_fine, const unsigned int fe_degree_coarse)
std::shared_ptr< const Utilities::MPI::Partitioner > partitioner_coarse
unsigned int n_components
void interpolate(LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
std::shared_ptr< const Utilities::MPI::Partitioner > partitioner_fine
internal::MatrixFreeFunctions::ConstraintInfo< dim, VectorizedArray< Number > > constraint_info
void restrict_and_add(LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
bool fine_element_is_continuous
void reinit(const DoFHandler< dim > &dof_handler_fine, const DoFHandler< dim > &dof_handler_coarse, const AffineConstraints< Number > &constraint_fine=AffineConstraints< Number >(), const AffineConstraints< Number > &constraint_coarse=AffineConstraints< Number >(), const unsigned int mg_level_fine=numbers::invalid_unsigned_int, const unsigned int mg_level_coarse=numbers::invalid_unsigned_int)
void reinit_geometric_transfer(const DoFHandler< dim > &dof_handler_fine, const DoFHandler< dim > &dof_handler_coarse, const AffineConstraints< Number > &constraint_fine=AffineConstraints< Number >(), const AffineConstraints< Number > &constraint_coarse=AffineConstraints< Number >(), const unsigned int mg_level_fine=numbers::invalid_unsigned_int, const unsigned int mg_level_coarse=numbers::invalid_unsigned_int)
AlignedVector< VectorizedArray< Number > > weights_compressed
std::vector< unsigned int > level_dof_indices_fine
std::vector< MGTransferScheme > schemes
void restrict_and_add(VectorType &dst, const VectorType &src) const
void enable_inplace_operations_if_possible(const std::shared_ptr< const Utilities::MPI::Partitioner > &partitioner_coarse, const std::shared_ptr< const Utilities::MPI::Partitioner > &partitioner_fine)
std::size_t memory_consumption() const
void interpolate(VectorType &dst, const VectorType &src) const
void prolongate_and_add(VectorType &dst, const VectorType &src) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
const std::shared_ptr< const Utilities::MPI::Partitioner > & get_partitioner() const
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
static const unsigned int invalid_unsigned_int
unsigned int n_coarse_cells
internal::MatrixFreeFunctions::ShapeInfo< VectorizedArray< Number > > shape_info_coarse
AlignedVector< VectorizedArray< Number > > prolongation_matrix_1d
AlignedVector< VectorizedArray< Number > > prolongation_matrix
unsigned int n_dofs_per_cell_coarse
unsigned int n_dofs_per_cell_fine
unsigned int degree_coarse
AlignedVector< VectorizedArray< Number > > restriction_matrix_1d
AlignedVector< VectorizedArray< Number > > restriction_matrix
const ::Triangulation< dim, spacedim > & tria