Reference documentation for deal.II version GIT relicensing-218-g1f873f3929 2024-03-28 15:00:02+00:00
\(\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\}}\)
Loading...
Searching...
No Matches
mg_transfer.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2001 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_mg_transfer_h
16#define dealii_mg_transfer_h
17
18#include <deal.II/base/config.h>
19
21
23
25
35
38
39#include <memory>
40
41
43
44
45namespace internal
46{
47 template <typename VectorType>
49 {
52
53 static const bool requires_distributed_sparsity_pattern = false;
54
55 template <typename SparsityPatternType, int dim, int spacedim>
56 static void
57 reinit(Matrix &matrix,
58 Sparsity &sparsity,
59 int level,
60 const SparsityPatternType &sp,
62 {
63 sparsity.copy_from(sp);
64 (void)level;
65 matrix.reinit(sparsity);
66 }
67 };
68
69#ifdef DEAL_II_WITH_TRILINOS
70 template <typename Number>
71 struct MatrixSelector<LinearAlgebra::distributed::Vector<Number>>
72 {
75
76 static const bool requires_distributed_sparsity_pattern = false;
77
78 template <typename SparsityPatternType, int dim, int spacedim>
79 static void
80 reinit(Matrix &matrix,
81 Sparsity &,
82 int level,
83 const SparsityPatternType &sp,
85 {
86 const MPI_Comm communicator = dh.get_communicator();
87
88 matrix.reinit(dh.locally_owned_mg_dofs(level + 1),
90 sp,
91 communicator,
92 true);
93 }
94 };
95
96 template <>
98 {
101
102 static const bool requires_distributed_sparsity_pattern = false;
103
104 template <typename SparsityPatternType, int dim, int spacedim>
105 static void
106 reinit(Matrix &matrix,
107 Sparsity &,
108 int level,
109 const SparsityPatternType &sp,
111 {
112 const MPI_Comm communicator = dh.get_communicator();
113
114 matrix.reinit(dh.locally_owned_mg_dofs(level + 1),
116 sp,
117 communicator,
118 true);
119 }
120 };
121
122# ifdef DEAL_II_WITH_MPI
123# ifdef DEAL_II_TRILINOS_WITH_TPETRA
124 template <typename Number>
126 {
129
130 static const bool requires_distributed_sparsity_pattern = false;
131
132 template <typename SparsityPatternType, int dim, int spacedim>
133 static void
134 reinit(Matrix &matrix,
135 Sparsity &,
136 int level,
137 const SparsityPatternType &sp,
139 {
140 const MPI_Comm communicator = dh.get_communicator();
141
142 matrix.reinit(dh.locally_owned_mg_dofs(level + 1),
144 sp,
145 communicator,
146 true);
147 }
148 };
149# endif
150
151 template <>
153 {
156
157 static const bool requires_distributed_sparsity_pattern = false;
158
159 template <typename SparsityPatternType, int dim, int spacedim>
160 static void
161 reinit(Matrix &matrix,
162 Sparsity &,
163 int level,
164 const SparsityPatternType &sp,
166 {
167 const MPI_Comm communicator = dh.get_communicator();
168
169 matrix.reinit(dh.locally_owned_mg_dofs(level + 1),
171 sp,
172 communicator,
173 true);
174 }
175 };
176# endif
177
178#else
179 // ! DEAL_II_WITH_TRILINOS
180 template <typename Number>
181 struct MatrixSelector<LinearAlgebra::distributed::Vector<Number>>
182 {
185
186 static const bool requires_distributed_sparsity_pattern = false;
187
188 template <typename SparsityPatternType, int dim, int spacedim>
189 static void
190 reinit(Matrix &,
191 Sparsity &,
192 int,
193 const SparsityPatternType &,
195 {
197 false,
199 "ERROR: MGTransferPrebuilt with LinearAlgebra::distributed::Vector currently "
200 "needs deal.II to be configured with Trilinos."));
201 }
202 };
203
204#endif
205
206#ifdef DEAL_II_WITH_PETSC
207 template <>
209 {
212
214
215 template <typename SparsityPatternType, int dim, int spacedim>
216 static void
217 reinit(Matrix &matrix,
218 Sparsity &,
219 int level,
220 const SparsityPatternType &sp,
222 {
223 const MPI_Comm communicator = dh.get_communicator();
224
225 // Reinit PETSc matrix
226 matrix.reinit(dh.locally_owned_mg_dofs(level + 1),
228 sp,
229 communicator);
230 }
231 };
232#endif
233} // namespace internal
234
235/*
236 * MGTransferBase is defined in mg_base.h
237 */
238
248template <typename VectorType>
249class MGLevelGlobalTransfer : public MGTransferBase<VectorType>
250{
251public:
255 void
256 clear();
257
264 template <int dim, class InVector, int spacedim>
265 void
268 const InVector &src) const;
269
277 template <int dim, class OutVector, int spacedim>
278 void
280 OutVector &dst,
281 const MGLevelObject<VectorType> &src) const;
282
288 template <int dim, class OutVector, int spacedim>
289 void
291 OutVector &dst,
292 const MGLevelObject<VectorType> &src) const;
293
309 void
310 set_component_to_block_map(const std::vector<unsigned int> &map);
311
315 std::size_t
316 memory_consumption() const;
317
321 void
322 print_indices(std::ostream &os) const;
323
324protected:
328 template <int dim, int spacedim>
329 void
331 const DoFHandler<dim, spacedim> &dof_handler);
332
336 std::vector<types::global_dof_index> sizes;
337
345 std::vector<
346 std::vector<std::pair<types::global_dof_index, types::global_dof_index>>>
348
356 std::vector<
357 std::vector<std::pair<types::global_dof_index, types::global_dof_index>>>
359
367 std::vector<
368 std::vector<std::pair<types::global_dof_index, types::global_dof_index>>>
370
378
383 std::vector<unsigned int> component_to_block_map;
384
389
390private:
394 template <int dim, int spacedim>
395 void
396 assert_built(const DoFHandler<dim, spacedim> &dof_handler) const;
397};
398
399
400
409template <typename Number>
410class MGLevelGlobalTransfer<LinearAlgebra::distributed::Vector<Number>>
411 : public MGTransferBase<LinearAlgebra::distributed::Vector<Number>>
412{
413public:
417 void
418 clear();
419
426 template <int dim, typename Number2, int spacedim>
427 void
431
439 template <int dim, typename Number2, int spacedim>
440 void
442 const DoFHandler<dim, spacedim> &dof_handler,
445
451 template <int dim, typename Number2, int spacedim>
452 void
454 const DoFHandler<dim, spacedim> &dof_handler,
457
473 void
474 set_component_to_block_map(const std::vector<unsigned int> &map);
475
479 std::size_t
480 memory_consumption() const;
481
485 void
486 print_indices(std::ostream &os) const;
487
488protected:
493 template <int dim, typename Number2, int spacedim>
494 void
498 const bool solution_transfer) const;
499
503 template <int dim, int spacedim>
504 void
506 const DoFHandler<dim, spacedim> &dof_handler);
507
511 std::vector<types::global_dof_index> sizes;
512
521 std::vector<Table<2, unsigned int>> copy_indices;
522
526 std::vector<Table<2, unsigned int>> solution_copy_indices;
527
535 std::vector<Table<2, unsigned int>> copy_indices_global_mine;
536
544 std::vector<Table<2, unsigned int>> copy_indices_level_mine;
545
549 std::vector<Table<2, unsigned int>> solution_copy_indices_level_mine;
550
558
566
571 std::vector<unsigned int> component_to_block_map;
572
577
584
590
597
601 std::function<void(const unsigned int,
604
605private:
609 template <int dim, int spacedim>
610 void
611 assert_built(const DoFHandler<dim, spacedim> &dof_handler) const;
612};
613
614
615
626template <typename VectorType>
628{
629public:
635
641
645 virtual ~MGTransferPrebuilt() override = default;
646
650 void
652
656 void
657 clear();
658
663 template <int dim, int spacedim>
664 void
665 build(const DoFHandler<dim, spacedim> &dof_handler);
666
678 virtual void
679 prolongate(const unsigned int to_level,
680 VectorType &dst,
681 const VectorType &src) const override;
682
698 virtual void
699 restrict_and_add(const unsigned int from_level,
700 VectorType &dst,
701 const VectorType &src) const override;
702
707
712
716 std::size_t
717 memory_consumption() const;
718
722 void
723 print_matrices(std::ostream &os) const;
724
725private:
729 std::vector<
730 std::shared_ptr<typename internal::MatrixSelector<VectorType>::Sparsity>>
732
738 std::vector<
739 std::shared_ptr<typename internal::MatrixSelector<VectorType>::Matrix>>
741
746 std::vector<std::vector<bool>> interface_dofs;
747};
748
749
754
755#endif
const IndexSet & locally_owned_mg_dofs(const unsigned int level) const
MPI_Comm get_communicator() const
void copy_from_mg(const DoFHandler< dim, spacedim > &dof_handler, LinearAlgebra::distributed::Vector< Number2 > &dst, const MGLevelObject< LinearAlgebra::distributed::Vector< Number > > &src) const
std::vector< Table< 2, unsigned int > > solution_copy_indices_level_mine
void copy_to_mg(const DoFHandler< dim, spacedim > &dof_handler, MGLevelObject< LinearAlgebra::distributed::Vector< Number > > &dst, const LinearAlgebra::distributed::Vector< Number2 > &src) const
MGLevelObject< LinearAlgebra::distributed::Vector< Number > > ghosted_level_vector
void set_component_to_block_map(const std::vector< unsigned int > &map)
void copy_from_mg_add(const DoFHandler< dim, spacedim > &dof_handler, LinearAlgebra::distributed::Vector< Number2 > &dst, const MGLevelObject< LinearAlgebra::distributed::Vector< Number > > &src) const
LinearAlgebra::distributed::Vector< Number > solution_ghosted_global_vector
void copy_to_mg(const DoFHandler< dim, spacedim > &dof_handler, MGLevelObject< LinearAlgebra::distributed::Vector< Number > > &dst, const LinearAlgebra::distributed::Vector< Number2 > &src, const bool solution_transfer) const
void assert_built(const DoFHandler< dim, spacedim > &dof_handler) const
std::vector< Table< 2, unsigned int > > copy_indices_global_mine
LinearAlgebra::distributed::Vector< Number > ghosted_global_vector
std::function< void(const unsigned int, LinearAlgebra::distributed::Vector< Number > &)> initialize_dof_vector
std::vector< Table< 2, unsigned int > > copy_indices_level_mine
void print_indices(std::ostream &os) const
void assert_built(const DoFHandler< dim, spacedim > &dof_handler) const
std::vector< std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > copy_indices_level_mine
std::size_t memory_consumption() const
void set_component_to_block_map(const std::vector< unsigned int > &map)
std::vector< std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > copy_indices
SmartPointer< const MGConstrainedDoFs > mg_constrained_dofs
void fill_and_communicate_copy_indices(const DoFHandler< dim, spacedim > &dof_handler)
std::vector< unsigned int > component_to_block_map
void copy_from_mg_add(const DoFHandler< dim, spacedim > &dof_handler, OutVector &dst, const MGLevelObject< VectorType > &src) const
std::vector< types::global_dof_index > sizes
std::vector< std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > copy_indices_global_mine
void copy_from_mg(const DoFHandler< dim, spacedim > &dof_handler, OutVector &dst, const MGLevelObject< VectorType > &src) const
void copy_to_mg(const DoFHandler< dim, spacedim > &dof_handler, MGLevelObject< VectorType > &dst, const InVector &src) const
virtual ~MGTransferPrebuilt() override=default
std::vector< std::vector< bool > > interface_dofs
virtual void restrict_and_add(const unsigned int from_level, VectorType &dst, const VectorType &src) const override
void build(const DoFHandler< dim, spacedim > &dof_handler)
std::vector< std::shared_ptr< typename internal::MatrixSelector< VectorType >::Sparsity > > prolongation_sparsities
std::size_t memory_consumption() const
MGTransferPrebuilt()=default
virtual void prolongate(const unsigned int to_level, VectorType &dst, const VectorType &src) const override
std::vector< std::shared_ptr< typename internal::MatrixSelector< VectorType >::Matrix > > prolongation_matrices
void print_matrices(std::ostream &os) const
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
void copy_from(const size_type n_rows, const size_type n_cols, const ForwardIterator begin, const ForwardIterator end)
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
unsigned int level
Definition grid_out.cc:4616
#define DeclException0(Exception0)
Definition exceptions.h:471
static ::ExceptionBase & ExcMatricesNotBuilt()
static ::ExceptionBase & ExcNoProlongation()
static ::ExceptionBase & ExcNotImplemented()
#define AssertThrow(cond, exc)
static void reinit(Matrix &matrix, Sparsity &, int level, const SparsityPatternType &sp, const DoFHandler< dim, spacedim > &dh)
Definition mg_transfer.h:80
static void reinit(Matrix &matrix, Sparsity &, int level, const SparsityPatternType &sp, const DoFHandler< dim, spacedim > &dh)
static void reinit(Matrix &matrix, Sparsity &, int level, const SparsityPatternType &sp, const DoFHandler< dim, spacedim > &dh)
static void reinit(Matrix &matrix, Sparsity &, int level, const SparsityPatternType &sp, const DoFHandler< dim, spacedim > &dh)
static void reinit(Matrix &matrix, Sparsity &, int level, const SparsityPatternType &sp, const DoFHandler< dim, spacedim > &dh)
static const bool requires_distributed_sparsity_pattern
Definition mg_transfer.h:53
static void reinit(Matrix &matrix, Sparsity &sparsity, int level, const SparsityPatternType &sp, const DoFHandler< dim, spacedim > &)
Definition mg_transfer.h:57
::SparseMatrix< typename VectorType::value_type > Matrix
Definition mg_transfer.h:51
::SparsityPattern Sparsity
Definition mg_transfer.h:50