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\}}\)
mg_transfer.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 2020 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_mg_transfer_h
17 #define dealii_mg_transfer_h
18 
19 #include <deal.II/base/config.h>
20 
22 
24 
26 
36 
39 
40 #include <memory>
41 
42 
44 
45 
46 namespace internal
47 {
48  template <typename VectorType>
50  {
53 
54  static const bool requires_distributed_sparsity_pattern = false;
55 
56  template <typename SparsityPatternType, typename DoFHandlerType>
57  static void
59  Sparsity & sparsity,
60  int level,
61  const SparsityPatternType &sp,
62  const DoFHandlerType &)
63  {
64  sparsity.copy_from(sp);
65  (void)level;
66  matrix.reinit(sparsity);
67  }
68  };
69 
70 #ifdef DEAL_II_WITH_TRILINOS
71  template <typename Number>
73  {
76 
77  static const bool requires_distributed_sparsity_pattern = false;
78 
79  template <typename SparsityPatternType, typename DoFHandlerType>
80  static void
82  Sparsity &,
83  int level,
84  const SparsityPatternType &sp,
85  DoFHandlerType & dh)
86  {
87  const ::parallel::TriangulationBase<DoFHandlerType::dimension,
88  DoFHandlerType::space_dimension>
89  *dist_tria = dynamic_cast<const ::parallel::TriangulationBase<
90  DoFHandlerType::dimension,
91  DoFHandlerType::space_dimension> *>(&(dh.get_triangulation()));
92  MPI_Comm communicator =
93  dist_tria != nullptr ? dist_tria->get_communicator() : MPI_COMM_SELF;
94 
95  matrix.reinit(dh.locally_owned_mg_dofs(level + 1),
96  dh.locally_owned_mg_dofs(level),
97  sp,
98  communicator,
99  true);
100  }
101  };
102 
103  template <>
105  {
108 
109  static const bool requires_distributed_sparsity_pattern = false;
110 
111  template <typename SparsityPatternType, typename DoFHandlerType>
112  static void
114  Sparsity &,
115  int level,
116  const SparsityPatternType &sp,
117  DoFHandlerType & dh)
118  {
119  const ::parallel::TriangulationBase<DoFHandlerType::dimension,
120  DoFHandlerType::space_dimension>
121  *dist_tria = dynamic_cast<const ::parallel::TriangulationBase<
122  DoFHandlerType::dimension,
123  DoFHandlerType::space_dimension> *>(&(dh.get_triangulation()));
124  MPI_Comm communicator =
125  dist_tria != nullptr ? dist_tria->get_communicator() : MPI_COMM_SELF;
126  matrix.reinit(dh.locally_owned_mg_dofs(level + 1),
127  dh.locally_owned_mg_dofs(level),
128  sp,
129  communicator,
130  true);
131  }
132  };
133 
134 # ifdef DEAL_II_WITH_MPI
135 # ifdef DEAL_II_TRILINOS_WITH_TPETRA
136  template <typename Number>
138  {
141 
142  static const bool requires_distributed_sparsity_pattern = false;
143 
144  template <typename SparsityPatternType, typename DoFHandlerType>
145  static void
147  Sparsity &,
148  int level,
149  const SparsityPatternType &sp,
150  DoFHandlerType & dh)
151  {
152  const ::parallel::TriangulationBase<DoFHandlerType::dimension,
153  DoFHandlerType::space_dimension>
154  *dist_tria = dynamic_cast<const ::parallel::TriangulationBase<
155  DoFHandlerType::dimension,
156  DoFHandlerType::space_dimension> *>(&(dh.get_triangulation()));
157  MPI_Comm communicator =
158  dist_tria != nullptr ? dist_tria->get_communicator() : MPI_COMM_SELF;
159  matrix.reinit(dh.locally_owned_mg_dofs(level + 1),
160  dh.locally_owned_mg_dofs(level),
161  sp,
162  communicator,
163  true);
164  }
165  };
166 # endif
167 
168  template <>
170  {
173 
174  static const bool requires_distributed_sparsity_pattern = false;
175 
176  template <typename SparsityPatternType, typename DoFHandlerType>
177  static void
179  Sparsity &,
180  int level,
181  const SparsityPatternType &sp,
182  DoFHandlerType & dh)
183  {
184  const ::parallel::TriangulationBase<DoFHandlerType::dimension,
185  DoFHandlerType::space_dimension>
186  *dist_tria = dynamic_cast<const ::parallel::TriangulationBase<
187  DoFHandlerType::dimension,
188  DoFHandlerType::space_dimension> *>(&(dh.get_triangulation()));
189  MPI_Comm communicator =
190  dist_tria != nullptr ? dist_tria->get_communicator() : MPI_COMM_SELF;
191  matrix.reinit(dh.locally_owned_mg_dofs(level + 1),
192  dh.locally_owned_mg_dofs(level),
193  sp,
194  communicator,
195  true);
196  }
197  };
198 # endif
199 
200 #else
201  // ! DEAL_II_WITH_TRILINOS
202  template <typename Number>
203  struct MatrixSelector<LinearAlgebra::distributed::Vector<Number>>
204  {
205  using Sparsity = ::SparsityPattern;
207 
208  static const bool requires_distributed_sparsity_pattern = false;
209 
210  template <typename SparsityPatternType, typename DoFHandlerType>
211  static void
212  reinit(Matrix &,
213  Sparsity &,
214  int,
215  const SparsityPatternType &,
216  const DoFHandlerType &)
217  {
218  AssertThrow(
219  false,
221  "ERROR: MGTransferPrebuilt with LinearAlgebra::distributed::Vector currently "
222  "needs deal.II to be configured with Trilinos."));
223  }
224  };
225 
226 #endif
227 
228 #ifdef DEAL_II_WITH_PETSC
229  template <>
231  {
234 
235  static const bool requires_distributed_sparsity_pattern = true;
236 
237  template <typename SparsityPatternType, typename DoFHandlerType>
238  static void
240  Sparsity &,
241  int level,
242  const SparsityPatternType &sp,
243  const DoFHandlerType & dh)
244  {
245  const ::parallel::TriangulationBase<DoFHandlerType::dimension,
246  DoFHandlerType::space_dimension>
247  *dist_tria = dynamic_cast<const ::parallel::TriangulationBase<
248  DoFHandlerType::dimension,
249  DoFHandlerType::space_dimension> *>(&(dh.get_triangulation()));
250  MPI_Comm communicator =
251  dist_tria != nullptr ? dist_tria->get_communicator() : MPI_COMM_SELF;
252  // Reinit PETSc matrix
253  matrix.reinit(dh.locally_owned_mg_dofs(level + 1),
254  dh.locally_owned_mg_dofs(level),
255  sp,
256  communicator);
257  }
258  };
259 #endif
260 } // namespace internal
261 
262 /*
263  * MGTransferBase is defined in mg_base.h
264  */
265 
268 
269 
270 
278 template <typename VectorType>
279 class MGLevelGlobalTransfer : public MGTransferBase<VectorType>
280 {
281 public:
285  void
286  clear();
287 
294  template <int dim, class InVector, int spacedim>
295  void
296  copy_to_mg(const DoFHandler<dim, spacedim> &dof_handler,
298  const InVector & src) const;
299 
307  template <int dim, class OutVector, int spacedim>
308  void
309  copy_from_mg(const DoFHandler<dim, spacedim> &dof_handler,
310  OutVector & dst,
311  const MGLevelObject<VectorType> &src) const;
312 
318  template <int dim, class OutVector, int spacedim>
319  void
320  copy_from_mg_add(const DoFHandler<dim, spacedim> &dof_handler,
321  OutVector & dst,
322  const MGLevelObject<VectorType> &src) const;
323 
339  void
340  set_component_to_block_map(const std::vector<unsigned int> &map);
341 
345  std::size_t
346  memory_consumption() const;
347 
351  void
352  print_indices(std::ostream &os) const;
353 
354 protected:
358  template <int dim, int spacedim>
359  void
361  const DoFHandler<dim, spacedim> &dof_handler);
362 
366  std::vector<types::global_dof_index> sizes;
367 
375  std::vector<
376  std::vector<std::pair<types::global_dof_index, types::global_dof_index>>>
378 
386  std::vector<
387  std::vector<std::pair<types::global_dof_index, types::global_dof_index>>>
389 
397  std::vector<
398  std::vector<std::pair<types::global_dof_index, types::global_dof_index>>>
400 
408 
413  std::vector<unsigned int> component_to_block_map;
414 
420 
421 private:
425  template <int dim, int spacedim>
426  void
427  assert_built(const DoFHandler<dim, spacedim> &dof_handler) const;
428 };
429 
430 
431 
443 template <typename Number>
445  : public MGTransferBase<LinearAlgebra::distributed::Vector<Number>>
446 {
447 public:
451  void
452  clear();
453 
460  template <int dim, typename Number2, int spacedim>
461  void
462  copy_to_mg(const DoFHandler<dim, spacedim> &dof_handler,
465 
473  template <int dim, typename Number2, int spacedim>
474  void
475  copy_from_mg(
476  const DoFHandler<dim, spacedim> & dof_handler,
479 
485  template <int dim, typename Number2, int spacedim>
486  void
488  const DoFHandler<dim, spacedim> & dof_handler,
491 
507  void
508  set_component_to_block_map(const std::vector<unsigned int> &map);
509 
513  std::size_t
514  memory_consumption() const;
515 
519  void
520  print_indices(std::ostream &os) const;
521 
522 protected:
527  template <int dim, typename Number2, int spacedim>
528  void
529  copy_to_mg(const DoFHandler<dim, spacedim> &dof_handler,
532  const bool solution_transfer) const;
533 
537  template <int dim, int spacedim>
538  void
540  const DoFHandler<dim, spacedim> &dof_handler);
541 
545  std::vector<types::global_dof_index> sizes;
546 
555  std::vector<Table<2, unsigned int>> copy_indices;
556 
560  std::vector<Table<2, unsigned int>> solution_copy_indices;
561 
569  std::vector<Table<2, unsigned int>> copy_indices_global_mine;
570 
574  std::vector<Table<2, unsigned int>> solution_copy_indices_global_mine;
575 
583  std::vector<Table<2, unsigned int>> copy_indices_level_mine;
584 
588  std::vector<Table<2, unsigned int>> solution_copy_indices_level_mine;
589 
597 
605 
610  std::vector<unsigned int> component_to_block_map;
611 
615  SmartPointer<
616  const MGConstrainedDoFs,
619 
626 
632 
639 
645 
646 private:
650  template <int dim, int spacedim>
651  void
652  assert_built(const DoFHandler<dim, spacedim> &dof_handler) const;
653 };
654 
655 
656 
670 template <typename VectorType>
671 class MGTransferPrebuilt : public MGLevelGlobalTransfer<VectorType>
672 {
673 public:
678  MGTransferPrebuilt() = default;
679 
685 
689  virtual ~MGTransferPrebuilt() override = default;
690 
694  void
696 
700  void
701  clear();
702 
707  template <int dim, int spacedim>
708  void
709  build(const DoFHandler<dim, spacedim> &dof_handler);
710 
716  template <int dim, int spacedim>
717  DEAL_II_DEPRECATED void
718  build_matrices(const DoFHandler<dim, spacedim> &dof_handler);
719 
731  virtual void
732  prolongate(const unsigned int to_level,
733  VectorType & dst,
734  const VectorType & src) const override;
735 
751  virtual void
752  restrict_and_add(const unsigned int from_level,
753  VectorType & dst,
754  const VectorType & src) const override;
755 
760 
765 
769  std::size_t
770  memory_consumption() const;
771 
775  void
776  print_matrices(std::ostream &os) const;
777 
778 private:
782  std::vector<
783  std::shared_ptr<typename internal::MatrixSelector<VectorType>::Sparsity>>
785 
791  std::vector<
792  std::shared_ptr<typename internal::MatrixSelector<VectorType>::Matrix>>
794 
799  std::vector<std::vector<bool>> interface_dofs;
800 };
801 
802 
807 
808 #endif
MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > >::copy_indices_global_mine
std::vector< Table< 2, unsigned int > > copy_indices_global_mine
Definition: mg_transfer.h:569
internal::MatrixSelector::Matrix
::SparseMatrix< typename VectorType::value_type > Matrix
Definition: mg_transfer.h:52
mg_base.h
sparse_matrix.h
TrilinosWrappers::MPI::Vector
Definition: trilinos_vector.h:400
LinearAlgebraDealII::Vector
Vector< double > Vector
Definition: generic_linear_algebra.h:43
MGTransferBase
Definition: mg_base.h:177
TrilinosWrappers::SparseMatrix
Definition: trilinos_sparse_matrix.h:515
LinearAlgebra
Definition: communication_pattern_base.h:30
MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > >::copy_indices_level_mine
std::vector< Table< 2, unsigned int > > copy_indices_level_mine
Definition: mg_transfer.h:583
LinearAlgebra::distributed::Vector< Number >
MGLevelGlobalTransfer
Definition: mg_transfer.h:279
MGConstrainedDoFs
Definition: mg_constrained_dofs.h:45
StandardExceptions::ExcNotImplemented
static ::ExceptionBase & ExcNotImplemented()
MGTransferPrebuilt::ExcNoProlongation
static ::ExceptionBase & ExcNoProlongation()
PETScWrappers::MPI::SparseMatrix
Definition: petsc_sparse_matrix.h:367
MGLevelGlobalTransfer::memory_consumption
std::size_t memory_consumption() const
Definition: mg_level_global_transfer.cc:143
MGLevelGlobalTransfer::mg_constrained_dofs
SmartPointer< const MGConstrainedDoFs, MGLevelGlobalTransfer< VectorType > > mg_constrained_dofs
Definition: mg_transfer.h:419
MGTransferPrebuilt::prolongate
virtual void prolongate(const unsigned int to_level, VectorType &dst, const VectorType &src) const override
Definition: mg_transfer_prebuilt.cc:81
MGLevelGlobalTransfer::print_indices
void print_indices(std::ostream &os) const
Definition: mg_level_global_transfer.cc:114
SparseMatrix
Definition: sparse_matrix.h:497
VectorType
petsc_sparse_matrix.h
MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > >::sizes
std::vector< types::global_dof_index > sizes
Definition: mg_transfer.h:545
trilinos_sparse_matrix.h
MPI_Comm
MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > >::copy_indices
std::vector< Table< 2, unsigned int > > copy_indices
Definition: mg_transfer.h:555
MGLevelGlobalTransfer::fill_and_communicate_copy_indices
void fill_and_communicate_copy_indices(const DoFHandler< dim, spacedim > &dof_handler)
Definition: mg_level_global_transfer.cc:51
block_sparsity_pattern.h
internal::MatrixSelector< LinearAlgebra::distributed::Vector< Number > >::reinit
static void reinit(Matrix &matrix, Sparsity &, int level, const SparsityPatternType &sp, DoFHandlerType &dh)
Definition: mg_transfer.h:81
MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > >::ghosted_global_vector
LinearAlgebra::distributed::Vector< Number > ghosted_global_vector
Definition: mg_transfer.h:625
DoFHandler
Definition: dof_handler.h:205
la_parallel_vector.h
TrilinosWrappers::SparsityPattern
Definition: trilinos_sparsity_pattern.h:279
MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > >::perform_renumbered_plain_copy
bool perform_renumbered_plain_copy
Definition: mg_transfer.h:604
MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > >::solution_copy_indices_level_mine
std::vector< Table< 2, unsigned int > > solution_copy_indices_level_mine
Definition: mg_transfer.h:588
MGLevelGlobalTransfer::copy_indices_level_mine
std::vector< std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > copy_indices_level_mine
Definition: mg_transfer.h:399
MGLevelGlobalTransfer::copy_indices_global_mine
std::vector< std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > copy_indices_global_mine
Definition: mg_transfer.h:388
level
unsigned int level
Definition: grid_out.cc:4355
internal::MatrixSelector::reinit
static void reinit(Matrix &matrix, Sparsity &sparsity, int level, const SparsityPatternType &sp, const DoFHandlerType &)
Definition: mg_transfer.h:58
internal::MatrixSelector<::TrilinosWrappers::MPI::Vector >::reinit
static void reinit(Matrix &matrix, Sparsity &, int level, const SparsityPatternType &sp, DoFHandlerType &dh)
Definition: mg_transfer.h:113
MGLevelGlobalTransfer::copy_from_mg_add
void copy_from_mg_add(const DoFHandler< dim, spacedim > &dof_handler, OutVector &dst, const MGLevelObject< VectorType > &src) const
MGTransferPrebuilt::~MGTransferPrebuilt
virtual ~MGTransferPrebuilt() override=default
MGTransferPrebuilt::prolongation_matrices
std::vector< std::shared_ptr< typename internal::MatrixSelector< VectorType >::Matrix > > prolongation_matrices
Definition: mg_transfer.h:793
MGTransferPrebuilt::build_matrices
void build_matrices(const DoFHandler< dim, spacedim > &dof_handler)
Definition: mg_transfer_prebuilt.cc:355
internal::MatrixSelector::requires_distributed_sparsity_pattern
static const bool requires_distributed_sparsity_pattern
Definition: mg_transfer.h:54
block_vector.h
MGTransferPrebuilt::initialize_constraints
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
Definition: mg_transfer_prebuilt.cc:59
MGTransferPrebuilt::MGTransferPrebuilt
MGTransferPrebuilt()=default
LinearAlgebra::EpetraWrappers::Vector
Definition: trilinos_epetra_vector.h:65
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > >::solution_copy_indices
std::vector< Table< 2, unsigned int > > solution_copy_indices
Definition: mg_transfer.h:560
LinearAlgebra::TpetraWrappers::Vector
Definition: trilinos_tpetra_vector.h:84
LAPACKSupport::matrix
@ matrix
Contents is actually a matrix.
Definition: lapack_support.h:60
MGTransferPrebuilt::memory_consumption
std::size_t memory_consumption() const
Definition: mg_transfer_prebuilt.cc:379
MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > >::component_to_block_map
std::vector< unsigned int > component_to_block_map
Definition: mg_transfer.h:610
DynamicSparsityPattern
Definition: dynamic_sparsity_pattern.h:323
internal::MatrixSelector::Sparsity
::SparsityPattern Sparsity
Definition: mg_transfer.h:51
SparsityPattern
Definition: sparsity_pattern.h:865
MGLevelGlobalTransfer::clear
void clear()
Definition: mg_level_global_transfer.cc:99
DEAL_II_DEPRECATED
#define DEAL_II_DEPRECATED
Definition: config.h:98
MGTransferPrebuilt
Definition: mg_transfer.h:671
MGTransferPrebuilt::build
void build(const DoFHandler< dim, spacedim > &dof_handler)
Definition: mg_transfer_prebuilt.cc:147
mg_level_object.h
MGTransferPrebuilt::clear
void clear()
Definition: mg_transfer_prebuilt.cc:69
internal::MatrixSelector<::PETScWrappers::MPI::Vector >::reinit
static void reinit(Matrix &matrix, Sparsity &, int level, const SparsityPatternType &sp, const DoFHandlerType &dh)
Definition: mg_transfer.h:239
PETScWrappers::MPI::Vector
Definition: petsc_vector.h:158
internal::MatrixSelector<::LinearAlgebra::TpetraWrappers::Vector< Number > >::reinit
static void reinit(Matrix &matrix, Sparsity &, int level, const SparsityPatternType &sp, DoFHandlerType &dh)
Definition: mg_transfer.h:146
mg_constrained_dofs.h
SmartPointer
Definition: smartpointer.h:68
internal::MatrixSelector<::LinearAlgebra::EpetraWrappers::Vector >::reinit
static void reinit(Matrix &matrix, Sparsity &, int level, const SparsityPatternType &sp, DoFHandlerType &dh)
Definition: mg_transfer.h:178
LinearAlgebraDealII::SparseMatrix
SparseMatrix< double > SparseMatrix
Definition: generic_linear_algebra.h:53
MGLevelGlobalTransfer::perform_plain_copy
bool perform_plain_copy
Definition: mg_transfer.h:407
MGLevelObject< VectorType >
MGLevelGlobalTransfer::set_component_to_block_map
void set_component_to_block_map(const std::vector< unsigned int > &map)
dof_handler.h
DeclException0
#define DeclException0(Exception0)
Definition: exceptions.h:473
MGTransferPrebuilt::prolongation_sparsities
std::vector< std::shared_ptr< typename internal::MatrixSelector< VectorType >::Sparsity > > prolongation_sparsities
Definition: mg_transfer.h:784
MGLevelGlobalTransfer::sizes
std::vector< types::global_dof_index > sizes
Definition: mg_transfer.h:366
affine_constraints.h
MGLevelGlobalTransfer::component_to_block_map
std::vector< unsigned int > component_to_block_map
Definition: mg_transfer.h:413
MGLevelGlobalTransfer::copy_to_mg
void copy_to_mg(const DoFHandler< dim, spacedim > &dof_handler, MGLevelObject< VectorType > &dst, const InVector &src) const
MGLevelGlobalTransfer::copy_indices
std::vector< std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > copy_indices
Definition: mg_transfer.h:377
MGTransferPrebuilt::ExcMatricesNotBuilt
static ::ExceptionBase & ExcMatricesNotBuilt()
MGTransferPrebuilt::interface_dofs
std::vector< std::vector< bool > > interface_dofs
Definition: mg_transfer.h:799
petsc_vector.h
vector_memory.h
parallel::TriangulationBase
Definition: tria_base.h:78
MGLevelGlobalTransfer::assert_built
void assert_built(const DoFHandler< dim, spacedim > &dof_handler) const
config.h
MGTransferPrebuilt::print_matrices
void print_matrices(std::ostream &os) const
Definition: mg_transfer_prebuilt.cc:365
internal
Definition: aligned_vector.h:369
internal::MatrixSelector
Definition: mg_transfer.h:49
MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > >::ghosted_level_vector
MGLevelObject< LinearAlgebra::distributed::Vector< Number > > ghosted_level_vector
Definition: mg_transfer.h:638
MGLevelGlobalTransfer::copy_from_mg
void copy_from_mg(const DoFHandler< dim, spacedim > &dof_handler, OutVector &dst, const MGLevelObject< VectorType > &src) const
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > >::solution_copy_indices_global_mine
std::vector< Table< 2, unsigned int > > solution_copy_indices_global_mine
Definition: mg_transfer.h:574
SparsityPattern::copy_from
void copy_from(const size_type n_rows, const size_type n_cols, const ForwardIterator begin, const ForwardIterator end)
AssertThrow
#define AssertThrow(cond, exc)
Definition: exceptions.h:1531
MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > >::perform_plain_copy
bool perform_plain_copy
Definition: mg_transfer.h:596
MGTransferPrebuilt::restrict_and_add
virtual void restrict_and_add(const unsigned int from_level, VectorType &dst, const VectorType &src) const override
Definition: mg_transfer_prebuilt.cc:103
MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > >::solution_ghosted_global_vector
LinearAlgebra::distributed::Vector< Number > solution_ghosted_global_vector
Definition: mg_transfer.h:631
tria_base.h
MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > >::mg_constrained_dofs
SmartPointer< const MGConstrainedDoFs, MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > > > mg_constrained_dofs
Definition: mg_transfer.h:618
MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > >::solution_ghosted_level_vector
MGLevelObject< LinearAlgebra::distributed::Vector< Number > > solution_ghosted_level_vector
Definition: mg_transfer.h:644