Reference documentation for deal.II version GIT b8135fa6eb 2023-03-25 15:55:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
petsc_solver.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2021 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.md at
12 // the top level directory of deal.II.
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>
26 
27 # include <petscksp.h>
28 
29 # include <memory>
30 
31 # ifdef DEAL_II_WITH_SLEPC
33 # endif
34 
36 
37 // Forward declarations
38 # ifndef DOXYGEN
39 # ifdef DEAL_II_WITH_SLEPC
40 namespace SLEPcWrappers
41 {
42  // forward declarations
43  class TransformationBase;
44 } // namespace SLEPcWrappers
45 # endif
46 # endif
47 
48 namespace PETScWrappers
49 {
50  // forward declarations
51 # ifndef DOXYGEN
52  class MatrixBase;
53  class VectorBase;
54  class PreconditionBase;
55 # endif
56 
57 
88  class SolverBase
89  {
90  public:
95 
99  virtual ~SolverBase() = default;
100 
108  void
109  solve(const MatrixBase & A,
110  VectorBase & x,
111  const VectorBase & b,
112  const PreconditionBase &preconditioner);
113 
114 
119  virtual void
120  reset();
121 
122 
127  void
128  set_prefix(const std::string &prefix);
129 
130 
134  SolverControl &
135  control() const;
136 
141  void
142  initialize(const PreconditionBase &preconditioner);
143 
144  protected:
152 
157  virtual void
158  set_solver_type(KSP &ksp) const = 0;
159 
166  std::string prefix_name;
167 
168  private:
175  static PetscErrorCode
176  convergence_test(KSP ksp,
177  const PetscInt iteration,
178  const PetscReal residual_norm,
179  KSPConvergedReason *reason,
180  void * solver_control);
181 
199  struct SolverData
200  {
204  ~SolverData();
205 
209  KSP ksp;
210  };
211 
216  std::unique_ptr<SolverData> solver_data;
217 
218 # ifdef DEAL_II_WITH_SLEPC
219  // Make the transformation class a friend, since it needs to set the KSP
220  // solver.
222 # endif
223  };
224 
225 
226 
234  {
235  public:
240  {
244  explicit AdditionalData(const double omega = 1);
245 
249  double omega;
250  };
251 
260  const AdditionalData &data = AdditionalData());
261 
268  DEAL_II_DEPRECATED_EARLY
270  const MPI_Comm & mpi_communicator,
271  const AdditionalData &data = AdditionalData());
272 
273  protected:
278 
283  virtual void
284  set_solver_type(KSP &ksp) const override;
285  };
286 
287 
288 
296  {
297  public:
302  {};
303 
312  const AdditionalData &data = AdditionalData());
313 
320  DEAL_II_DEPRECATED_EARLY
322  const MPI_Comm & mpi_communicator,
323  const AdditionalData &data = AdditionalData());
324 
325  protected:
330 
335  virtual void
336  set_solver_type(KSP &ksp) const override;
337  };
338 
339 
340 
346  class SolverCG : public SolverBase
347  {
348  public:
353  {};
354 
362  SolverCG(SolverControl &cn, const AdditionalData &data = AdditionalData());
363 
370  DEAL_II_DEPRECATED_EARLY
371  SolverCG(SolverControl & cn,
372  const MPI_Comm & mpi_communicator,
373  const AdditionalData &data = AdditionalData());
374 
375  protected:
380 
385  virtual void
386  set_solver_type(KSP &ksp) const override;
387  };
388 
389 
390 
396  class SolverBiCG : public SolverBase
397  {
398  public:
403  {};
404 
413  const AdditionalData &data = AdditionalData());
414 
421  DEAL_II_DEPRECATED_EARLY
423  const MPI_Comm & mpi_communicator,
424  const AdditionalData &data = AdditionalData());
425 
426  protected:
431 
436  virtual void
437  set_solver_type(KSP &ksp) const override;
438  };
439 
440 
441 
447  class SolverGMRES : public SolverBase
448  {
449  public:
454  {
459  AdditionalData(const unsigned int restart_parameter = 30,
460  const bool right_preconditioning = false);
461 
465  unsigned int restart_parameter;
466 
471  };
472 
481  const AdditionalData &data = AdditionalData());
482 
489  DEAL_II_DEPRECATED_EARLY
491  const MPI_Comm & mpi_communicator,
492  const AdditionalData &data = AdditionalData());
493 
494  protected:
499 
504  virtual void
505  set_solver_type(KSP &ksp) const override;
506  };
507 
508 
509 
516  class SolverBicgstab : public SolverBase
517  {
518  public:
523  {};
524 
533  const AdditionalData &data = AdditionalData());
534 
541  DEAL_II_DEPRECATED_EARLY
543  const MPI_Comm & mpi_communicator,
544  const AdditionalData &data = AdditionalData());
545 
546  protected:
551 
556  virtual void
557  set_solver_type(KSP &ksp) const override;
558  };
559 
560 
561 
568  class SolverCGS : public SolverBase
569  {
570  public:
575  {};
576 
585 
592  DEAL_II_DEPRECATED_EARLY
594  const MPI_Comm & mpi_communicator,
595  const AdditionalData &data = AdditionalData());
596 
597  protected:
602 
607  virtual void
608  set_solver_type(KSP &ksp) const override;
609  };
610 
611 
612 
618  class SolverTFQMR : public SolverBase
619  {
620  public:
625  {};
626 
635  const AdditionalData &data = AdditionalData());
636 
643  DEAL_II_DEPRECATED_EARLY
645  const MPI_Comm & mpi_communicator,
646  const AdditionalData &data = AdditionalData());
647 
648  protected:
653 
658  virtual void
659  set_solver_type(KSP &ksp) const override;
660  };
661 
662 
663 
674  class SolverTCQMR : public SolverBase
675  {
676  public:
681  {};
682 
691  const AdditionalData &data = AdditionalData());
692 
699  DEAL_II_DEPRECATED_EARLY
701  const MPI_Comm & mpi_communicator,
702  const AdditionalData &data = AdditionalData());
703 
704  protected:
709 
714  virtual void
715  set_solver_type(KSP &ksp) const override;
716  };
717 
718 
719 
725  class SolverCR : public SolverBase
726  {
727  public:
732  {};
733 
741  SolverCR(SolverControl &cn, const AdditionalData &data = AdditionalData());
742 
749  DEAL_II_DEPRECATED_EARLY
750  SolverCR(SolverControl & cn,
751  const MPI_Comm & mpi_communicator,
752  const AdditionalData &data = AdditionalData());
753 
754  protected:
759 
764  virtual void
765  set_solver_type(KSP &ksp) const override;
766  };
767 
768 
769 
776  class SolverLSQR : public SolverBase
777  {
778  public:
783  {};
784 
793  const AdditionalData &data = AdditionalData());
794 
801  DEAL_II_DEPRECATED_EARLY
803  const MPI_Comm & mpi_communicator,
804  const AdditionalData &data = AdditionalData());
805 
806  protected:
811 
816  virtual void
817  set_solver_type(KSP &ksp) const override;
818  };
819 
820 
832  class SolverPreOnly : public SolverBase
833  {
834  public:
839  {};
840 
849  const AdditionalData &data = AdditionalData());
850 
857  DEAL_II_DEPRECATED_EARLY
859  const MPI_Comm & mpi_communicator,
860  const AdditionalData &data = AdditionalData());
861 
862  protected:
867 
872  virtual void
873  set_solver_type(KSP &ksp) const override;
874  };
875 
900  {
901  public:
906  {};
907 
912  const AdditionalData &data = AdditionalData());
913 
920  DEAL_II_DEPRECATED_EARLY
922  const MPI_Comm & mpi_communicator,
923  const AdditionalData &data = AdditionalData());
924 
928  void
929  solve(const MatrixBase &A, VectorBase &x, const VectorBase &b);
930 
936  void
937  set_symmetric_mode(const bool flag);
938 
939  protected:
944 
945  virtual void
946  set_solver_type(KSP &ksp) const override;
947 
948  private:
955  static PetscErrorCode
956  convergence_test(KSP ksp,
957  const PetscInt iteration,
958  const PetscReal residual_norm,
959  KSPConvergedReason *reason,
960  void * solver_control);
961 
968  {
973 
974  KSP ksp;
975  PC pc;
976  };
977 
978  std::unique_ptr<SolverDataMUMPS> solver_data;
979 
985  };
986 } // namespace PETScWrappers
987 
989 
990 #endif // DEAL_II_WITH_PETSC
991 
992 #endif
virtual void set_solver_type(KSP &ksp) const =0
SolverControl & control() const
SolverControl & solver_control
Definition: petsc_solver.h:151
void initialize(const PreconditionBase &preconditioner)
SolverBase(SolverControl &cn)
Definition: petsc_solver.cc:45
static PetscErrorCode convergence_test(KSP ksp, const PetscInt iteration, const PetscReal residual_norm, KSPConvergedReason *reason, void *solver_control)
virtual ~SolverBase()=default
void set_prefix(const std::string &prefix)
void solve(const MatrixBase &A, VectorBase &x, const VectorBase &b, const PreconditionBase &preconditioner)
Definition: petsc_solver.cc:52
std::unique_ptr< SolverData > solver_data
Definition: petsc_solver.h:216
SolverBiCG(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
Definition: petsc_solver.h:430
virtual void set_solver_type(KSP &ksp) const override
virtual void set_solver_type(KSP &ksp) const override
const AdditionalData additional_data
Definition: petsc_solver.h:550
SolverBicgstab(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
Definition: petsc_solver.h:601
virtual void set_solver_type(KSP &ksp) const override
SolverCGS(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
Definition: petsc_solver.h:379
SolverCG(SolverControl &cn, const AdditionalData &data=AdditionalData())
virtual void set_solver_type(KSP &ksp) const override
virtual void set_solver_type(KSP &ksp) const override
const AdditionalData additional_data
Definition: petsc_solver.h:758
SolverCR(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
Definition: petsc_solver.h:329
virtual void set_solver_type(KSP &ksp) const override
SolverChebychev(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
Definition: petsc_solver.h:498
virtual void set_solver_type(KSP &ksp) const override
SolverGMRES(SolverControl &cn, const AdditionalData &data=AdditionalData())
virtual void set_solver_type(KSP &ksp) const override
SolverLSQR(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
Definition: petsc_solver.h:810
SolverPreOnly(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
Definition: petsc_solver.h:866
virtual void set_solver_type(KSP &ksp) const override
virtual void set_solver_type(KSP &ksp) const override
const AdditionalData additional_data
Definition: petsc_solver.h:277
SolverRichardson(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
Definition: petsc_solver.h:708
virtual void set_solver_type(KSP &ksp) const override
SolverTCQMR(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
Definition: petsc_solver.h:652
virtual void set_solver_type(KSP &ksp) const override
SolverTFQMR(SolverControl &cn, const AdditionalData &data=AdditionalData())
void solve(const MatrixBase &A, VectorBase &x, const VectorBase &b)
virtual void set_solver_type(KSP &ksp) const override
const AdditionalData additional_data
Definition: petsc_solver.h:943
void set_symmetric_mode(const bool flag)
static PetscErrorCode convergence_test(KSP ksp, const PetscInt iteration, const PetscReal residual_norm, KSPConvergedReason *reason, void *solver_control)
SparseDirectMUMPS(SolverControl &cn, const AdditionalData &data=AdditionalData())
std::unique_ptr< SolverDataMUMPS > solver_data
Definition: petsc_solver.h:978
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:474
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:475
static const char A
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
AdditionalData(const unsigned int restart_parameter=30, const bool right_preconditioning=false)