deal.II version GIT relicensing-2013-g7f3fb24d6c 2024-10-21 19:30:00+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_mpi_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_mpi_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, typename MemorySpace>
126 ::LinearAlgebra::TpetraWrappers::Vector<Number, MemorySpace>>
127 {
130
131 static const bool requires_distributed_sparsity_pattern = false;
132
133 template <typename SparsityPatternType, int dim, int spacedim>
134 static void
135 reinit(Matrix &matrix,
136 Sparsity &,
137 int level,
138 const SparsityPatternType &sp,
140 {
141 const MPI_Comm communicator = dh.get_mpi_communicator();
142
143 matrix.reinit(dh.locally_owned_mg_dofs(level + 1),
145 sp,
146 communicator,
147 true);
148 }
149 };
150# endif
151
152 template <>
154 {
157
158 static const bool requires_distributed_sparsity_pattern = false;
159
160 template <typename SparsityPatternType, int dim, int spacedim>
161 static void
162 reinit(Matrix &matrix,
163 Sparsity &,
164 int level,
165 const SparsityPatternType &sp,
167 {
168 const MPI_Comm communicator = dh.get_mpi_communicator();
169
170 matrix.reinit(dh.locally_owned_mg_dofs(level + 1),
172 sp,
173 communicator,
174 true);
175 }
176 };
177# endif
178
179#else
180 // ! DEAL_II_WITH_TRILINOS
181 template <typename Number>
182 struct MatrixSelector<LinearAlgebra::distributed::Vector<Number>>
183 {
186
187 static const bool requires_distributed_sparsity_pattern = false;
188
189 template <typename SparsityPatternType, int dim, int spacedim>
190 static void
191 reinit(Matrix &,
192 Sparsity &,
193 int,
194 const SparsityPatternType &,
196 {
198 false,
200 "ERROR: MGTransferPrebuilt with LinearAlgebra::distributed::Vector currently "
201 "needs deal.II to be configured with Trilinos."));
202 }
203 };
204
205#endif
206
207#ifdef DEAL_II_WITH_PETSC
208 template <>
210 {
213
215
216 template <typename SparsityPatternType, int dim, int spacedim>
217 static void
218 reinit(Matrix &matrix,
219 Sparsity &,
220 int level,
221 const SparsityPatternType &sp,
223 {
224 const MPI_Comm communicator = dh.get_mpi_communicator();
225
226 // Reinit PETSc matrix
227 matrix.reinit(dh.locally_owned_mg_dofs(level + 1),
229 sp,
230 communicator);
231 }
232 };
233#endif
234} // namespace internal
235
236/*
237 * MGTransferBase is defined in mg_base.h
238 */
239
249template <typename VectorType>
250class MGLevelGlobalTransfer : public MGTransferBase<VectorType>
251{
252public:
256 void
257 clear();
258
265 template <int dim, class InVector, int spacedim>
266 void
269 const InVector &src) const;
270
278 template <int dim, class OutVector, int spacedim>
279 void
281 OutVector &dst,
282 const MGLevelObject<VectorType> &src) const;
283
289 template <int dim, class OutVector, int spacedim>
290 void
292 OutVector &dst,
293 const MGLevelObject<VectorType> &src) const;
294
310 void
311 set_component_to_block_map(const std::vector<unsigned int> &map);
312
316 std::size_t
317 memory_consumption() const;
318
322 void
323 print_indices(std::ostream &os) const;
324
325protected:
329 template <int dim, int spacedim>
330 void
332 const DoFHandler<dim, spacedim> &dof_handler);
333
337 std::vector<types::global_dof_index> sizes;
338
346 std::vector<
347 std::vector<std::pair<types::global_dof_index, types::global_dof_index>>>
349
357 std::vector<
358 std::vector<std::pair<types::global_dof_index, types::global_dof_index>>>
360
368 std::vector<
369 std::vector<std::pair<types::global_dof_index, types::global_dof_index>>>
371
379
384 std::vector<unsigned int> component_to_block_map;
385
390
391private:
395 template <int dim, int spacedim>
396 void
397 assert_built(const DoFHandler<dim, spacedim> &dof_handler) const;
398};
399
400
401
410template <typename Number>
411class MGLevelGlobalTransfer<LinearAlgebra::distributed::Vector<Number>>
412 : public MGTransferBase<LinearAlgebra::distributed::Vector<Number>>
413{
414public:
418 void
419 clear();
420
427 template <int dim, typename Number2, int spacedim>
428 void
432
440 template <int dim, typename Number2, int spacedim>
441 void
443 const DoFHandler<dim, spacedim> &dof_handler,
446
452 template <int dim, typename Number2, int spacedim>
453 void
455 const DoFHandler<dim, spacedim> &dof_handler,
458
474 void
475 set_component_to_block_map(const std::vector<unsigned int> &map);
476
480 std::size_t
481 memory_consumption() const;
482
486 void
487 print_indices(std::ostream &os) const;
488
489protected:
494 template <int dim, typename Number2, int spacedim>
495 void
499 const bool solution_transfer) const;
500
504 template <int dim, int spacedim>
505 void
507 const DoFHandler<dim, spacedim> &dof_handler);
508
512 std::vector<types::global_dof_index> sizes;
513
522 std::vector<Table<2, unsigned int>> copy_indices;
523
527 std::vector<Table<2, unsigned int>> solution_copy_indices;
528
536 std::vector<Table<2, unsigned int>> copy_indices_global_mine;
537
545 std::vector<Table<2, unsigned int>> copy_indices_level_mine;
546
550 std::vector<Table<2, unsigned int>> solution_copy_indices_level_mine;
551
559
567
572 std::vector<unsigned int> component_to_block_map;
573
578
585
591
598
602 std::function<void(const unsigned int,
605
606private:
610 template <int dim, int spacedim>
611 void
612 assert_built(const DoFHandler<dim, spacedim> &dof_handler) const;
613};
614
615
616
627template <typename VectorType>
629{
630public:
636
642
646 virtual ~MGTransferPrebuilt() override = default;
647
651 void
653
657 void
658 clear();
659
664 template <int dim, int spacedim>
665 void
666 build(const DoFHandler<dim, spacedim> &dof_handler);
667
679 virtual void
680 prolongate(const unsigned int to_level,
681 VectorType &dst,
682 const VectorType &src) const override;
683
699 virtual void
700 restrict_and_add(const unsigned int from_level,
701 VectorType &dst,
702 const VectorType &src) const override;
703
708
713
717 std::size_t
718 memory_consumption() const;
719
723 void
724 print_matrices(std::ostream &os) const;
725
726private:
730 std::vector<
731 std::shared_ptr<typename internal::MatrixSelector<VectorType>::Sparsity>>
733
739 std::vector<
740 std::shared_ptr<typename internal::MatrixSelector<VectorType>::Matrix>>
742
747 std::vector<std::vector<bool>> interface_dofs;
748};
749
750
755
756#endif
const IndexSet & locally_owned_mg_dofs(const unsigned int level) const
MPI_Comm get_mpi_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
ObserverPointer< const MGConstrainedDoFs > mg_constrained_dofs
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)
ObserverPointer< const MGConstrainedDoFs > mg_constrained_dofs
std::vector< std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > copy_indices
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:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
unsigned int level
Definition grid_out.cc:4626
#define DeclException0(Exception0)
Definition exceptions.h:466
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