Reference documentation for deal.II version 8.4.2
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 
71  void vmult (VectorBase &dst,
72  const VectorBase &src) const;
73 
74 
78  const PC &get_pc () const;
79 
80  protected:
84  PC pc;
85 
89  Mat matrix;
90 
95  void create_pc ();
96 
102  operator Mat () const;
103 
108  friend class SolverBase;
109  };
110 
111 
112 
125  {
126  public:
132  {};
133 
139 
140 
146  const AdditionalData &additional_data = AdditionalData());
147 
152  PreconditionJacobi (const MPI_Comm communicator,
153  const AdditionalData &additional_data = AdditionalData());
154 
161  void initialize (const MatrixBase &matrix,
162  const AdditionalData &additional_data = AdditionalData());
163 
164  protected:
169 
175  void initialize();
176  };
177 
178 
179 
199  {
200  public:
206  {};
207 
213 
219  const AdditionalData &additional_data = AdditionalData());
220 
225  PreconditionBlockJacobi (const MPI_Comm communicator,
226  const AdditionalData &additional_data = AdditionalData());
227 
228 
235  void initialize (const MatrixBase &matrix,
236  const AdditionalData &additional_data = AdditionalData());
237 
238  protected:
243 
249  void initialize();
250 
251  };
252 
253 
254 
267  {
268  public:
274  {
278  AdditionalData (const double omega = 1);
279 
283  double omega;
284  };
285 
290  PreconditionSOR ();
291 
297  const AdditionalData &additional_data = AdditionalData());
298 
305  void initialize (const MatrixBase &matrix,
306  const AdditionalData &additional_data = AdditionalData());
307 
308  protected:
313  };
314 
315 
316 
329  {
330  public:
336  {
340  AdditionalData (const double omega = 1);
341 
345  double omega;
346  };
347 
352  PreconditionSSOR ();
353 
358  PreconditionSSOR (const MatrixBase &matrix,
359  const AdditionalData &additional_data = AdditionalData());
360 
367  void initialize (const MatrixBase &matrix,
368  const AdditionalData &additional_data = AdditionalData());
369 
370  protected:
375  };
376 
377 
378 
391  {
392  public:
398  {
402  AdditionalData (const double omega = 1);
403 
407  double omega;
408  };
409 
415 
420  PreconditionEisenstat (const MatrixBase &matrix,
421  const AdditionalData &additional_data = AdditionalData());
422 
429  void initialize (const MatrixBase &matrix,
430  const AdditionalData &additional_data = AdditionalData());
431 
432  protected:
437  };
438 
439 
440 
453  {
454  public:
460  {
464  AdditionalData (const unsigned int levels = 0);
465 
469  unsigned int levels;
470  };
471 
476  PreconditionICC ();
477 
482  PreconditionICC (const MatrixBase &matrix,
483  const AdditionalData &additional_data = AdditionalData());
484 
491  void initialize (const MatrixBase &matrix,
492  const AdditionalData &additional_data = AdditionalData());
493 
494  protected:
499  };
500 
501 
502 
515  {
516  public:
522  {
526  AdditionalData (const unsigned int levels = 0);
527 
531  unsigned int levels;
532  };
533 
538  PreconditionILU ();
539 
544  PreconditionILU (const MatrixBase &matrix,
545  const AdditionalData &additional_data = AdditionalData());
546 
553  void initialize (const MatrixBase &matrix,
554  const AdditionalData &additional_data = AdditionalData());
555 
556  protected:
561  };
562 
563 
564 
578  {
579  public:
585  {
590  AdditionalData (const double pivoting = 1.e-6,
591  const double zero_pivot = 1.e-12,
592  const double damping = 0.0);
593 
599  double pivoting;
600 
605  double zero_pivot;
606 
611  double damping;
612  };
613 
618  PreconditionLU ();
619 
624  PreconditionLU (const MatrixBase &matrix,
625  const AdditionalData &additional_data = AdditionalData());
626 
633  void initialize (const MatrixBase &matrix,
634  const AdditionalData &additional_data = AdditionalData());
635 
636  protected:
641  };
642 
643 
644 
645 
658  {
659  public:
665  {
671  const bool symmetric_operator = false,
672  const double strong_threshold = 0.25,
673  const double max_row_sum = 0.9,
674  const unsigned int aggressive_coarsening_num_levels = 0,
675  const bool output_details = false
676  );
677 
685 
692 
704  double max_row_sum;
705 
712 
718  };
719 
725 
730  PreconditionBoomerAMG (const MatrixBase &matrix,
731  const AdditionalData &additional_data = AdditionalData());
732 
737  PreconditionBoomerAMG (const MPI_Comm communicator,
738  const AdditionalData &additional_data = AdditionalData());
739 
740 
747  void initialize (const MatrixBase &matrix,
748  const AdditionalData &additional_data = AdditionalData());
749 
750  protected:
755 
761  void initialize();
762 
763  };
764 
765 
766 
788  {
789  public:
795  {
800  const unsigned int symmetric = 1,
801  const unsigned int n_levels = 1,
802  const double threshold = 0.1,
803  const double filter = 0.05,
804  const bool output_details = false
805  );
806 
818  unsigned int symmetric;
819 
826  unsigned int n_levels;
827 
838  double threshold;
839 
850  double filter;
851 
857  };
858 
859 
860 
866 
871  PreconditionParaSails (const MatrixBase &matrix,
872  const AdditionalData &additional_data = AdditionalData());
873 
880  void initialize (const MatrixBase &matrix,
881  const AdditionalData &additional_data = AdditionalData());
882 
883  private:
888  };
889 
890 
891 
899  {
900  public:
906  {};
907 
912  PreconditionNone ();
913 
919  PreconditionNone (const MatrixBase &matrix,
920  const AdditionalData &additional_data = AdditionalData());
921 
929  void initialize (const MatrixBase &matrix,
930  const AdditionalData &additional_data = AdditionalData());
931 
932  private:
937  };
938 }
939 
940 
941 
942 DEAL_II_NAMESPACE_CLOSE
943 
944 
945 #endif // DEAL_II_WITH_PETSC
946 
947 /*---------------------------- petsc_precondition.h ---------------------------*/
948 
949 #endif
950 /*---------------------------- petsc_precondition.h ---------------------------*/
void vmult(VectorBase &dst, const VectorBase &src) const