22 template <
int dim,
typename VectorType>
31 for (
unsigned int l = transfer_matrices_.
min_level();
43 template <
int dim,
typename VectorType>
52 for (
unsigned int l = transfer_matrices_.
min_level();
58 const_cast<dealii::TrilinosWrappers::SparseMatrix *
>(
59 &transfer_matrices_[
l]);
66 template <
int dim,
typename VectorType>
69 const unsigned int to_level,
71 const VectorType &src)
const
73 dst =
typename VectorType::value_type(0.0);
79 template <
int dim,
typename VectorType>
82 const unsigned int to_level,
84 const VectorType &src)
const
87 ExcMessage(
"Transfer matrix has not been initialized."));
93 template <
int dim,
typename VectorType>
96 const unsigned int from_level,
98 const VectorType &src)
const
106 ExcMessage(
"Matrix has not been initialized."));
112 template <
int dim,
typename VectorType>
117 const VectorType &src)
const
129 if constexpr (std::is_same_v<VectorType, TrilinosWrappers::MPI::Vector>)
131 else if constexpr (std::is_same_v<
134 dst[dst.
max_level()].copy_locally_owned_data_from(src);
141 template <
int dim,
typename VectorType>
149 if constexpr (std::is_same_v<VectorType, TrilinosWrappers::MPI::Vector>)
151 else if constexpr (std::is_same_v<
154 dst.copy_locally_owned_data_from(src[src.
max_level()]);
unsigned int max_level() const
unsigned int min_level() const
unsigned int n_levels() const
MGLevelObject< SmartPointer< TrilinosWrappers::SparseMatrix > > transfer_matrices
void prolongate(const unsigned int to_level, VectorType &dst, const VectorType &src) const override
void copy_from_mg(const DoFHandler< dim > &dof_handler, VectorType &dst, const MGLevelObject< VectorType > &src) const
MGTransferAgglomeration(const MGLevelObject< TrilinosWrappers::SparseMatrix * > &transfer_matrices, const std::vector< DoFHandler< dim > * > &dof_handlers)
void restrict_and_add(const unsigned int from_level, VectorType &dst, const VectorType &src) const override
std::vector< const DoFHandler< dim > * > dof_handlers
void prolongate_and_add(const unsigned int to_level, VectorType &dst, const VectorType &src) const override
void copy_to_mg(const DoFHandler< dim > &dof_handler, MGLevelObject< VectorType > &dst, const VectorType &src) const
#define DEAL_II_NOT_IMPLEMENTED()
#define Assert(cond, exc)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)