Reference documentation for deal.II version 9.4.1
\(\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\}}\)
Loading...
Searching...
No Matches
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
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
40namespace SLEPcWrappers
41{
42 // forward declarations
43 class TransformationBase;
44} // namespace SLEPcWrappers
45# endif
46# endif
47
48namespace PETScWrappers
49{
50 // forward declarations
51# ifndef DOXYGEN
52 class MatrixBase;
53 class VectorBase;
54 class PreconditionBase;
55# endif
56
57
113 {
114 public:
125
129 virtual ~SolverBase() = default;
130
138 void
139 solve(const MatrixBase & A,
140 VectorBase & x,
141 const VectorBase & b,
142 const PreconditionBase &preconditioner);
143
144
149 virtual void
150 reset();
151
152
157 void
158 set_prefix(const std::string &prefix);
159
160
165 control() const;
166
171 void
172 initialize(const PreconditionBase &preconditioner);
173
174 protected:
182
187
192 virtual void
193 set_solver_type(KSP &ksp) const = 0;
194
201 std::string prefix_name;
202
203 private:
210 static PetscErrorCode
211 convergence_test(KSP ksp,
212 const PetscInt iteration,
213 const PetscReal residual_norm,
214 KSPConvergedReason *reason,
215 void * solver_control);
216
235 {
239 ~SolverData();
240
244 KSP ksp;
245 };
246
251 std::unique_ptr<SolverData> solver_data;
252
253# ifdef DEAL_II_WITH_SLEPC
254 // Make the transformation class a friend, since it needs to set the KSP
255 // solver.
257# endif
258 };
259
260
261
269 {
270 public:
275 {
279 explicit AdditionalData(const double omega = 1);
280
284 double omega;
285 };
286
304 const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
305 const AdditionalData &data = AdditionalData());
306
307 protected:
312
317 virtual void
318 set_solver_type(KSP &ksp) const override;
319 };
320
321
322
330 {
331 public:
336 {};
337
355 const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
356 const AdditionalData &data = AdditionalData());
357
358 protected:
363
368 virtual void
369 set_solver_type(KSP &ksp) const override;
370 };
371
372
373
379 class SolverCG : public SolverBase
380 {
381 public:
386 {};
387
405 const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
406 const AdditionalData &data = AdditionalData());
407
408 protected:
413
418 virtual void
419 set_solver_type(KSP &ksp) const override;
420 };
421
422
423
429 class SolverBiCG : public SolverBase
430 {
431 public:
436 {};
437
455 const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
456 const AdditionalData &data = AdditionalData());
457
458 protected:
463
468 virtual void
469 set_solver_type(KSP &ksp) const override;
470 };
471
472
473
479 class SolverGMRES : public SolverBase
480 {
481 public:
486 {
491 AdditionalData(const unsigned int restart_parameter = 30,
492 const bool right_preconditioning = false);
493
497 unsigned int restart_parameter;
498
503 };
504
522 const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
523 const AdditionalData &data = AdditionalData());
524
525 protected:
530
535 virtual void
536 set_solver_type(KSP &ksp) const override;
537 };
538
539
540
548 {
549 public:
554 {};
555
573 const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
574 const AdditionalData &data = AdditionalData());
575
576 protected:
581
586 virtual void
587 set_solver_type(KSP &ksp) const override;
588 };
589
596 class SolverCGS : public SolverBase
597 {
598 public:
603 {};
604
622 const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
623 const AdditionalData &data = AdditionalData());
624
625 protected:
630
635 virtual void
636 set_solver_type(KSP &ksp) const override;
637 };
638
639
640
646 class SolverTFQMR : public SolverBase
647 {
648 public:
653 {};
654
672 const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
673 const AdditionalData &data = AdditionalData());
674
675 protected:
680
685 virtual void
686 set_solver_type(KSP &ksp) const override;
687 };
688
689
690
701 class SolverTCQMR : public SolverBase
702 {
703 public:
708 {};
709
727 const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
728 const AdditionalData &data = AdditionalData());
729
730 protected:
735
740 virtual void
741 set_solver_type(KSP &ksp) const override;
742 };
743
744
745
751 class SolverCR : public SolverBase
752 {
753 public:
758 {};
759
777 const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
778 const AdditionalData &data = AdditionalData());
779
780 protected:
785
790 virtual void
791 set_solver_type(KSP &ksp) const override;
792 };
793
794
795
802 class SolverLSQR : public SolverBase
803 {
804 public:
809 {};
810
828 const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
829 const AdditionalData &data = AdditionalData());
830
831 protected:
836
841 virtual void
842 set_solver_type(KSP &ksp) const override;
843 };
844
845
858 {
859 public:
864 {};
865
883 const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
884 const AdditionalData &data = AdditionalData());
885
886 protected:
891
896 virtual void
897 set_solver_type(KSP &ksp) const override;
898 };
899
924 {
925 public:
930 {};
935 const MPI_Comm & mpi_communicator = PETSC_COMM_SELF,
936 const AdditionalData &data = AdditionalData());
937
941 void
942 solve(const MatrixBase &A, VectorBase &x, const VectorBase &b);
943
949 void
950 set_symmetric_mode(const bool flag);
951
952 protected:
957
958 virtual void
959 set_solver_type(KSP &ksp) const override;
960
961 private:
968 static PetscErrorCode
969 convergence_test(KSP ksp,
970 const PetscInt iteration,
971 const PetscReal residual_norm,
972 KSPConvergedReason *reason,
973 void * solver_control);
974
981 {
986
987 KSP ksp;
988 PC pc;
989 };
990
991 std::unique_ptr<SolverDataMUMPS> solver_data;
992
998 };
999} // namespace PETScWrappers
1000
1002
1003# endif // DEAL_II_WITH_PETSC
1004
1005/*---------------------------- petsc_solver.h ---------------------------*/
1006
1007#endif
1008/*---------------------------- petsc_solver.h ---------------------------*/
virtual void set_solver_type(KSP &ksp) const =0
SolverControl & control() const
SolverControl & solver_control
Definition: petsc_solver.h:181
void initialize(const PreconditionBase &preconditioner)
static PetscErrorCode convergence_test(KSP ksp, const PetscInt iteration, const PetscReal residual_norm, KSPConvergedReason *reason, void *solver_control)
virtual ~SolverBase()=default
const MPI_Comm mpi_communicator
Definition: petsc_solver.h:186
void set_prefix(const std::string &prefix)
void solve(const MatrixBase &A, VectorBase &x, const VectorBase &b, const PreconditionBase &preconditioner)
Definition: petsc_solver.cc:53
std::unique_ptr< SolverData > solver_data
Definition: petsc_solver.h:251
const AdditionalData additional_data
Definition: petsc_solver.h:462
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:580
const AdditionalData additional_data
Definition: petsc_solver.h:629
virtual void set_solver_type(KSP &ksp) const override
const AdditionalData additional_data
Definition: petsc_solver.h:412
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:784
const AdditionalData additional_data
Definition: petsc_solver.h:362
virtual void set_solver_type(KSP &ksp) const override
const AdditionalData additional_data
Definition: petsc_solver.h:529
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:835
const AdditionalData additional_data
Definition: petsc_solver.h:890
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:311
const AdditionalData additional_data
Definition: petsc_solver.h:734
virtual void set_solver_type(KSP &ksp) const override
const AdditionalData additional_data
Definition: petsc_solver.h:679
virtual void set_solver_type(KSP &ksp) const override
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:956
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)
std::unique_ptr< SolverDataMUMPS > solver_data
Definition: petsc_solver.h:991
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443