Reference documentation for deal.II version 8.5.1
petsc_precondition.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2016 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 
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 
205  {
206  public:
212  {};
213 
219 
226 
231  PreconditionBlockJacobi (const MPI_Comm communicator,
233 
234 
241  void initialize (const MatrixBase &matrix,
243 
244  protected:
249 
255  void initialize();
256 
257  };
258 
259 
260 
271  {
272  public:
278  {
282  AdditionalData (const double omega = 1);
283 
287  double omega;
288  };
289 
294  PreconditionSOR ();
295 
302 
309  void initialize (const MatrixBase &matrix,
311 
312  protected:
317  };
318 
319 
320 
331  {
332  public:
338  {
342  AdditionalData (const double omega = 1);
343 
347  double omega;
348  };
349 
354  PreconditionSSOR ();
355 
362 
369  void initialize (const MatrixBase &matrix,
371 
372  protected:
377  };
378 
379 
380 
394  {
395  public:
401  {
405  AdditionalData (const double omega = 1);
406 
410  double omega;
411  };
412 
418 
425 
432  void initialize (const MatrixBase &matrix,
434 
435  protected:
440  };
441 
442 
443 
454  {
455  public:
461  {
465  AdditionalData (const unsigned int levels = 0);
466 
470  unsigned int levels;
471  };
472 
477  PreconditionICC ();
478 
485 
492  void initialize (const MatrixBase &matrix,
494 
495  protected:
500  };
501 
502 
503 
514  {
515  public:
521  {
525  AdditionalData (const unsigned int levels = 0);
526 
530  unsigned int levels;
531  };
532 
537  PreconditionILU ();
538 
545 
552  void initialize (const MatrixBase &matrix,
554 
555  protected:
560  };
561 
562 
563 
575  {
576  public:
582  {
587  AdditionalData (const double pivoting = 1.e-6,
588  const double zero_pivot = 1.e-12,
589  const double damping = 0.0);
590 
596  double pivoting;
597 
602  double zero_pivot;
603 
608  double damping;
609  };
610 
615  PreconditionLU ();
616 
623 
630  void initialize (const MatrixBase &matrix,
632 
633  protected:
638  };
639 
640 
641 
642 
655  {
656  public:
662  {
668  const bool symmetric_operator = false,
669  const double strong_threshold = 0.25,
670  const double max_row_sum = 0.9,
671  const unsigned int aggressive_coarsening_num_levels = 0,
672  const bool output_details = false
673  );
674 
682 
689 
701  double max_row_sum;
702 
709 
715  };
716 
722 
729 
734  PreconditionBoomerAMG (const MPI_Comm communicator,
736 
737 
744  void initialize (const MatrixBase &matrix,
746 
747  protected:
752 
758  void initialize();
759 
760  };
761 
762 
763 
785  {
786  public:
792  {
797  const unsigned int symmetric = 1,
798  const unsigned int n_levels = 1,
799  const double threshold = 0.1,
800  const double filter = 0.05,
801  const bool output_details = false
802  );
803 
815  unsigned int symmetric;
816 
823  unsigned int n_levels;
824 
835  double threshold;
836 
847  double filter;
848 
854  };
855 
856 
857 
863 
870 
877  void initialize (const MatrixBase &matrix,
879 
880  private:
885  };
886 
887 
888 
896  {
897  public:
903  {};
904 
909  PreconditionNone ();
910 
918 
926  void initialize (const MatrixBase &matrix,
928 
929  private:
934  };
935 }
936 
937 
938 
939 DEAL_II_NAMESPACE_CLOSE
940 
941 
942 #endif // DEAL_II_WITH_PETSC
943 
944 /*---------------------------- petsc_precondition.h ---------------------------*/
945 
946 #endif
947 /*---------------------------- 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)