16#ifndef dealii_mg_transfer_matrix_free_h
17#define dealii_mg_transfer_matrix_free_h
55template <
int dim,
typename Number>
106 const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
107 &external_partitioners =
108 std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>());
116 const std::function<
void(
const unsigned int,
136 const unsigned int to_level,
142 const unsigned int to_level,
166 const unsigned int from_level,
184 template <
typename BlockVectorType2,
int spacedim>
189 const BlockVectorType2 & src)
const;
245 std::vector<std::vector<std::pair<unsigned int, unsigned int>>>
297 template <
int degree>
300 const unsigned int to_level,
307 template <
int degree>
321template <
int dim,
typename Number,
typename TransferType>
323 :
public MGTransferBase<LinearAlgebra::distributed::BlockVector<Number>>
346 const unsigned int to_level,
352 const unsigned int to_level,
376 const unsigned int from_level,
390 template <
typename BlockVectorType2,
int spacedim>
395 const BlockVectorType2 & src)
const;
400 template <
typename BlockVectorType2,
int spacedim>
405 const BlockVectorType2 & src)
const;
410 template <
typename BlockVectorType2,
int spacedim>
414 BlockVectorType2 & dst,
421 template <
typename BlockVectorType2,
int spacedim>
425 BlockVectorType2 & dst,
440 virtual const TransferType &
465template <
int dim,
typename Number>
469 MGTransferMatrixFree<dim, Number>>
488 const std::vector<MGConstrainedDoFs> &mg_constrained_dofs);
506 const std::vector<MGConstrainedDoFs> &mg_constrained_dofs);
551template <
int dim,
typename Number>
552template <
typename BlockVectorType2,
int spacedim>
557 const BlockVectorType2 & src)
const
559 const unsigned int min_level = dst.min_level();
560 const unsigned int max_level = dst.max_level();
570 dst[
level].locally_owned_size() !=
572 dst[
level].reinit(this->vector_partitioners[
level]);
575 this->copy_to_mg(dof_handler, dst, src,
true);
582 dst[max_level].update_ghost_values();
589 if (dst[
level].get_partitioner().get() ==
590 this->vector_partitioners[
level].get())
594 ghosted_fine.
reinit(this->vector_partitioners[
level]);
597 input = &ghosted_fine;
603 std::vector<types::global_dof_index> dof_indices(fe.
n_dofs_per_cell());
605 if (cell->is_locally_owned_on_level())
610 if (cell->is_active())
613 std::fill(dof_values_coarse.begin(), dof_values_coarse.end(), 0.);
614 for (
unsigned int child = 0; child < cell->n_children(); ++child)
616 cell->child(child)->get_mg_dof_indices(dof_indices);
619 this->mg_constrained_dofs,
level, dof_indices);
622 dof_values_fine(i) = (*input)(dof_indices[i]);
624 .vmult(tmp, dof_values_fine);
627 dof_values_coarse[i] += tmp[i];
628 else if (tmp(i) != 0.)
629 dof_values_coarse[i] = tmp[i];
631 cell->get_mg_dof_indices(dof_indices);
635 dst[
level - 1](dof_indices[i]) = dof_values_coarse[i];
644template <
int dim,
typename Number,
typename TransferType>
645template <
typename BlockVectorType2,
int spacedim>
650 const BlockVectorType2 & src)
const
654 "This object was initialized with support for usage with one "
655 "DoFHandler for each block, but this method assumes that "
656 "the same DoFHandler is used for all the blocks!"));
657 const std::vector<const DoFHandler<dim, spacedim> *> mg_dofs(src.n_blocks(),
660 copy_to_mg(mg_dofs, dst, src);
665template <
int dim,
typename Number,
typename TransferType>
666template <
typename BlockVectorType2,
int spacedim>
671 const BlockVectorType2 & src)
const
673 const unsigned int n_blocks = src.n_blocks();
679 const unsigned int min_level = dst.min_level();
680 const unsigned int max_level = dst.max_level();
684 dst[
level].reinit(n_blocks);
688 min_level, max_level);
692 const unsigned int data_block = same_for_all ? 0 :
b;
693 get_matrix_free_transfer(data_block)
694 .copy_to_mg(*dof_handler[b], dst_non_block, src.block(b));
696 for (
unsigned int l = min_level;
l <= max_level; ++
l)
697 dst[l].block(b) = dst_non_block[
l];
704template <
int dim,
typename Number,
typename TransferType>
705template <
typename BlockVectorType2,
int spacedim>
709 BlockVectorType2 & dst,
715 "This object was initialized with support for usage with one "
716 "DoFHandler for each block, but this method assumes that "
717 "the same DoFHandler is used for all the blocks!"));
718 const std::vector<const DoFHandler<dim, spacedim> *> mg_dofs(dst.n_blocks(),
721 copy_from_mg(mg_dofs, dst, src);
724template <
int dim,
typename Number,
typename TransferType>
725template <
typename BlockVectorType2,
int spacedim>
729 BlockVectorType2 & dst,
733 const unsigned int n_blocks = dst.n_blocks();
739 const unsigned int min_level = src.min_level();
740 const unsigned int max_level = src.max_level();
742 for (
unsigned int l = min_level;
l <= max_level; ++
l)
747 min_level, max_level);
751 for (
unsigned int l = min_level;
l <= max_level; ++
l)
753 src_non_block[
l].reinit(src[l].block(b));
754 src_non_block[
l] = src[
l].block(b);
756 const unsigned int data_block = same_for_all ? 0 :
b;
757 get_matrix_free_transfer(data_block)
758 .copy_from_mg(*dof_handler[b], dst.block(b), src_non_block);
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
const IndexSet & locally_owned_mg_dofs(const unsigned int level) const
const Triangulation< dim, spacedim > & get_triangulation() const
types::global_dof_index n_dofs() const
unsigned int n_dofs_per_cell() const
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
bool restriction_is_additive(const unsigned int index) const
size_type n_elements() const
bool is_element(const size_type index) const
void update_ghost_values() const
void copy_locally_owned_data_from(const Vector< Number2, MemorySpace > &src)
void reinit(const size_type size, const bool omit_zeroing_entries=false)
SmartPointer< const MGConstrainedDoFs, MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > > > mg_constrained_dofs
std::function< void(const unsigned int, LinearAlgebra::distributed::Vector< Number > &)> initialize_dof_vector
virtual void prolongate(const unsigned int to_level, LinearAlgebra::distributed::BlockVector< Number > &dst, const LinearAlgebra::distributed::BlockVector< Number > &src) const override
void copy_to_mg(const std::vector< const DoFHandler< dim, spacedim > * > &dof_handler, MGLevelObject< LinearAlgebra::distributed::BlockVector< Number > > &dst, const BlockVectorType2 &src) const
virtual void prolongate_and_add(const unsigned int to_level, LinearAlgebra::distributed::BlockVector< Number > &dst, const LinearAlgebra::distributed::BlockVector< Number > &src) const override
void copy_to_mg(const DoFHandler< dim, spacedim > &dof_handler, MGLevelObject< LinearAlgebra::distributed::BlockVector< Number > > &dst, const BlockVectorType2 &src) const
virtual const TransferType & get_matrix_free_transfer(const unsigned int b) const =0
void copy_from_mg(const DoFHandler< dim, spacedim > &dof_handler, BlockVectorType2 &dst, const MGLevelObject< LinearAlgebra::distributed::BlockVector< Number > > &src) const
virtual void restrict_and_add(const unsigned int from_level, LinearAlgebra::distributed::BlockVector< Number > &dst, const LinearAlgebra::distributed::BlockVector< Number > &src) const override
MGTransferBlockMatrixFreeBase(const bool same_for_all)
void copy_from_mg(const std::vector< const DoFHandler< dim, spacedim > * > &dof_handler, BlockVectorType2 &dst, const MGLevelObject< LinearAlgebra::distributed::BlockVector< Number > > &src) const
static const bool supports_dof_handler_vector
std::vector< MGTransferMatrixFree< dim, Number > > matrix_free_transfer_vector
std::size_t memory_consumption() const
MGTransferBlockMatrixFree()=default
virtual ~MGTransferBlockMatrixFree() override=default
const MGTransferMatrixFree< dim, Number > & get_matrix_free_transfer(const unsigned int b) const override
void build(const DoFHandler< dim, dim > &dof_handler)
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
void do_restrict_add(const unsigned int from_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
void do_prolongate_add(const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
virtual ~MGTransferMatrixFree() override=default
std::vector< std::vector< unsigned int > > level_dof_indices
unsigned int n_components
virtual void prolongate(const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const override
unsigned int n_child_cell_dofs
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
std::vector< std::vector< std::vector< unsigned short > > > dirichlet_indices
void interpolate_to_mg(const DoFHandler< dim, spacedim > &dof_handler, MGLevelObject< LinearAlgebra::distributed::Vector< Number > > &dst, const BlockVectorType2 &src) const
MGLevelObject< std::shared_ptr< const Utilities::MPI::Partitioner > > vector_partitioners
void build(const DoFHandler< dim, dim > &dof_handler, const std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > &external_partitioners=std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > >())
AlignedVector< VectorizedArray< Number > > evaluation_data
std::vector< std::vector< std::pair< unsigned int, unsigned int > > > parent_child_connect
virtual void prolongate_and_add(const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const override
std::size_t memory_consumption() const
std::vector< AlignedVector< VectorizedArray< Number > > > weights_on_refined
std::vector< unsigned int > n_owned_level_cells
AlignedVector< VectorizedArray< Number > > prolongation_matrix_1d
bool element_is_continuous
virtual void restrict_and_add(const unsigned int from_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const override
virtual unsigned int n_global_levels() const
void update_ghost_values() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
#define DeclException0(Exception0)
static ::ExceptionBase & ExcNoProlongation()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
std::enable_if_t< IsBlockVector< VectorType >::value, void > collect_sizes(VectorType &vector)
std::enable_if_t< IsBlockVector< VectorType >::value, unsigned int > n_blocks(const VectorType &vector)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
void resolve_identity_constraints(const MGConstrainedDoFs *mg_constrained_dofs, const unsigned int level, std::vector< types::global_dof_index > &dof_indices)