16#ifndef dealii_mg_transfer_matrix_free_h
17#define dealii_mg_transfer_matrix_free_h
53template <
int dim,
typename Number>
104 const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
105 &external_partitioners =
106 std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>());
124 const unsigned int to_level,
130 const unsigned int to_level,
154 const unsigned int from_level,
172 template <
typename Number2,
int spacedim>
233 std::vector<std::vector<std::pair<unsigned int, unsigned int>>>
285 template <
int degree>
288 const unsigned int to_level,
295 template <
int degree>
316template <
int dim,
typename Number>
318 :
public MGTransferBase<LinearAlgebra::distributed::BlockVector<Number>>
337 const std::vector<MGConstrainedDoFs> &mg_constrained_dofs);
355 const std::vector<MGConstrainedDoFs> &mg_constrained_dofs);
391 const unsigned int to_level,
397 const unsigned int to_level,
421 const unsigned int from_level,
435 template <
typename Number2,
int spacedim>
445 template <
typename Number2,
int spacedim>
455 template <
typename Number2,
int spacedim>
466 template <
typename Number2,
int spacedim>
507template <
int dim,
typename Number>
508template <
typename Number2,
int spacedim>
515 const unsigned int min_level = dst.min_level();
516 const unsigned int max_level = dst.max_level();
526 dst[
level].locally_owned_size() !=
528 dst[
level].reinit(this->vector_partitioners[
level]);
531 this->copy_to_mg(dof_handler, dst, src,
true);
538 dst[max_level].update_ghost_values();
544 if (dst[
level].get_partitioner().get() ==
545 this->vector_partitioners[
level].get())
549 ghosted_fine.
reinit(this->vector_partitioners[
level]);
552 input = &ghosted_fine;
558 std::vector<types::global_dof_index> dof_indices(fe.
n_dofs_per_cell());
560 if (cell->is_locally_owned_on_level())
565 if (cell->is_active())
568 std::fill(dof_values_coarse.begin(), dof_values_coarse.end(), 0.);
569 for (
unsigned int child = 0; child < cell->n_children(); ++child)
571 cell->child(child)->get_mg_dof_indices(dof_indices);
573 dof_values_fine(i) = (*input)(dof_indices[i]);
575 .vmult(tmp, dof_values_fine);
578 dof_values_coarse[i] += tmp[i];
579 else if (tmp(i) != 0.)
580 dof_values_coarse[i] = tmp[i];
582 cell->get_mg_dof_indices(dof_indices);
586 dst[
level - 1](dof_indices[i]) = dof_values_coarse[i];
595template <
int dim,
typename Number>
596template <
typename Number2,
int spacedim>
606 "This object was initialized with support for usage with one "
607 "DoFHandler for each block, but this method assumes that "
608 "the same DoFHandler is used for all the blocks!"));
609 const std::vector<const DoFHandler<dim, spacedim> *> mg_dofs(src.
n_blocks(),
612 copy_to_mg(mg_dofs, dst, src);
617template <
int dim,
typename Number>
618template <
typename Number2,
int spacedim>
631 const unsigned int min_level = dst.min_level();
632 const unsigned int max_level = dst.max_level();
641 for (
unsigned int i = 1; i <
n_blocks; ++i)
644 &(dof_handler[0]->get_triangulation())) == tria),
645 ExcMessage(
"The DoFHandler use different Triangulations!"));
648 do_reinit.
resize(min_level, max_level);
651 do_reinit[
level] =
false;
654 do_reinit[
level] =
true;
665 do_reinit[
level] =
true;
673 if (do_reinit[
level])
683 dst[
level].collect_sizes();
692 min_level, max_level);
696 for (
unsigned int l = min_level;
l <= max_level; ++
l)
697 dst_non_block[
l].
reinit(dst[
l].block(
b));
698 const unsigned int data_block = same_for_all ? 0 :
b;
699 matrix_free_transfer_vector[data_block].copy_to_mg(*dof_handler[
b],
703 for (
unsigned int l = min_level;
l <= max_level; ++
l)
704 dst[
l].block(
b) = dst_non_block[
l];
708template <
int dim,
typename Number>
709template <
typename Number2,
int spacedim>
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>
725template <
typename Number2,
int spacedim>
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));
756 const unsigned int data_block = same_for_all ? 0 :
b;
757 matrix_free_transfer_vector[data_block].copy_from_mg(*dof_handler[
b],
const FiniteElement< dim, spacedim > & get_fe(const unsigned int 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
MPI_Comm get_communicator() 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
SmartPointer< const MGConstrainedDoFs, MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > > > mg_constrained_dofs
void resize(const unsigned int new_minlevel, const unsigned int new_maxlevel, Args &&... args)
void copy_to_mg(const DoFHandler< dim, spacedim > &dof_handler, MGLevelObject< LinearAlgebra::distributed::BlockVector< Number > > &dst, const LinearAlgebra::distributed::BlockVector< Number2 > &src) const
std::vector< MGTransferMatrixFree< dim, Number > > matrix_free_transfer_vector
virtual void prolongate_and_add(const unsigned int to_level, LinearAlgebra::distributed::BlockVector< Number > &dst, const LinearAlgebra::distributed::BlockVector< Number > &src) const override
virtual void restrict_and_add(const unsigned int from_level, LinearAlgebra::distributed::BlockVector< Number > &dst, const LinearAlgebra::distributed::BlockVector< Number > &src) const override
static const bool supports_dof_handler_vector
virtual void prolongate(const unsigned int to_level, LinearAlgebra::distributed::BlockVector< Number > &dst, const LinearAlgebra::distributed::BlockVector< Number > &src) const override
std::size_t memory_consumption() const
MGTransferBlockMatrixFree()=default
void copy_to_mg(const std::vector< const DoFHandler< dim, spacedim > * > &dof_handler, MGLevelObject< LinearAlgebra::distributed::BlockVector< Number > > &dst, const LinearAlgebra::distributed::BlockVector< Number2 > &src) const
void copy_from_mg(const DoFHandler< dim, spacedim > &dof_handler, LinearAlgebra::distributed::BlockVector< Number2 > &dst, const MGLevelObject< LinearAlgebra::distributed::BlockVector< Number > > &src) const
virtual ~MGTransferBlockMatrixFree() override=default
void build(const DoFHandler< dim, dim > &dof_handler)
void copy_from_mg(const std::vector< const DoFHandler< dim, spacedim > * > &dof_handler, LinearAlgebra::distributed::BlockVector< Number2 > &dst, const MGLevelObject< LinearAlgebra::distributed::BlockVector< Number > > &src) const
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
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
void interpolate_to_mg(const DoFHandler< dim, spacedim > &dof_handler, MGLevelObject< LinearAlgebra::distributed::Vector< Number > > &dst, const LinearAlgebra::distributed::Vector< Number2 > &src) const
virtual unsigned int n_global_levels() 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()
unsigned int n_blocks() const
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
BlockType & block(const unsigned int i)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
void update_ghost_values() const
size_type locally_owned_size() const
virtual size_type size() const override
void copy_locally_owned_data_from(const Vector< Number2, MemorySpace > &src)
void reinit(const size_type size, const bool omit_zeroing_entries=false)
void update_ghost_values() const
std::enable_if< IsBlockVector< VectorType >::value, unsignedint >::type 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 reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)