Reference documentation for deal.II version 9.2.0
\(\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 - 2019 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 PreconditionerBase;
55 # endif
56 
57 
114  {
115  public:
126 
130  virtual ~SolverBase() = default;
131 
139  void
140  solve(const MatrixBase & A,
141  VectorBase & x,
142  const VectorBase & b,
143  const PreconditionerBase &preconditioner);
144 
145 
150  virtual void
151  reset();
152 
153 
158  void
159  set_prefix(const std::string &prefix);
160 
161 
165  SolverControl &
166  control() const;
167 
172  void
173  initialize(const PreconditionerBase &preconditioner);
174 
175  protected:
183 
188 
193  virtual void
194  set_solver_type(KSP &ksp) const = 0;
195 
202  std::string prefix_name;
203 
204  private:
211  static PetscErrorCode
212  convergence_test(KSP ksp,
213  const PetscInt iteration,
214  const PetscReal residual_norm,
215  KSPConvergedReason *reason,
216  void * solver_control);
217 
235  struct SolverData
236  {
240  ~SolverData();
241 
245  KSP ksp;
246  };
247 
252  std::unique_ptr<SolverData> solver_data;
253 
254 # ifdef DEAL_II_WITH_SLEPC
255  // Make the transformation class a friend, since it needs to set the KSP
256  // solver.
258 # endif
259  };
260 
261 
262 
271  {
272  public:
277  {
281  explicit AdditionalData(const double omega = 1);
282 
286  double omega;
287  };
288 
306  const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
307  const AdditionalData &data = AdditionalData());
308 
309  protected:
314 
319  virtual void
320  set_solver_type(KSP &ksp) const override;
321  };
322 
323 
324 
333  {
334  public:
339  {};
340 
358  const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
359  const AdditionalData &data = AdditionalData());
360 
361  protected:
366 
371  virtual void
372  set_solver_type(KSP &ksp) const override;
373  };
374 
375 
376 
383  class SolverCG : public SolverBase
384  {
385  public:
390  {};
391 
408  SolverCG(SolverControl & cn,
409  const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
410  const AdditionalData &data = AdditionalData());
411 
412  protected:
417 
422  virtual void
423  set_solver_type(KSP &ksp) const override;
424  };
425 
426 
427 
434  class SolverBiCG : public SolverBase
435  {
436  public:
441  {};
442 
460  const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
461  const AdditionalData &data = AdditionalData());
462 
463  protected:
468 
473  virtual void
474  set_solver_type(KSP &ksp) const override;
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
542  set_solver_type(KSP &ksp) const override;
543  };
544 
545 
546 
554  class SolverBicgstab : public SolverBase
555  {
556  public:
561  {};
562 
580  const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
581  const AdditionalData &data = AdditionalData());
582 
583  protected:
588 
593  virtual void
594  set_solver_type(KSP &ksp) const override;
595  };
596 
604  class SolverCGS : public SolverBase
605  {
606  public:
611  {};
612 
630  const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
631  const AdditionalData &data = AdditionalData());
632 
633  protected:
638 
643  virtual void
644  set_solver_type(KSP &ksp) const override;
645  };
646 
647 
648 
655  class SolverTFQMR : public SolverBase
656  {
657  public:
662  {};
663 
681  const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
682  const AdditionalData &data = AdditionalData());
683 
684  protected:
689 
694  virtual void
695  set_solver_type(KSP &ksp) const override;
696  };
697 
698 
699 
711  class SolverTCQMR : public SolverBase
712  {
713  public:
718  {};
719 
737  const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
738  const AdditionalData &data = AdditionalData());
739 
740  protected:
745 
750  virtual void
751  set_solver_type(KSP &ksp) const override;
752  };
753 
754 
755 
762  class SolverCR : public SolverBase
763  {
764  public:
769  {};
770 
787  SolverCR(SolverControl & cn,
788  const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
789  const AdditionalData &data = AdditionalData());
790 
791  protected:
796 
801  virtual void
802  set_solver_type(KSP &ksp) const override;
803  };
804 
805 
806 
814  class SolverLSQR : public SolverBase
815  {
816  public:
821  {};
822 
840  const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
841  const AdditionalData &data = AdditionalData());
842 
843  protected:
848 
853  virtual void
854  set_solver_type(KSP &ksp) const override;
855  };
856 
857 
870  class SolverPreOnly : public SolverBase
871  {
872  public:
877  {};
878 
896  const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
897  const AdditionalData &data = AdditionalData());
898 
899  protected:
904 
909  virtual void
910  set_solver_type(KSP &ksp) const override;
911  };
912 
938  {
939  public:
944  {};
949  const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
950  const AdditionalData &data = AdditionalData());
951 
955  void
956  solve(const MatrixBase &A, VectorBase &x, const VectorBase &b);
957 
963  void
964  set_symmetric_mode(const bool flag);
965 
966  protected:
971 
972  virtual void
973  set_solver_type(KSP &ksp) const override;
974 
975  private:
982  static PetscErrorCode
983  convergence_test(KSP ksp,
984  const PetscInt iteration,
985  const PetscReal residual_norm,
986  KSPConvergedReason *reason,
987  void * solver_control);
988 
995  {
1000 
1001  KSP ksp;
1002  PC pc;
1003  };
1004 
1005  std::unique_ptr<SolverDataMUMPS> solver_data;
1006 
1012  };
1013 } // namespace PETScWrappers
1014 
1016 
1017 # endif // DEAL_II_WITH_PETSC
1018 
1019 /*---------------------------- petsc_solver.h ---------------------------*/
1020 
1021 #endif
1022 /*---------------------------- petsc_solver.h ---------------------------*/
PETScWrappers::SolverGMRES::additional_data
const AdditionalData additional_data
Definition: petsc_solver.h:535
PETScWrappers::SolverCG::AdditionalData
Definition: petsc_solver.h:389
PETScWrappers::SolverBicgstab::SolverBicgstab
SolverBicgstab(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
Definition: petsc_solver.cc:448
PETScWrappers::SolverCGS::set_solver_type
virtual void set_solver_type(KSP &ksp) const override
Definition: petsc_solver.cc:481
PETScWrappers::SolverCG
Definition: petsc_solver.h:383
PETScWrappers::SolverGMRES::set_solver_type
virtual void set_solver_type(KSP &ksp) const override
Definition: petsc_solver.cc:396
PETScWrappers::SolverTCQMR::set_solver_type
virtual void set_solver_type(KSP &ksp) const override
Definition: petsc_solver.cc:529
PETScWrappers::SolverTCQMR::additional_data
const AdditionalData additional_data
Definition: petsc_solver.h:744
PETScWrappers::SolverChebychev::additional_data
const AdditionalData additional_data
Definition: petsc_solver.h:365
PETScWrappers::SolverBase
Definition: petsc_solver.h:113
PETScWrappers
Definition: petsc_block_sparse_matrix.h:36
PETScWrappers::SolverBase::~SolverBase
virtual ~SolverBase()=default
PETScWrappers::SolverPreOnly::SolverPreOnly
SolverPreOnly(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
Definition: petsc_solver.cc:592
PETScWrappers::SolverGMRES::AdditionalData::restart_parameter
unsigned int restart_parameter
Definition: petsc_solver.h:503
PETScWrappers::SparseDirectMUMPS::solver_data
std::unique_ptr< SolverDataMUMPS > solver_data
Definition: petsc_solver.h:1005
PETScWrappers::SolverCR::set_solver_type
virtual void set_solver_type(KSP &ksp) const override
Definition: petsc_solver.cc:553
PETScWrappers::SolverCG::additional_data
const AdditionalData additional_data
Definition: petsc_solver.h:416
PETScWrappers::SolverBase::set_solver_type
virtual void set_solver_type(KSP &ksp) const =0
PETScWrappers::SparseDirectMUMPS::solve
void solve(const MatrixBase &A, VectorBase &x, const VectorBase &b)
Definition: petsc_solver.cc:671
PETScWrappers::SolverBase::SolverData::~SolverData
~SolverData()
Definition: petsc_solver.cc:38
PETScWrappers::SolverBiCG
Definition: petsc_solver.h:434
PETScWrappers::SolverTCQMR::SolverTCQMR
SolverTCQMR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
Definition: petsc_solver.cc:520
PETScWrappers::SolverChebychev::SolverChebychev
SolverChebychev(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
Definition: petsc_solver.cc:306
PETScWrappers::SolverBicgstab::additional_data
const AdditionalData additional_data
Definition: petsc_solver.h:587
PETScWrappers::SolverBiCG::set_solver_type
virtual void set_solver_type(KSP &ksp) const override
Definition: petsc_solver.cc:363
PETScWrappers::SparseDirectMUMPS::set_symmetric_mode
void set_symmetric_mode(const bool flag)
Definition: petsc_solver.cc:884
MPI_Comm
PETScWrappers::SolverLSQR::AdditionalData
Definition: petsc_solver.h:820
exceptions.h
PETScWrappers::SolverCGS::SolverCGS
SolverCGS(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
Definition: petsc_solver.cc:472
PETScWrappers::SolverRichardson::AdditionalData
Definition: petsc_solver.h:276
PETScWrappers::SolverCG::set_solver_type
virtual void set_solver_type(KSP &ksp) const override
Definition: petsc_solver.cc:339
PETScWrappers::SolverBase::prefix_name
std::string prefix_name
Definition: petsc_solver.h:202
PETScWrappers::SolverBiCG::AdditionalData
Definition: petsc_solver.h:440
PETScWrappers::SolverTFQMR::additional_data
const AdditionalData additional_data
Definition: petsc_solver.h:688
PETScWrappers::SolverBase::control
SolverControl & control() const
Definition: petsc_solver.cc:168
PETScWrappers::SolverCR::SolverCR
SolverCR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
Definition: petsc_solver.cc:544
PETScWrappers::SolverGMRES::AdditionalData::AdditionalData
AdditionalData(const unsigned int restart_parameter=30, const bool right_preconditioning=false)
Definition: petsc_solver.cc:378
PETScWrappers::SolverTFQMR::SolverTFQMR
SolverTFQMR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
Definition: petsc_solver.cc:496
SLEPcWrappers::TransformationBase
Definition: slepc_spectral_transformation.h:74
PETScWrappers::SolverBase::convergence_test
static PetscErrorCode convergence_test(KSP ksp, const PetscInt iteration, const PetscReal residual_norm, KSPConvergedReason *reason, void *solver_control)
Definition: petsc_solver.cc:175
PETScWrappers::SolverBase::SolverData
Definition: petsc_solver.h:235
PETScWrappers::SolverRichardson::additional_data
const AdditionalData additional_data
Definition: petsc_solver.h:313
PETScWrappers::SolverGMRES::AdditionalData::right_preconditioning
bool right_preconditioning
Definition: petsc_solver.h:508
PETScWrappers::SolverCGS::AdditionalData
Definition: petsc_solver.h:610
PETScWrappers::SolverTCQMR
Definition: petsc_solver.h:711
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
Physics::Elasticity::Kinematics::b
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
PETScWrappers::SparseDirectMUMPS::SolverDataMUMPS::~SolverDataMUMPS
~SolverDataMUMPS()
Definition: petsc_solver.cc:626
PETScWrappers::MatrixBase
Definition: petsc_matrix_base.h:284
PETScWrappers::SolverPreOnly::AdditionalData
Definition: petsc_solver.h:876
PETScWrappers::SolverCR::additional_data
const AdditionalData additional_data
Definition: petsc_solver.h:795
PETScWrappers::SolverTFQMR
Definition: petsc_solver.h:655
PETScWrappers::SolverTFQMR::AdditionalData
Definition: petsc_solver.h:661
PETScWrappers::SparseDirectMUMPS::set_solver_type
virtual void set_solver_type(KSP &ksp) const override
Definition: petsc_solver.cc:644
SLEPcWrappers
Definition: slepc_solver.h:137
PETScWrappers::SparseDirectMUMPS::SolverDataMUMPS
Definition: petsc_solver.h:994
PETScWrappers::SparseDirectMUMPS::SolverDataMUMPS::ksp
KSP ksp
Definition: petsc_solver.h:1001
PETScWrappers::SolverBiCG::additional_data
const AdditionalData additional_data
Definition: petsc_solver.h:467
PETScWrappers::SolverRichardson::AdditionalData::omega
double omega
Definition: petsc_solver.h:286
PETScWrappers::SolverBase::solver_control
SolverControl & solver_control
Definition: petsc_solver.h:182
PETScWrappers::SolverCR
Definition: petsc_solver.h:762
PETScWrappers::SolverRichardson
Definition: petsc_solver.h:270
PETScWrappers::SparseDirectMUMPS::SparseDirectMUMPS
SparseDirectMUMPS(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
Definition: petsc_solver.cc:634
PETScWrappers::SparseDirectMUMPS
Definition: petsc_solver.h:937
PETScWrappers::SparseDirectMUMPS::AdditionalData
Definition: petsc_solver.h:943
LAPACKSupport::A
static const char A
Definition: lapack_support.h:155
PETScWrappers::SolverChebychev::AdditionalData
Definition: petsc_solver.h:338
PETScWrappers::SolverTCQMR::AdditionalData
Definition: petsc_solver.h:717
PETScWrappers::SolverCGS::additional_data
const AdditionalData additional_data
Definition: petsc_solver.h:637
slepc_spectral_transformation.h
PETScWrappers::SolverGMRES
Definition: petsc_solver.h:485
PETScWrappers::SolverBase::reset
virtual void reset()
Definition: petsc_solver.cc:161
PETScWrappers::SolverBase::set_prefix
void set_prefix(const std::string &prefix)
Definition: petsc_solver.cc:154
PETScWrappers::SolverRichardson::set_solver_type
virtual void set_solver_type(KSP &ksp) const override
Definition: petsc_solver.cc:267
PETScWrappers::SolverBicgstab::set_solver_type
virtual void set_solver_type(KSP &ksp) const override
Definition: petsc_solver.cc:457
PETScWrappers::SolverBase::solve
void solve(const MatrixBase &A, VectorBase &x, const VectorBase &b, const PreconditionerBase &preconditioner)
Definition: petsc_solver.cc:53
PETScWrappers::SolverChebychev::set_solver_type
virtual void set_solver_type(KSP &ksp) const override
Definition: petsc_solver.cc:315
PETScWrappers::SparseDirectMUMPS::additional_data
const AdditionalData additional_data
Definition: petsc_solver.h:970
PETScWrappers::SolverLSQR::SolverLSQR
SolverLSQR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
Definition: petsc_solver.cc:568
PETScWrappers::SparseDirectMUMPS::convergence_test
static PetscErrorCode convergence_test(KSP ksp, const PetscInt iteration, const PetscReal residual_norm, KSPConvergedReason *reason, void *solver_control)
Definition: petsc_solver.cc:847
PETScWrappers::SolverBase::mpi_communicator
const MPI_Comm mpi_communicator
Definition: petsc_solver.h:187
PETScWrappers::SolverBase::SolverData::ksp
KSP ksp
Definition: petsc_solver.h:245
PETScWrappers::SolverBicgstab
Definition: petsc_solver.h:554
PETScWrappers::SolverCG::SolverCG
SolverCG(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
Definition: petsc_solver.cc:330
config.h
PETScWrappers::PreconditionerBase
Definition: petsc_precondition.h:57
PETScWrappers::SparseDirectMUMPS::SolverDataMUMPS::pc
PC pc
Definition: petsc_solver.h:1002
PETScWrappers::SolverPreOnly::set_solver_type
virtual void set_solver_type(KSP &ksp) const override
Definition: petsc_solver.cc:601
PETScWrappers::SolverRichardson::SolverRichardson
SolverRichardson(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
Definition: petsc_solver.cc:258
PETScWrappers::SolverBase::solver_data
std::unique_ptr< SolverData > solver_data
Definition: petsc_solver.h:252
SolverControl
Definition: solver_control.h:67
PETScWrappers::SolverPreOnly::additional_data
const AdditionalData additional_data
Definition: petsc_solver.h:903
PETScWrappers::SolverBase::initialize
void initialize(const PreconditionerBase &preconditioner)
Definition: petsc_solver.cc:213
PETScWrappers::SolverTFQMR::set_solver_type
virtual void set_solver_type(KSP &ksp) const override
Definition: petsc_solver.cc:505
PETScWrappers::SolverPreOnly
Definition: petsc_solver.h:870
PETScWrappers::SolverBiCG::SolverBiCG
SolverBiCG(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
Definition: petsc_solver.cc:354
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
PETScWrappers::SolverLSQR
Definition: petsc_solver.h:814
PETScWrappers::SolverLSQR::additional_data
const AdditionalData additional_data
Definition: petsc_solver.h:847
PETScWrappers::SolverGMRES::AdditionalData
Definition: petsc_solver.h:491
PETScWrappers::VectorBase
Definition: petsc_vector_base.h:240
PETScWrappers::SolverCR::AdditionalData
Definition: petsc_solver.h:768
PETScWrappers::SolverGMRES::SolverGMRES
SolverGMRES(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
Definition: petsc_solver.cc:387
PETScWrappers::SolverCGS
Definition: petsc_solver.h:604
PETScWrappers::SolverRichardson::AdditionalData::AdditionalData
AdditionalData(const double omega=1)
Definition: petsc_solver.cc:252
PETScWrappers::SolverChebychev
Definition: petsc_solver.h:332
PETScWrappers::SolverBicgstab::AdditionalData
Definition: petsc_solver.h:560
PETScWrappers::SolverBase::SolverBase
SolverBase(SolverControl &cn, const MPI_Comm &mpi_communicator)
Definition: petsc_solver.cc:45
PETScWrappers::SolverLSQR::set_solver_type
virtual void set_solver_type(KSP &ksp) const override
Definition: petsc_solver.cc:577
solver_control.h
PETScWrappers::SparseDirectMUMPS::symmetric_mode
bool symmetric_mode
Definition: petsc_solver.h:1011