Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
petsc_precondition.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2018 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 
22 # ifdef DEAL_II_WITH_PETSC
23 
24 # include <deal.II/lac/exceptions.h>
25 
26 # include <petscpc.h>
27 
28 DEAL_II_NAMESPACE_OPEN
29 
30 
31 
32 namespace PETScWrappers
33 {
34  // forward declarations
35  class MatrixBase;
36  class VectorBase;
37  class SolverBase;
38 
39 
57  {
58  public:
63 
67  virtual ~PreconditionerBase();
68 
73  void
74  clear();
75 
79  void
80  vmult(VectorBase &dst, const VectorBase &src) const;
81 
82 
86  const PC &
87  get_pc() const;
88 
89  protected:
93  PC pc;
94 
98  Mat matrix;
99 
104  void
105  create_pc();
106 
112  operator Mat() const;
113 
118  friend class SolverBase;
119  };
120 
121 
122 
135  {
136  public:
142  {};
143 
148  PreconditionJacobi() = default;
149 
150 
156  const MatrixBase & matrix,
158 
164  const MPI_Comm communicator,
166 
173  void
174  initialize(const MatrixBase & matrix,
176 
177  protected:
182 
188  void
189  initialize();
190  };
191 
192 
193 
219  {
220  public:
226  {};
227 
232  PreconditionBlockJacobi() = default;
233 
239  const MatrixBase & matrix,
241 
247  const MPI_Comm communicator,
249 
250 
257  void
258  initialize(const MatrixBase & matrix,
260 
261  protected:
266 
272  void
273  initialize();
274  };
275 
276 
277 
288  {
289  public:
295  {
299  AdditionalData(const double omega = 1);
300 
304  double omega;
305  };
306 
311  PreconditionSOR() = default;
312 
319 
326  void
327  initialize(const MatrixBase & matrix,
329 
330  protected:
335  };
336 
337 
338 
349  {
350  public:
356  {
360  AdditionalData(const double omega = 1);
361 
365  double omega;
366  };
367 
372  PreconditionSSOR() = default;
373 
380 
387  void
388  initialize(const MatrixBase & matrix,
390 
391  protected:
396  };
397 
398 
399 
413  {
414  public:
420  {
424  AdditionalData(const double omega = 1);
425 
429  double omega;
430  };
431 
436  PreconditionEisenstat() = default;
437 
443  const MatrixBase & matrix,
445 
452  void
453  initialize(const MatrixBase & matrix,
455 
456  protected:
461  };
462 
463 
464 
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 
536  {
537  public:
543  {
547  AdditionalData(const unsigned int levels = 0);
548 
552  unsigned int levels;
553  };
554 
559  PreconditionILU() = default;
560 
567 
574  void
575  initialize(const MatrixBase & matrix,
577 
578  protected:
583  };
584 
585 
586 
602  {
603  public:
609  {
614  AdditionalData(const double pivoting = 1.e-6,
615  const double zero_pivot = 1.e-12,
616  const double damping = 0.0);
617 
623  double pivoting;
624 
629  double zero_pivot;
630 
635  double damping;
636  };
637 
642  PreconditionLU() = default;
643 
650 
657  void
658  initialize(const MatrixBase & matrix,
660 
661  protected:
666  };
667 
668 
669 
682  {
683  public:
689  {
694  AdditionalData(const bool symmetric_operator = false,
695  const double strong_threshold = 0.25,
696  const double max_row_sum = 0.9,
697  const unsigned int aggressive_coarsening_num_levels = 0,
698  const bool output_details = false);
699 
707 
714 
726  double max_row_sum;
727 
734 
740  };
741 
746  PreconditionBoomerAMG() = default;
747 
753  const MatrixBase & matrix,
755 
761  const MPI_Comm communicator,
763 
764 
771  void
772  initialize(const MatrixBase & matrix,
774 
775  protected:
780 
786  void
787  initialize();
788  };
789 
790 
791 
813  {
814  public:
820  {
824  AdditionalData(const unsigned int symmetric = 1,
825  const unsigned int n_levels = 1,
826  const double threshold = 0.1,
827  const double filter = 0.05,
828  const bool output_details = false);
829 
841  unsigned int symmetric;
842 
849  unsigned int n_levels;
850 
861  double threshold;
862 
873  double filter;
874 
880  };
881 
882 
883 
888  PreconditionParaSails() = default;
889 
895  const MatrixBase & matrix,
897 
904  void
905  initialize(const MatrixBase & matrix,
907 
908  private:
913  };
914 
915 
916 
924  {
925  public:
931  {};
932 
937  PreconditionNone() = default;
938 
946 
954  void
955  initialize(const MatrixBase & matrix,
957 
958  private:
963  };
964 } // namespace PETScWrappers
965 
966 
967 
968 DEAL_II_NAMESPACE_CLOSE
969 
970 
971 # endif // DEAL_II_WITH_PETSC
972 
973 #endif
974 /*--------------------------- petsc_precondition.h --------------------------*/
void vmult(VectorBase &dst, const VectorBase &src) const
AdditionalData(const double pivoting=1.e-6, const double zero_pivot=1.e-12, const double damping=0.0)
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())
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)
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
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)