|
Reference documentation for deal.II version 9.2.0
|
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\)
\(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\)
\(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\)
\(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Go to the documentation of this file.
16 #ifndef dealii_mg_transfer_matrix_free_h
17 #define dealii_mg_transfer_matrix_free_h
56 template <
int dim,
typename Number>
107 const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
108 &external_partitioners =
109 std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>());
127 const unsigned int to_level,
151 const unsigned int from_level,
165 template <
typename Number2,
int spacedim>
226 std::vector<std::vector<std::pair<unsigned int, unsigned int>>>
278 template <
int degree>
281 const unsigned int to_level,
288 template <
int degree>
312 template <
int dim,
typename Number>
314 :
public MGTransferBase<LinearAlgebra::distributed::BlockVector<Number>>
333 const std::vector<MGConstrainedDoFs> &mg_constrained_dofs);
351 const std::vector<MGConstrainedDoFs> &mg_constrained_dofs);
387 const unsigned int to_level,
411 const unsigned int from_level,
425 template <
typename Number2,
int spacedim>
435 template <
typename Number2,
int spacedim>
445 template <
typename Number2,
int spacedim>
456 template <
typename Number2,
int spacedim>
497 template <
int dim,
typename Number>
498 template <
typename Number2,
int spacedim>
505 const unsigned int min_level = dst.min_level();
506 const unsigned int max_level = dst.max_level();
524 relevant_dofs[
level]);
525 if (dst[
level].size() !=
527 dst[
level].local_size() !=
530 relevant_dofs[
level],
537 this->copy_to_mg(dof_handler, dst, src,
true);
542 dst[max_level].update_ghost_values();
549 relevant_dofs[
level],
551 ghosted_vector = dst[
level];
552 ghosted_vector.update_ghost_values();
557 std::vector<types::global_dof_index> dof_indices(fe.
dofs_per_cell);
561 for (; cell != endc; ++cell)
562 if (cell->is_locally_owned_on_level())
567 if (cell->is_active())
570 std::fill(dof_values_coarse.begin(), dof_values_coarse.end(), 0.);
571 for (
unsigned int child = 0; child < cell->n_children(); ++child)
573 cell->child(child)->get_mg_dof_indices(dof_indices);
575 dof_values_fine(i) = ghosted_vector(dof_indices[i]);
577 .
vmult(tmp, dof_values_fine);
580 dof_values_coarse[i] += tmp[i];
581 else if (tmp(i) != 0.)
582 dof_values_coarse[i] = tmp[i];
584 cell->get_mg_dof_indices(dof_indices);
586 dst[
level - 1](dof_indices[i]) = dof_values_coarse[i];
590 dst[
level - 1].update_ghost_values();
596 template <
int dim,
typename Number>
597 template <
typename Number2,
int spacedim>
607 "This object was initialized with support for usage with one "
608 "DoFHandler for each block, but this method assumes that "
609 "the same DoFHandler is used for all the blocks!"));
610 const std::vector<const DoFHandler<dim, spacedim> *> mg_dofs(src.
n_blocks(),
613 copy_to_mg(mg_dofs, dst, src);
616 template <
int dim,
typename Number>
617 template <
typename Number2,
int spacedim>
630 const unsigned int min_level = dst.min_level();
631 const unsigned int max_level = dst.max_level();
640 for (
unsigned int i = 1; i <
n_blocks; ++i)
643 &(dof_handler[0]->get_triangulation())) == tria),
644 ExcMessage(
"The DoFHandler use different Triangulations!"));
647 do_reinit.
resize(min_level, max_level);
650 do_reinit[
level] =
false;
653 do_reinit[
level] =
true;
664 do_reinit[
level] =
true;
672 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];
708 template <
int dim,
typename Number>
709 template <
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);
724 template <
int dim,
typename Number>
725 template <
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 IndexSet & locally_owned_mg_dofs(const unsigned int level) const
void vmult(Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const
std::vector< AlignedVector< VectorizedArray< Number > > > weights_on_refined
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
void interpolate_to_mg(const DoFHandler< dim, spacedim > &dof_handler, MGLevelObject< LinearAlgebra::distributed::Vector< Number >> &dst, const LinearAlgebra::distributed::Vector< Number2 > &src) const
static const bool supports_dof_handler_vector
void do_restrict_add(const unsigned int from_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< 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
void copy_from_mg(const DoFHandler< dim, spacedim > &dof_handler, LinearAlgebra::distributed::BlockVector< Number2 > &dst, const MGLevelObject< LinearAlgebra::distributed::BlockVector< Number >> &src) const
virtual void prolongate(const unsigned int to_level, LinearAlgebra::distributed::BlockVector< Number > &dst, const LinearAlgebra::distributed::BlockVector< Number > &src) const override
AlignedVector< VectorizedArray< Number > > prolongation_matrix_1d
std::size_t memory_consumption() const
virtual ~MGTransferMatrixFree() override=default
size_type n_elements() const
std::enable_if< IsBlockVector< VectorType >::value, unsigned int >::type n_blocks(const VectorType &vector)
static ::ExceptionBase & ExcNoProlongation()
cell_iterator begin(const unsigned int level=0) const
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 >>())
std::vector< std::vector< unsigned int > > level_dof_indices
unsigned int n_child_cell_dofs
static ::ExceptionBase & ExcMessage(std::string arg1)
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
const Triangulation< dim, spacedim > & get_triangulation() const
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
virtual ~MGTransferBlockMatrixFree() override=default
const unsigned int dofs_per_cell
std::vector< std::vector< std::pair< unsigned int, unsigned int > > > parent_child_connect
#define DEAL_II_NAMESPACE_OPEN
unsigned int n_blocks() const
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
MGTransferBlockMatrixFree()=default
size_type local_size() const
#define AssertDimension(dim1, dim2)
std::vector< std::vector< std::vector< unsigned short > > > dirichlet_indices
MGLevelObject< std::shared_ptr< const Utilities::MPI::Partitioner > > vector_partitioners
unsigned int n_components
cell_iterator end() const
virtual void restrict_and_add(const unsigned int from_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const override
virtual const MPI_Comm & get_communicator() const
void resize(const unsigned int new_minlevel, const unsigned int new_maxlevel)
#define Assert(cond, exc)
bool restriction_is_additive(const unsigned int index) const
#define DeclException0(Exception0)
void reinit(const size_type size, const bool omit_zeroing_entries=false)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
virtual size_type size() const override
#define DEAL_II_NAMESPACE_CLOSE
std::size_t memory_consumption() const
BlockType & block(const unsigned int i)
std::vector< unsigned int > n_owned_level_cells
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
#define AssertThrow(cond, exc)
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
void do_prolongate_add(const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
bool element_is_continuous
virtual void prolongate(const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const override
AlignedVector< VectorizedArray< Number > > evaluation_data
void build(const DoFHandler< dim, dim > &dof_handler)
SmartPointer< const MGConstrainedDoFs, MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > > > mg_constrained_dofs