PolyDEAL
 
Loading...
Searching...
No Matches
multigrid_amg.h
Go to the documentation of this file.
1// -----------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception OR LGPL-2.1-or-later
4// Copyright (C) XXXX - YYYY by the polyDEAL authors
5//
6// This file is part of the polyDEAL library.
7//
8// Detailed license information governing the source code
9// can be found in LICENSE.md at the top level directory.
10//
11// -----------------------------------------------------------------------------
12
13
14#ifndef multigrid_amg_h
15#define multigrid_amg_h
16
17
18#include <deal.II/base/config.h>
19
21
22#include <deal.II/fe/fe_dgq.h>
23
25
27
29
31
32namespace dealii
33{
41 template <int dim, typename OperatorType, typename Number = double>
43 {
44 private:
49 MPI_Comm communicator;
50
54 const OperatorType *mf_operator;
55
59 std::vector<TrilinosWrappers::SparseMatrix *> transfer_matrices;
60
64 std::vector<LinearOperatorMG<VectorType, VectorType>> level_operators;
65
70
71 public:
77 const OperatorType &mf_operator_,
78 const std::vector<TrilinosWrappers::SparseMatrix *> transfers_)
79 {
80 // Only DGQ discretizations are supported.
81 Assert(dynamic_cast<const FE_DGQ<dim> *>(
82 &mf_operator_.get_matrix_free()->get_dof_handler().get_fe()) !=
83 nullptr,
84 ExcNotImplemented());
86
87 // Check parallel layout is identical on every level
88 for (unsigned int l = 0; l < transfers_.size(); ++l)
89 Assert((mf_operator_.get_matrix_free()->get_locally_owned_set() ==
90 transfers_[l]->locally_owned_range_indices()),
91 ExcInternalError());
92
93 transfer_matrices.resize(transfers_.size());
94 level_operators.resize(transfers_.size());
95 // get communicator from first Trilinos matrix
96 communicator = transfers_[0]->get_mpi_communicator();
97
98 for (unsigned int l = 0; l < transfers_.size(); ++l)
99 {
100 // Set the pointer to the correct matrix
101 transfer_matrices[l] = transfers_[l];
102
103 // Define vmult-type lambdas for each linear operator.
104 level_operators[l].vmult = [this, l](VectorType &dst,
105 const VectorType &src) {
106 transfer_matrices[l]->vmult(dst, src);
107 };
108 level_operators[l].vmult_add = [this, l](VectorType &dst,
109 const VectorType &src) {
110 transfer_matrices[l]->vmult_add(dst, src);
111 };
112 level_operators[l].Tvmult = [this, l](VectorType &dst,
113 const VectorType &src) {
114 transfer_matrices[l]->Tvmult(dst, src);
115 };
116 level_operators[l].Tvmult_add = [this, l](VectorType &dst,
117 const VectorType &src) {
118 transfer_matrices[l]->Tvmult_add(dst, src);
119 };
120
121 // Inform each linear operator about the parallel layout. Use the
122 // given trilinos matrices.
123 level_operators[l].reinit_domain_vector = [this, l](VectorType &v,
124 bool) {
125 v.reinit(transfer_matrices[l]->locally_owned_domain_indices(),
127 };
128
129 level_operators[l].reinit_range_vector = [this, l](VectorType &v,
130 bool) {
131 v.reinit(transfer_matrices[l]->locally_owned_range_indices(),
133 };
134 }
135
136 // Do the same for the matrix-free object.
137 // First, set the pointer
138 mf_operator = &mf_operator_;
139
140 // Then, populate the corresponding lambdas (std::functions)
141 mf_linear_operator.vmult = [this](VectorType &dst,
142 const VectorType &src) {
143 mf_operator->vmult(dst, src);
144 };
145 mf_linear_operator.vmult_add = [this](VectorType &dst,
146 const VectorType &src) {
147 mf_operator->vmult_add(dst, src);
148 };
149 mf_linear_operator.Tvmult = [this](VectorType &dst,
150 const VectorType &src) {
151 mf_operator->Tvmult(dst, src);
152 };
153 mf_linear_operator.Tvmult_add = [this](VectorType &dst,
154 const VectorType &src) {
155 mf_operator->Tvmult_add(dst, src);
156 };
157
158 mf_linear_operator.reinit_domain_vector = [&](VectorType &v, bool) {};
159 mf_linear_operator.reinit_range_vector = [&](VectorType &v, bool) {
160 v.reinit(mf_operator->get_matrix_free()
161 ->get_dof_handler()
162 .locally_owned_dofs(),
164 };
165 }
166
175 void
178 {
179 Assert(mg_matrices.n_levels() > 1,
180 ExcMessage("Vector of matrices set to invalid size."));
181 Assert(!mf_linear_operator.is_null_operator, ExcInternalError());
182
184 const unsigned int min_level = mg_matrices.min_level();
185 const unsigned int max_level = mg_matrices.max_level();
186
187 mg_matrices[min_level] = mf_linear_operator; // finest level
188
189 // do the same, but using transfers to define level matrices
190 for (unsigned int l = min_level + 1; l < max_level; l++)
191 mg_matrices[l] = transpose_operator(level_operators[l - 1]) *
193 }
194 };
195
196
197
201 template <int dim, typename MatrixType, typename Number = double>
203 {
204 private:
209 MPI_Comm communicator;
210
214 // const MatrixType &fine_operator;
215
219 std::vector<const MatrixType *> transfer_matrices;
220
221
222
223 public:
228 AmgProjector(const std::vector<TrilinosWrappers::SparseMatrix> &transfers_)
229 {
231 // get communicator from first Trilinos matrix
232 communicator = transfers_[0].get_mpi_communicator();
233 transfer_matrices.resize(transfers_.size());
234
235 // Store pointers to fine operator and transfers
236 for (unsigned int l = 0; l < transfers_.size(); ++l)
237 {
238 // std::cout << "transferring matrix l= " << l << std::endl;
239 transfer_matrices[l] = &transfers_[l];
240 }
241 }
242
243
250 {
252
253 // get communicator from first Trilinos matrix
254 communicator = transfers_[0].get_mpi_communicator();
255 transfer_matrices.resize(transfers_.n_levels());
256
257 // Store pointers to fine operator and transfers
258 for (unsigned int l = 0; l < transfer_matrices.size(); ++l)
259 {
260 // std::cout << "transferring matrix l= " << l << std::endl;
261 transfer_matrices[l] = &transfers_[l];
262 }
263 }
264
265
266
274 void
276 MGLevelObject<std::unique_ptr<MatrixType>> &mg_matrices)
277 {
278 Assert(mg_matrices.n_levels() > 1,
279 ExcMessage("Vector of matrices set to invalid size."));
281
282 const unsigned int min_level = mg_matrices.min_level();
283 const unsigned int max_level = mg_matrices.max_level();
284
285 // mg_matrices[max_level].copy_from(fine_operator); // finest level
286 // MPI_Barrier(communicator);
287 // std::cout << "copied finest operator " << max_level << std::endl;
288
289 // do the same, but using transfers to define level matrices
290 // std::cout << "min level = " << min_level << std::endl;
291 for (unsigned int l = max_level; l-- > min_level;)
292 {
293 // Set parallel layout of intermediate operators AP
294 MatrixType level_operator;
295
296 mg_matrices[l] = std::make_unique<MatrixType>();
297
298 // First, compute AP
299 mg_matrices[l + 1]->mmult(
300 level_operator,
301 *transfer_matrices[l]); // result stored in level_operators[l]
302 // Multiply by the transpose
303 transfer_matrices[l]->Tmmult(*mg_matrices[l], level_operator);
304 }
305 }
306
307
308
317 void
319 MGLevelObject<std::unique_ptr<MatrixType>> &mg_matrices,
323 &multigrid_matrices_lo)
324 {
325 Assert(mg_matrices.n_levels() > 1,
326 ExcMessage("Vector of matrices set to invalid size."));
328
329 const unsigned int min_level = mg_matrices.min_level();
330 const unsigned int max_level = mg_matrices.max_level();
331
332 for (unsigned int l = max_level; l-- > min_level;)
333 {
334 // Set parallel layout of intermediate operators AP
335 MatrixType level_operator;
336
337 mg_matrices[l] = std::make_unique<MatrixType>();
338
339 // First, compute AP
340 mg_matrices[l + 1]->mmult(
341 level_operator,
342 *transfer_matrices[l]); // result stored in level_operators[l]
343 // Multiply by the transpose
344 transfer_matrices[l]->Tmmult(*mg_matrices[l], level_operator);
345
346 // Wrap every matrix into a linear operator now.
347 multigrid_matrices_lo[l] =
348 linear_operator_mg<VectorType, VectorType>(*mg_matrices[l]);
349 multigrid_matrices_lo[l].vmult =
350 [&mg_matrices, l](VectorType &dst, const VectorType &src) {
351 mg_matrices[l]->vmult(dst, src);
352 };
353 multigrid_matrices_lo[l].n_rows = mg_matrices[l]->m();
354 multigrid_matrices_lo[l].n_cols = mg_matrices[l]->n();
355 }
356 }
357
358
359 /*
360 * Same as above, but with a vector of matrices.
361 */
362 void
364 MGLevelObject<MatrixType> &mg_matrices,
368 &multigrid_matrices_lo)
369 {
370 Assert(mg_matrices.n_levels() > 1,
371 ExcMessage("Vector of matrices set to invalid size."));
373
374 const unsigned int min_level = mg_matrices.min_level();
375 const unsigned int max_level = mg_matrices.max_level();
376
377 for (unsigned int l = max_level; l-- > min_level;)
378 {
379 // Set intermediate operators "AP"
380 MatrixType level_operator;
381
382 // First, compute AP
383 mg_matrices[l + 1].mmult(
384 level_operator,
385 *transfer_matrices[l]); // result stored in level_operators[l]
386 // Multiply by the transpose
387 transfer_matrices[l]->Tmmult(mg_matrices[l], level_operator);
388 // Wrap every matrix into a linear operator now.
389 multigrid_matrices_lo[l] =
390 linear_operator_mg<VectorType, VectorType>(mg_matrices[l]);
391 multigrid_matrices_lo[l].vmult =
392 [&mg_matrices, l](VectorType &dst, const VectorType &src) {
393 mg_matrices[l].vmult(dst, src);
394 };
395 multigrid_matrices_lo[l].n_rows = mg_matrices[l].m();
396 multigrid_matrices_lo[l].n_cols = mg_matrices[l].n();
397 }
398 }
399
400
401
402 void
404 {
405 Assert(mg_matrices.n_levels() > 1,
406 ExcMessage("Vector of matrices set to invalid size."));
408
409 const unsigned int min_level = mg_matrices.min_level();
410 const unsigned int max_level = mg_matrices.max_level();
411
412 // mg_matrices[max_level].copy_from(fine_operator); // finest level
413 // MPI_Barrier(communicator);
414 // std::cout << "copied finest operator " << max_level << std::endl;
415
416 // do the same, but using transfers to define level matrices
417 // std::cout << "min level = " << min_level << std::endl;
418 for (unsigned int l = max_level; l-- > min_level;)
419 {
420 // Set parallel layout of intermediate operators AP
421 MatrixType level_operator(
422 transfer_matrices[l]->locally_owned_range_indices(),
423 transfer_matrices[l]->locally_owned_domain_indices(),
425
426 // First, compute AP
427 mg_matrices[l + 1].mmult(
428 level_operator,
429 *transfer_matrices[l]); // result stored in level_operators[l]
430 // Multiply by the transpose
431 transfer_matrices[l]->Tmmult(mg_matrices[l], level_operator);
432 }
433 }
434 };
435
436
437
444 template <int dim, typename VectorType>
445 class MGTransferAgglomeration : public MGTransferBase<VectorType>
446 {
447 public:
454 const std::vector<DoFHandler<dim> *> &dof_handlers);
455
462 const std::vector<DoFHandler<dim> *> &dof_handlers);
463
464
469 void
470 prolongate(const unsigned int to_level,
471 VectorType &dst,
472 const VectorType &src) const override;
473
474
478 void
479 prolongate_and_add(const unsigned int to_level,
480 VectorType &dst,
481 const VectorType &src) const override;
482
483
487 void
488 restrict_and_add(const unsigned int from_level,
489 VectorType &dst,
490 const VectorType &src) const override;
491
492
499 void
500 copy_to_mg(const DoFHandler<dim> &dof_handler,
502 const VectorType &src) const;
503
504
508 void
509 copy_from_mg(const DoFHandler<dim> &dof_handler,
510 VectorType &dst,
511 const MGLevelObject<VectorType> &src) const;
512
513
514 private:
520
524 std::vector<const DoFHandler<dim> *> dof_handlers;
525 };
526
527} // namespace dealii
528
529
530#endif
void reinit(const size_type size, const bool omit_zeroing_entries=false)
unsigned int max_level() const
unsigned int min_level() const
unsigned int n_levels() const
void compute_level_matrices(MGLevelObject< MatrixType > &mg_matrices)
std::vector< const MatrixType * > transfer_matrices
void compute_level_matrices_as_linear_operators(MGLevelObject< MatrixType > &mg_matrices, MGLevelObject< LinearOperatorMG< LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< Number > > > &multigrid_matrices_lo)
void compute_level_matrices_as_linear_operators(MGLevelObject< std::unique_ptr< MatrixType > > &mg_matrices, MGLevelObject< LinearOperatorMG< LinearAlgebra::distributed::Vector< Number >, LinearAlgebra::distributed::Vector< Number > > > &multigrid_matrices_lo)
LinearAlgebra::distributed::Vector< Number > VectorType
AmgProjector(const MGLevelObject< TrilinosWrappers::SparseMatrix > &transfers_)
void compute_level_matrices(MGLevelObject< std::unique_ptr< MatrixType > > &mg_matrices)
AmgProjector(const std::vector< TrilinosWrappers::SparseMatrix > &transfers_)
void prolongate(const unsigned int to_level, VectorType &dst, const VectorType &src) const override
void copy_from_mg(const DoFHandler< dim > &dof_handler, VectorType &dst, const MGLevelObject< VectorType > &src) const
MGTransferAgglomeration(const MGLevelObject< TrilinosWrappers::SparseMatrix * > &transfer_matrices, const std::vector< DoFHandler< dim > * > &dof_handlers)
MGLevelObject< ObserverPointer< TrilinosWrappers::SparseMatrix > > transfer_matrices
void restrict_and_add(const unsigned int from_level, VectorType &dst, const VectorType &src) const override
std::vector< const DoFHandler< dim > * > dof_handlers
void prolongate_and_add(const unsigned int to_level, VectorType &dst, const VectorType &src) const override
void copy_to_mg(const DoFHandler< dim > &dof_handler, MGLevelObject< VectorType > &dst, const VectorType &src) const
std::vector< TrilinosWrappers::SparseMatrix * > transfer_matrices
const OperatorType * mf_operator
std::vector< LinearOperatorMG< VectorType, VectorType > > level_operators
LinearAlgebra::distributed::Vector< Number > VectorType
LinearOperatorMG< VectorType, VectorType > mf_linear_operator
void compute_level_matrices(MGLevelObject< LinearOperatorMG< VectorType, VectorType > > &mg_matrices)
MatrixFreeProjector(const OperatorType &mf_operator_, const std::vector< TrilinosWrappers::SparseMatrix * > transfers_)
#define Assert(cond, exc)
LinearOperator< Domain, Range, Payload > transpose_operator(const LinearOperator< Range, Domain, Payload > &op)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)