Reference documentation for deal.II version 9.0.0
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 at
12 // the top level of the deal.II distribution.
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 # include <petscpc.h>
26 
27 DEAL_II_NAMESPACE_OPEN
28 
29 
30 
31 namespace PETScWrappers
32 {
33  // forward declarations
34  class MatrixBase;
35  class VectorBase;
36  class SolverBase;
37 
38 
56  {
57  public:
62 
66  virtual ~PreconditionerBase ();
67 
72  void clear ();
73 
77  void vmult (VectorBase &dst,
78  const VectorBase &src) const;
79 
80 
84  const PC &get_pc () const;
85 
86  protected:
90  PC pc;
91 
95  Mat matrix;
96 
101  void create_pc ();
102 
108  operator Mat () const;
109 
114  friend class SolverBase;
115  };
116 
117 
118 
131  {
132  public:
138  {};
139 
144  PreconditionJacobi () = default;
145 
146 
153 
158  PreconditionJacobi (const MPI_Comm communicator,
160 
167  void initialize (const MatrixBase &matrix,
169 
170  protected:
175 
181  void initialize();
182  };
183 
184 
185 
211  {
212  public:
218  {};
219 
224  PreconditionBlockJacobi () = default;
225 
232 
237  PreconditionBlockJacobi (const MPI_Comm communicator,
239 
240 
247  void initialize (const MatrixBase &matrix,
249 
250  protected:
255 
261  void initialize();
262 
263  };
264 
265 
266 
277  {
278  public:
284  {
288  AdditionalData (const double omega = 1);
289 
293  double omega;
294  };
295 
300  PreconditionSOR () = default;
301 
308 
315  void initialize (const MatrixBase &matrix,
317 
318  protected:
323  };
324 
325 
326 
337  {
338  public:
344  {
348  AdditionalData (const double omega = 1);
349 
353  double omega;
354  };
355 
360  PreconditionSSOR () = default;
361 
368 
375  void initialize (const MatrixBase &matrix,
377 
378  protected:
383  };
384 
385 
386 
400  {
401  public:
407  {
411  AdditionalData (const double omega = 1);
412 
416  double omega;
417  };
418 
423  PreconditionEisenstat () = default;
424 
431 
438  void initialize (const MatrixBase &matrix,
440 
441  protected:
446  };
447 
448 
449 
460  {
461  public:
467  {
471  AdditionalData (const unsigned int levels = 0);
472 
476  unsigned int levels;
477  };
478 
483  PreconditionICC () = default;
484 
491 
498  void initialize (const MatrixBase &matrix,
500 
501  protected:
506  };
507 
508 
509 
520  {
521  public:
527  {
531  AdditionalData (const unsigned int levels = 0);
532 
536  unsigned int levels;
537  };
538 
543  PreconditionILU () = default;
544 
551 
558  void initialize (const MatrixBase &matrix,
560 
561  protected:
566  };
567 
568 
569 
585  {
586  public:
592  {
597  AdditionalData (const double pivoting = 1.e-6,
598  const double zero_pivot = 1.e-12,
599  const double damping = 0.0);
600 
606  double pivoting;
607 
612  double zero_pivot;
613 
618  double damping;
619  };
620 
625  PreconditionLU () = default;
626 
633 
640  void initialize (const MatrixBase &matrix,
642 
643  protected:
648  };
649 
650 
651 
652 
665  {
666  public:
672  {
678  const bool symmetric_operator = false,
679  const double strong_threshold = 0.25,
680  const double max_row_sum = 0.9,
681  const unsigned int aggressive_coarsening_num_levels = 0,
682  const bool output_details = false
683  );
684 
692 
699 
711  double max_row_sum;
712 
719 
725  };
726 
731  PreconditionBoomerAMG () = default;
732 
739 
744  PreconditionBoomerAMG (const MPI_Comm communicator,
746 
747 
754  void initialize (const MatrixBase &matrix,
756 
757  protected:
762 
768  void initialize();
769 
770  };
771 
772 
773 
795  {
796  public:
802  {
807  const unsigned int symmetric = 1,
808  const unsigned int n_levels = 1,
809  const double threshold = 0.1,
810  const double filter = 0.05,
811  const bool output_details = false
812  );
813 
825  unsigned int symmetric;
826 
833  unsigned int n_levels;
834 
845  double threshold;
846 
857  double filter;
858 
864  };
865 
866 
867 
872  PreconditionParaSails () = default;
873 
880 
887  void initialize (const MatrixBase &matrix,
889 
890  private:
895  };
896 
897 
898 
906  {
907  public:
913  {};
914 
919  PreconditionNone () = default;
920 
928 
936  void initialize (const MatrixBase &matrix,
938 
939  private:
944  };
945 }
946 
947 
948 
949 DEAL_II_NAMESPACE_CLOSE
950 
951 
952 #endif // DEAL_II_WITH_PETSC
953 
954 /*---------------------------- petsc_precondition.h ---------------------------*/
955 
956 #endif
957 /*---------------------------- 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)