Reference documentation for deal.II version 9.3.2
\(\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\}}\)
petsc_precondition.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2021 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_petsc_precondition_h
17 # define dealii_petsc_precondition_h
18 
19 
20 # include <deal.II/base/config.h>
21 
23 
24 # ifdef DEAL_II_WITH_PETSC
25 
26 # include <deal.II/lac/exceptions.h>
27 
28 # include <petscpc.h>
29 
31 
32 
33 
34 namespace PETScWrappers
35 {
36  // forward declarations
37 # ifndef DOXYGEN
38  class MatrixBase;
39  class VectorBase;
40  class SolverBase;
41 # endif
42 
59  {
60  public:
65 
69  virtual ~PreconditionBase();
70 
75  void
76  clear();
77 
81  void
82  vmult(VectorBase &dst, const VectorBase &src) const;
83 
87  void
88  Tvmult(VectorBase &dst, const VectorBase &src) const;
89 
90 
94  const PC &
95  get_pc() const;
96 
97  protected:
101  PC pc;
102 
106  Mat matrix;
107 
112  void
113  create_pc();
114 
120  operator Mat() const;
121 
122  // Make the solver class a friend, since it needs to call the conversion
123  // operator.
124  friend class SolverBase;
125  };
126 
127 
128 
140  {
141  public:
147  {};
148 
153  PreconditionJacobi() = default;
154 
155 
161  const MatrixBase & matrix,
163 
169  const MPI_Comm & communicator,
171 
178  void
179  initialize(const MatrixBase & matrix,
181 
182  protected:
187 
193  void
194  initialize();
195  };
196 
197 
198 
223  {
224  public:
230  {};
231 
237 
243  const MatrixBase & matrix,
245 
251  const MPI_Comm & communicator,
253 
254 
261  void
262  initialize(const MatrixBase & matrix,
264 
265  protected:
270 
276  void
277  initialize();
278  };
279 
280 
281 
291  {
292  public:
298  {
302  AdditionalData(const double omega = 1);
303 
307  double omega;
308  };
309 
314  PreconditionSOR() = default;
315 
322 
329  void
330  initialize(const MatrixBase & matrix,
332 
333  protected:
338  };
339 
340 
341 
351  {
352  public:
358  {
362  AdditionalData(const double omega = 1);
363 
367  double omega;
368  };
369 
374  PreconditionSSOR() = default;
375 
382 
389  void
390  initialize(const MatrixBase & matrix,
392 
393  protected:
398  };
399 
400 
401 
414  {
415  public:
421  {
425  AdditionalData(const double omega = 1);
426 
430  double omega;
431  };
432 
438 
444  const MatrixBase & matrix,
446 
453  void
454  initialize(const MatrixBase & matrix,
456 
457  protected:
462  };
463 
464 
465 
475  {
476  public:
482  {
486  AdditionalData(const unsigned int levels = 0);
487 
491  unsigned int levels;
492  };
493 
498  PreconditionICC() = default;
499 
506 
513  void
514  initialize(const MatrixBase & matrix,
516 
517  protected:
522  };
523 
524 
525 
535  {
536  public:
542  {
546  AdditionalData(const unsigned int levels = 0);
547 
551  unsigned int levels;
552  };
553 
558  PreconditionILU() = default;
559 
566 
573  void
574  initialize(const MatrixBase & matrix,
576 
577  protected:
582  };
583 
584 
585 
600  {
601  public:
607  {
612  AdditionalData(const double pivoting = 1.e-6,
613  const double zero_pivot = 1.e-12,
614  const double damping = 0.0);
615 
621  double pivoting;
622 
627  double zero_pivot;
628 
633  double damping;
634  };
635 
640  PreconditionLU() = default;
641 
648 
655  void
656  initialize(const MatrixBase & matrix,
658 
659  protected:
664  };
665 
666 
667 
679  {
680  public:
686  {
690  enum class RelaxationType
691  {
692  Jacobi,
695  SORJacobi,
702  CG,
703  Chebyshev,
704  FCFJacobi,
706  None
707  };
708 
714  const bool symmetric_operator = false,
715  const double strong_threshold = 0.25,
716  const double max_row_sum = 0.9,
717  const unsigned int aggressive_coarsening_num_levels = 0,
718  const bool output_details = false,
723  const unsigned int n_sweeps_coarse = 1,
724  const double tol = 0.0,
725  const unsigned int max_iter = 1,
726  const bool w_cycle = false);
727 
734 
741 
753  double max_row_sum;
754 
761 
767 
772 
777 
782 
786  unsigned int n_sweeps_coarse;
787 
791  double tol;
792 
796  unsigned int max_iter;
797 
802  bool w_cycle;
803  };
804 
810 
816  const MatrixBase & matrix,
818 
824  const MPI_Comm & communicator,
826 
827 
834  void
835  initialize(const MatrixBase & matrix,
837 
838  protected:
843 
849  void
850  initialize();
851  };
852 
853 
854 
875  {
876  public:
882  {
886  AdditionalData(const unsigned int symmetric = 1,
887  const unsigned int n_levels = 1,
888  const double threshold = 0.1,
889  const double filter = 0.05,
890  const bool output_details = false);
891 
903  unsigned int symmetric;
904 
911  unsigned int n_levels;
912 
923  double threshold;
924 
935  double filter;
936 
942  };
943 
944 
945 
951 
957  const MatrixBase & matrix,
959 
966  void
967  initialize(const MatrixBase & matrix,
969 
970  private:
975  };
976 
977 
978 
985  {
986  public:
992  {};
993 
998  PreconditionNone() = default;
999 
1007 
1015  void
1016  initialize(const MatrixBase & matrix,
1018 
1019  private:
1024  };
1025 
1031 } // namespace PETScWrappers
1032 
1033 
1034 
1036 
1037 
1038 # endif // DEAL_II_WITH_PETSC
1039 
1040 #endif
1041 /*--------------------------- petsc_precondition.h --------------------------*/
void Tvmult(VectorBase &dst, const VectorBase &src) const
void vmult(VectorBase &dst, const VectorBase &src) const
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())
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:396
#define DEAL_II_DEPRECATED_EARLY
Definition: config.h:167
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:397
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)