Reference documentation for deal.II version 9.0.0
petsc_solver.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_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 
27 # include <petscksp.h>
28 
29 # include <memory>
30 
31 #ifdef DEAL_II_WITH_SLEPC
32 #include <deal.II/lac/slepc_spectral_transformation.h>
33 #endif
34 
35 DEAL_II_NAMESPACE_OPEN
36 
37 #ifdef DEAL_II_WITH_SLEPC
38 namespace SLEPcWrappers
39 {
40  // forward declarations
41  class TransformationBase;
42 }
43 #endif
44 
45 namespace PETScWrappers
46 {
47  // forward declarations
48  class MatrixBase;
49  class VectorBase;
50  class PreconditionerBase;
51 
52 
108  {
109  public:
120  const MPI_Comm &mpi_communicator);
121 
125  virtual ~SolverBase () = default;
126 
134  void
135  solve (const MatrixBase &A,
136  VectorBase &x,
137  const VectorBase &b,
138  const PreconditionerBase &preconditioner);
139 
140 
145  virtual void reset();
146 
147 
152  void set_prefix(const std::string &prefix);
153 
154 
158  SolverControl &control() const;
159 
164  void initialize(const PreconditionerBase &preconditioner);
165 
166  protected:
167 
175 
179  const MPI_Comm mpi_communicator;
180 
185  virtual void set_solver_type (KSP &ksp) const = 0;
186 
193  std::string prefix_name;
194 
195  private:
202  static
203  PetscErrorCode convergence_test (KSP ksp,
204  const PetscInt iteration,
205  const PetscReal residual_norm,
206  KSPConvergedReason *reason,
207  void *solver_control);
208 
226  struct SolverData
227  {
231  ~SolverData ();
232 
236  KSP ksp;
237  };
238 
243  std::unique_ptr<SolverData> solver_data;
244 
245 #ifdef DEAL_II_WITH_SLEPC
246 
251 #endif
252  };
253 
254 
255 
264  {
265  public:
270  {
274  explicit
275  AdditionalData (const double omega = 1);
276 
280  double omega;
281  };
282 
300  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
301  const AdditionalData &data = AdditionalData());
302 
303  protected:
308 
313  virtual void set_solver_type (KSP &ksp) const;
314  };
315 
316 
317 
326  {
327  public:
332  {};
333 
351  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
352  const AdditionalData &data = AdditionalData());
353 
354  protected:
359 
364  virtual void set_solver_type (KSP &ksp) const;
365  };
366 
367 
368 
375  class SolverCG : public SolverBase
376  {
377  public:
382  {};
383 
400  SolverCG (SolverControl &cn,
401  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
402  const AdditionalData &data = AdditionalData());
403 
404  protected:
409 
414  virtual void set_solver_type (KSP &ksp) const;
415  };
416 
417 
418 
425  class SolverBiCG : public SolverBase
426  {
427  public:
432  {};
433 
451  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
452  const AdditionalData &data = AdditionalData());
453 
454  protected:
459 
464  virtual void set_solver_type (KSP &ksp) const;
465  };
466 
467 
468 
475  class SolverGMRES : public SolverBase
476  {
477  public:
482  {
487  AdditionalData (const unsigned int restart_parameter = 30,
488  const bool right_preconditioning = false);
489 
493  unsigned int restart_parameter;
494 
499  };
500 
518  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
519  const AdditionalData &data = AdditionalData());
520 
521  protected:
526 
531  virtual void set_solver_type (KSP &ksp) const;
532  };
533 
534 
535 
543  class SolverBicgstab : public SolverBase
544  {
545  public:
550  {};
551 
569  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
570  const AdditionalData &data = AdditionalData());
571 
572  protected:
577 
582  virtual void set_solver_type (KSP &ksp) const;
583  };
584 
592  class SolverCGS : public SolverBase
593  {
594  public:
599  {};
600 
618  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
619  const AdditionalData &data = AdditionalData());
620 
621  protected:
626 
631  virtual void set_solver_type (KSP &ksp) const;
632  };
633 
634 
635 
642  class SolverTFQMR : public SolverBase
643  {
644  public:
649  {};
650 
668  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
669  const AdditionalData &data = AdditionalData());
670 
671  protected:
676 
681  virtual void set_solver_type (KSP &ksp) const;
682  };
683 
684 
685 
686 
698  class SolverTCQMR : public SolverBase
699  {
700  public:
705  {};
706 
724  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
725  const AdditionalData &data = AdditionalData());
726 
727  protected:
732 
737  virtual void set_solver_type (KSP &ksp) const;
738  };
739 
740 
741 
748  class SolverCR : public SolverBase
749  {
750  public:
755  {};
756 
773  SolverCR (SolverControl &cn,
774  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
775  const AdditionalData &data = AdditionalData());
776 
777  protected:
782 
787  virtual void set_solver_type (KSP &ksp) const;
788  };
789 
790 
791 
799  class SolverLSQR : public SolverBase
800  {
801  public:
806  {};
807 
825  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
826  const AdditionalData &data = AdditionalData());
827 
828  protected:
833 
838  virtual void set_solver_type (KSP &ksp) const;
839  };
840 
841 
854  class SolverPreOnly : public SolverBase
855  {
856  public:
861  {};
862 
880  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
881  const AdditionalData &data = AdditionalData());
882 
883  protected:
888 
893  virtual void set_solver_type (KSP &ksp) const;
894  };
895 
921  {
922  public:
927  {};
932  const MPI_Comm &mpi_communicator = PETSC_COMM_SELF,
933  const AdditionalData &data = AdditionalData());
934 
938  void solve (const MatrixBase &A,
939  VectorBase &x,
940  const VectorBase &b);
941 
947  void set_symmetric_mode (const bool flag);
948 
949  protected:
954 
955  virtual void set_solver_type (KSP &ksp) const;
956 
957  private:
964  static
965  PetscErrorCode convergence_test (KSP ksp,
966  const PetscInt iteration,
967  const PetscReal residual_norm,
968  KSPConvergedReason *reason,
969  void *solver_control);
970 
977  {
981  ~SolverDataMUMPS ();
982 
983  KSP ksp;
984  PC pc;
985  };
986 
987  std::unique_ptr<SolverDataMUMPS> solver_data;
988 
994  };
995 }
996 
997 DEAL_II_NAMESPACE_CLOSE
998 
999 #endif // DEAL_II_WITH_PETSC
1000 
1001 /*---------------------------- petsc_solver.h ---------------------------*/
1002 
1003 #endif
1004 /*---------------------------- petsc_solver.h ---------------------------*/
void solve(const MatrixBase &A, VectorBase &x, const VectorBase &b, const PreconditionerBase &preconditioner)
Definition: petsc_solver.cc:54
virtual void set_solver_type(KSP &ksp) const =0
const AdditionalData additional_data
Definition: petsc_solver.h:408
const AdditionalData additional_data
Definition: petsc_solver.h:832
const AdditionalData additional_data
Definition: petsc_solver.h:307
virtual void set_solver_type(KSP &ksp) const
const AdditionalData additional_data
Definition: petsc_solver.h:358
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:731
virtual void set_solver_type(KSP &ksp) const
const AdditionalData additional_data
Definition: petsc_solver.h:525
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)
const MPI_Comm mpi_communicator
Definition: petsc_solver.h:179
std::unique_ptr< SolverData > solver_data
Definition: petsc_solver.h:243
SolverBase(SolverControl &cn, const MPI_Comm &mpi_communicator)
Definition: petsc_solver.cc:44
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:887
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:174
const AdditionalData additional_data
Definition: petsc_solver.h:458
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:576
void initialize(const PreconditionerBase &preconditioner)
virtual ~SolverBase()=default
SolverBicgstab(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
Definition: petsc_solver.h:953
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:625
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:781
SolverGMRES(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
Definition: petsc_solver.h:675
virtual void set_solver_type(KSP &ksp) const