Reference documentation for deal.II version 9.6.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\}}\)
Loading...
Searching...
No Matches
petsc_precondition.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2004 - 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_petsc_precondition_h
16#define dealii_petsc_precondition_h
17
18
19#include <deal.II/base/config.h>
20
21#include <deal.II/base/point.h>
23
24#ifdef DEAL_II_WITH_PETSC
25
27
28# include <petscpc.h>
29
30# include <functional>
31
33
34
35
36namespace PETScWrappers
37{
38 // forward declarations
39# ifndef DOXYGEN
40 class MatrixBase;
41 class VectorBase;
42# endif
43
60 {
61 public:
65 explicit PreconditionBase(const MPI_Comm mpi_communicator);
66
72
76 virtual ~PreconditionBase();
77
82 void
83 clear();
84
88 void
89 vmult(VectorBase &dst, const VectorBase &src) const;
90
94 void
95 Tvmult(VectorBase &dst, const VectorBase &src) const;
96
101 void
102 setup();
103
107 const PC &
108 get_pc() const;
109
114 get_mpi_communicator() const;
115
116 protected:
120 PC pc;
121
126 void
128
132 void
134 };
135
136
137
149 {
150 public:
156 {};
157
163
169 const MatrixBase &matrix,
171
177 const MPI_Comm communicator,
179
186 void
187 initialize(const MatrixBase &matrix,
189
190 protected:
195
201 void
202 initialize();
203 };
204
205
206
231 {
232 public:
238 {};
239
245
251 const MatrixBase &matrix,
253
259 const MPI_Comm communicator,
261
262
269 void
270 initialize(const MatrixBase &matrix,
272
273 protected:
278
284 void
285 initialize();
286 };
287
288
289
299 {
300 public:
306 {
310 AdditionalData(const double omega = 1);
311
315 double omega;
316 };
317
323
328 PreconditionSOR(const MatrixBase &matrix,
330
337 void
338 initialize(const MatrixBase &matrix,
340
341 protected:
346 };
347
348
349
359 {
360 public:
366 {
370 explicit AdditionalData(const double omega = 1);
371
375 double omega;
376 };
377
383
388 PreconditionSSOR(const MatrixBase &matrix,
390
397 void
398 initialize(const MatrixBase &matrix,
400
401 protected:
406 };
407
417 {
418 public:
424 {
428 AdditionalData(const unsigned int levels = 0);
429
433 unsigned int levels;
434 };
435
441
446 PreconditionICC(const MatrixBase &matrix,
448
455 void
456 initialize(const MatrixBase &matrix,
458
459 protected:
464 };
465
466
467
477 {
478 public:
484 {
488 AdditionalData(const unsigned int levels = 0);
489
493 unsigned int levels;
494 };
495
501
506 PreconditionILU(const MatrixBase &matrix,
508
515 void
516 initialize(const MatrixBase &matrix,
518
519 protected:
524 };
525
526
527
542 {
543 public:
549 {
554 AdditionalData(const double pivoting = 1.e-6,
555 const double zero_pivot = 1.e-12,
556 const double damping = 0.0);
557
563 double pivoting;
564
570
575 double damping;
576 };
577
583
588 PreconditionLU(const MatrixBase &matrix,
590
597 void
598 initialize(const MatrixBase &matrix,
600
601 protected:
606 };
607
608
609
621 {
622 public:
628 {
650
656 const bool symmetric_operator = false,
657 const double strong_threshold = 0.25,
658 const double max_row_sum = 0.9,
659 const unsigned int aggressive_coarsening_num_levels = 0,
660 const bool output_details = false,
661 const RelaxationType relaxation_type_up = RelaxationType::SORJacobi,
662 const RelaxationType relaxation_type_down = RelaxationType::SORJacobi,
663 const RelaxationType relaxation_type_coarse =
665 const unsigned int n_sweeps_coarse = 1,
666 const double tol = 0.0,
667 const unsigned int max_iter = 1,
668 const bool w_cycle = false);
669
676
683
696
703
709
714
719
724
728 unsigned int n_sweeps_coarse;
729
733 double tol;
734
738 unsigned int max_iter;
739
745 };
746
752
758 const MatrixBase &matrix,
760
766 const MPI_Comm communicator,
768
769
776 void
777 initialize(const MatrixBase &matrix,
779
780 protected:
785
791 void
792 initialize();
793 };
794
795
796
817 {
818 public:
824 {
828 AdditionalData(const unsigned int symmetric = 1,
829 const unsigned int n_levels = 1,
830 const double threshold = 0.1,
831 const double filter = 0.05,
832 const bool output_details = false);
833
845 unsigned int symmetric;
846
853 unsigned int n_levels;
854
865 double threshold;
866
877 double filter;
878
884 };
885
886
887
893
899 const MatrixBase &matrix,
901
908 void
909 initialize(const MatrixBase &matrix,
911
912 private:
917 };
918
919
920
927 {
928 public:
934 {};
935
941
947 PreconditionNone(const MatrixBase &matrix,
949
957 void
958 initialize(const MatrixBase &matrix,
960
961 private:
966 };
967
991 template <int dim>
993 {
994 public:
1000 {
1005 AdditionalData(const bool use_vertices = true,
1006 const bool use_edges = false,
1007 const bool use_faces = false,
1008 const bool symmetric = false,
1009 const std::vector<Point<dim>> coords = {});
1010
1016
1023
1030
1035
1040 std::vector<Point<dim>> coords;
1041 };
1042
1048
1053 PreconditionBDDC(const MatrixBase &matrix,
1055
1060 PreconditionBDDC(const MPI_Comm communicator,
1062
1069 void
1070 initialize(const MatrixBase &matrix,
1072
1073 protected:
1078
1084 void
1085 initialize();
1086 };
1087
1089 {
1090 public:
1096
1100 PreconditionShell(const MatrixBase &matrix);
1101
1105 PreconditionShell(const MPI_Comm communicator);
1106
1115 std::function<void(VectorBase &dst, const VectorBase &src)> vmult;
1116
1125 std::function<void(VectorBase &dst, const VectorBase &src)> vmultT;
1126
1127 protected:
1132 void
1133 initialize(const MPI_Comm comm);
1134
1139 void
1140 initialize(const MatrixBase &matrix);
1141
1142 private:
1146 static PetscErrorCode
1147 pcapply(PC pc, Vec src, Vec dst);
1148
1152 static PetscErrorCode
1153 pcapply_transpose(PC pc, Vec src, Vec dst);
1154
1158 static PetscErrorCode
1159 pcsetup(PC pc);
1160 };
1161} // namespace PETScWrappers
1162
1163
1165
1166
1167#endif // DEAL_II_WITH_PETSC
1168
1169#endif
void Tvmult(VectorBase &dst, const VectorBase &src) const
void vmult(VectorBase &dst, const VectorBase &src) const
void create_pc_with_mat(const MatrixBase &)
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
std::function< void(VectorBase &dst, const VectorBase &src)> vmult
static PetscErrorCode pcapply_transpose(PC pc, Vec src, Vec dst)
static PetscErrorCode pcsetup(PC pc)
void initialize(const MPI_Comm comm)
std::function< void(VectorBase &dst, const VectorBase &src)> vmultT
static PetscErrorCode pcapply(PC pc, Vec src, Vec dst)
Definition point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
*braid_SplitCommworld & comm
AdditionalData(const bool use_vertices=true, const bool use_edges=false, const bool use_faces=false, const bool symmetric=false, const std::vector< Point< dim > > coords={})
AdditionalData(const bool symmetric_operator=false, const double strong_threshold=0.25, const double max_row_sum=0.9, const unsigned int aggressive_coarsening_num_levels=0, const bool output_details=false, const RelaxationType relaxation_type_up=RelaxationType::SORJacobi, const RelaxationType relaxation_type_down=RelaxationType::SORJacobi, const RelaxationType relaxation_type_coarse=RelaxationType::GaussianElimination, const unsigned int n_sweeps_coarse=1, const double tol=0.0, const unsigned int max_iter=1, const bool w_cycle=false)
AdditionalData(const double pivoting=1.e-6, const double zero_pivot=1.e-12, const double damping=0.0)
AdditionalData(const unsigned int symmetric=1, const unsigned int n_levels=1, const double threshold=0.1, const double filter=0.05, const bool output_details=false)