Reference documentation for deal.II version 8.5.1
petsc_solver.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2017 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_solver_h
17 #define dealii__petsc_solver_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 <deal.II/lac/solver_control.h>
26 # include <deal.II/base/std_cxx11/shared_ptr.h>
27 
28 # include <petscksp.h>
29 
30 #ifdef DEAL_II_WITH_SLEPC
31 #include <deal.II/lac/slepc_spectral_transformation.h>
32 #endif
33 
34 DEAL_II_NAMESPACE_OPEN
35 
36 #ifdef DEAL_II_WITH_SLEPC
37 namespace SLEPcWrappers
38 {
39  // forward declarations
40  class TransformationBase;
41 }
42 #endif
43 
44 namespace PETScWrappers
45 {
46  // forward declarations
47  class MatrixBase;
48  class VectorBase;
49  class PreconditionerBase;
50 
51 
107  {
108  public:
119  const MPI_Comm &mpi_communicator);
120 
124  virtual ~SolverBase ();
125 
133  void
134  solve (const MatrixBase &A,
135  VectorBase &x,
136  const VectorBase &b,
137  const PreconditionerBase &preconditioner);
138 
139 
144  virtual void reset();
145 
146 
151  void set_prefix(const std::string &prefix);
152 
153 
157  SolverControl &control() const;
158 
163  void initialize(const PreconditionerBase &preconditioner);
164 
172  int,
173  << "An error with error number " << arg1
174  << " occurred while calling a PETSc function") DEAL_II_DEPRECATED;
175 
176  protected:
177 
185 
189  const MPI_Comm mpi_communicator;
190 
195  virtual void set_solver_type (KSP &ksp) const = 0;
196 
203  std::string prefix_name;
204 
205  private:
212  static
213  PetscErrorCode convergence_test (KSP ksp,
214  const PetscInt iteration,
215  const PetscReal residual_norm,
216  KSPConvergedReason *reason,
217  void *solver_control);
218 
236  struct SolverData
237  {
241  ~SolverData ();
242 
246  KSP ksp;
247  };
248 
253  std_cxx11::shared_ptr<SolverData> solver_data;
254 
255 #ifdef DEAL_II_WITH_SLEPC
256 
261 #endif
262  };
263 
264 
265 
274  {
275  public:
280  {
284  explicit
285  AdditionalData (const double omega = 1);
286 
290  double omega;
291  };
292 
310  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
311  const AdditionalData &data = AdditionalData());
312 
313  protected:
318 
323  virtual void set_solver_type (KSP &ksp) const;
324  };
325 
326 
327 
336  {
337  public:
342  {};
343 
361  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
362  const AdditionalData &data = AdditionalData());
363 
364  protected:
369 
374  virtual void set_solver_type (KSP &ksp) const;
375  };
376 
377 
378 
385  class SolverCG : public SolverBase
386  {
387  public:
392  {};
393 
410  SolverCG (SolverControl &cn,
411  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
412  const AdditionalData &data = AdditionalData());
413 
414  protected:
419 
424  virtual void set_solver_type (KSP &ksp) const;
425  };
426 
427 
428 
435  class SolverBiCG : public SolverBase
436  {
437  public:
442  {};
443 
461  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
462  const AdditionalData &data = AdditionalData());
463 
464  protected:
469 
474  virtual void set_solver_type (KSP &ksp) const;
475  };
476 
477 
478 
485  class SolverGMRES : public SolverBase
486  {
487  public:
492  {
497  AdditionalData (const unsigned int restart_parameter = 30,
498  const bool right_preconditioning = false);
499 
503  unsigned int restart_parameter;
504 
509  };
510 
528  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
529  const AdditionalData &data = AdditionalData());
530 
531  protected:
536 
541  virtual void set_solver_type (KSP &ksp) const;
542  };
543 
544 
545 
553  class SolverBicgstab : public SolverBase
554  {
555  public:
560  {};
561 
579  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
580  const AdditionalData &data = AdditionalData());
581 
582  protected:
587 
592  virtual void set_solver_type (KSP &ksp) const;
593  };
594 
602  class SolverCGS : public SolverBase
603  {
604  public:
609  {};
610 
628  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
629  const AdditionalData &data = AdditionalData());
630 
631  protected:
636 
641  virtual void set_solver_type (KSP &ksp) const;
642  };
643 
644 
645 
652  class SolverTFQMR : public SolverBase
653  {
654  public:
659  {};
660 
678  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
679  const AdditionalData &data = AdditionalData());
680 
681  protected:
686 
691  virtual void set_solver_type (KSP &ksp) const;
692  };
693 
694 
695 
696 
708  class SolverTCQMR : public SolverBase
709  {
710  public:
715  {};
716 
734  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
735  const AdditionalData &data = AdditionalData());
736 
737  protected:
742 
747  virtual void set_solver_type (KSP &ksp) const;
748  };
749 
750 
751 
758  class SolverCR : public SolverBase
759  {
760  public:
765  {};
766 
783  SolverCR (SolverControl &cn,
784  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
785  const AdditionalData &data = AdditionalData());
786 
787  protected:
792 
797  virtual void set_solver_type (KSP &ksp) const;
798  };
799 
800 
801 
809  class SolverLSQR : public SolverBase
810  {
811  public:
816  {};
817 
835  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
836  const AdditionalData &data = AdditionalData());
837 
838  protected:
843 
848  virtual void set_solver_type (KSP &ksp) const;
849  };
850 
851 
864  class SolverPreOnly : public SolverBase
865  {
866  public:
871  {};
872 
890  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
891  const AdditionalData &data = AdditionalData());
892 
893  protected:
898 
903  virtual void set_solver_type (KSP &ksp) const;
904  };
905 
931  {
932  public:
937  {};
942  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
943  const AdditionalData &data = AdditionalData());
944 
948  void solve (const MatrixBase &A,
949  VectorBase &x,
950  const VectorBase &b);
951 
957  void set_symmetric_mode (const bool flag);
958 
959  protected:
964 
965  virtual void set_solver_type (KSP &ksp) const;
966 
967  private:
974  static
975  PetscErrorCode convergence_test (KSP ksp,
976  const PetscInt iteration,
977  const PetscReal residual_norm,
978  KSPConvergedReason *reason,
979  void *solver_control);
980 
987  {
991  ~SolverDataMUMPS ();
992 
993  KSP ksp;
994  PC pc;
995  };
996 
997  std_cxx11::shared_ptr<SolverDataMUMPS> solver_data;
998 
1004  };
1005 }
1006 
1007 DEAL_II_NAMESPACE_CLOSE
1008 
1009 #endif // DEAL_II_WITH_PETSC
1010 
1011 /*---------------------------- petsc_solver.h ---------------------------*/
1012 
1013 #endif
1014 /*---------------------------- petsc_solver.h ---------------------------*/
void solve(const MatrixBase &A, VectorBase &x, const VectorBase &b, const PreconditionerBase &preconditioner)
Definition: petsc_solver.cc:58
virtual void set_solver_type(KSP &ksp) const =0
const AdditionalData additional_data
Definition: petsc_solver.h:418
const AdditionalData additional_data
Definition: petsc_solver.h:842
const AdditionalData additional_data
Definition: petsc_solver.h:317
virtual void set_solver_type(KSP &ksp) const
const AdditionalData additional_data
Definition: petsc_solver.h:368
AdditionalData(const unsigned int restart_parameter=30, const bool right_preconditioning=false)
SolverRichardson(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
Definition: petsc_solver.h:741
virtual void set_solver_type(KSP &ksp) const
std_cxx11::shared_ptr< SolverData > solver_data
Definition: petsc_solver.h:253
const AdditionalData additional_data
Definition: petsc_solver.h:535
virtual void set_solver_type(KSP &ksp) const
SparseDirectMUMPS(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
void set_prefix(const std::string &prefix)
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:564
const MPI_Comm mpi_communicator
Definition: petsc_solver.h:189
SolverBase(SolverControl &cn, const MPI_Comm &mpi_communicator)
Definition: petsc_solver.cc:43
static PetscErrorCode convergence_test(KSP ksp, const PetscInt iteration, const PetscReal residual_norm, KSPConvergedReason *reason, void *solver_control)
virtual void set_solver_type(KSP &ksp) const
const AdditionalData additional_data
Definition: petsc_solver.h:897
virtual void set_solver_type(KSP &ksp) const
virtual void set_solver_type(KSP &ksp) const
SolverCR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
SolverControl & solver_control
Definition: petsc_solver.h:184
const AdditionalData additional_data
Definition: petsc_solver.h:468
SolverBiCG(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
SolverLSQR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
SolverChebychev(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
SolverCGS(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
virtual void set_solver_type(KSP &ksp) const
const AdditionalData additional_data
Definition: petsc_solver.h:586
void initialize(const PreconditionerBase &preconditioner)
SolverBicgstab(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
Definition: petsc_solver.h:963
void solve(const MatrixBase &A, VectorBase &x, const VectorBase &b)
virtual void set_solver_type(KSP &ksp) const
virtual void set_solver_type(KSP &ksp) const
SolverTCQMR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
SolverControl & control() const
virtual void set_solver_type(KSP &ksp) const
const AdditionalData additional_data
Definition: petsc_solver.h:635
SolverTFQMR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
void set_symmetric_mode(const bool flag)
virtual void set_solver_type(KSP &ksp) const
virtual void set_solver_type(KSP &ksp) const
static PetscErrorCode convergence_test(KSP ksp, const PetscInt iteration, const PetscReal residual_norm, KSPConvergedReason *reason, void *solver_control)
SolverCG(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
SolverPreOnly(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
Definition: petsc_solver.h:791
SolverGMRES(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
Definition: petsc_solver.h:685
virtual void set_solver_type(KSP &ksp) const