Reference documentation for deal.II version 9.6.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\}}\)
Loading...
Searching...
No Matches
petsc_solver.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2004 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_petsc_solver_h
16#define dealii_petsc_solver_h
17
18
19#include <deal.II/base/config.h>
20
21#ifdef DEAL_II_WITH_PETSC
22
24
27
28# include <petscksp.h>
29
30# ifdef DEAL_II_WITH_SLEPC
32# endif
33
35
36// Forward declarations
37# ifndef DOXYGEN
38# ifdef DEAL_II_WITH_SLEPC
39namespace SLEPcWrappers
40{
41 // forward declarations
42 class TransformationBase;
43} // namespace SLEPcWrappers
44# endif
45# endif
46
47namespace PETScWrappers
48{
49 // forward declarations
50# ifndef DOXYGEN
51 class MatrixBase;
52 class VectorBase;
53 class PreconditionBase;
54# endif
55
56
92 {
93 public:
98 SolverBase();
99
104
108 virtual ~SolverBase();
109
117 void
118 solve(const MatrixBase &A,
119 VectorBase &x,
120 const VectorBase &b,
121 const PreconditionBase &preconditioner);
122
127 virtual void
128 reset();
129
134 void
135 set_prefix(const std::string &prefix);
136
141 control() const;
142
147 void
148 initialize(const PreconditionBase &preconditioner);
149
153 KSP
154 petsc_ksp();
155
161 operator KSP() const;
162
163 protected:
167 KSP ksp;
168
176
180 void
182
187 virtual void
188 set_solver_type(KSP &ksp) const;
189
196 void
198
205 std::string prefix_name;
206
207 private:
214 static PetscErrorCode
216 const PetscInt iteration,
217 const PetscReal residual_norm,
218 KSPConvergedReason *reason,
219 void *solver_control);
220
221
222# ifdef DEAL_II_WITH_SLEPC
223 // Make the transformation class a friend, since it needs to set the KSP
224 // solver.
226# endif
227 };
228
229
230
238 {
239 public:
244 {
248 explicit AdditionalData(const double omega = 1);
249
253 double omega;
254 };
255
264 const AdditionalData &data = AdditionalData());
265
274 const MPI_Comm mpi_communicator,
275 const AdditionalData &data = AdditionalData());
276
277 protected:
282
287 virtual void
288 set_solver_type(KSP &ksp) const override;
289 };
290
291
292
300 {
301 public:
306 {};
307
316 const AdditionalData &data = AdditionalData());
317
326 const MPI_Comm mpi_communicator,
327 const AdditionalData &data = AdditionalData());
328
329 protected:
334
339 virtual void
340 set_solver_type(KSP &ksp) const override;
341 };
342
343
344
350 class SolverCG : public SolverBase
351 {
352 public:
357 {};
358
367
376 const MPI_Comm mpi_communicator,
377 const AdditionalData &data = AdditionalData());
378
379 protected:
384
389 virtual void
390 set_solver_type(KSP &ksp) const override;
391 };
392
393
394
400 class SolverBiCG : public SolverBase
401 {
402 public:
407 {};
408
417 const AdditionalData &data = AdditionalData());
418
427 const MPI_Comm mpi_communicator,
428 const AdditionalData &data = AdditionalData());
429
430 protected:
435
440 virtual void
441 set_solver_type(KSP &ksp) const override;
442 };
443
444
445
451 class SolverGMRES : public SolverBase
452 {
453 public:
458 {
463 AdditionalData(const unsigned int restart_parameter = 30,
464 const bool right_preconditioning = false);
465
469 unsigned int restart_parameter;
470
475 };
476
485 const AdditionalData &data = AdditionalData());
486
495 const MPI_Comm mpi_communicator,
496 const AdditionalData &data = AdditionalData());
497
498 protected:
503
508 virtual void
509 set_solver_type(KSP &ksp) const override;
510 };
511
512
513
521 {
522 public:
527 {};
528
537 const AdditionalData &data = AdditionalData());
538
547 const MPI_Comm mpi_communicator,
548 const AdditionalData &data = AdditionalData());
549
550 protected:
555
560 virtual void
561 set_solver_type(KSP &ksp) const override;
562 };
563
564
565
572 class SolverCGS : public SolverBase
573 {
574 public:
579 {};
580
589
598 const MPI_Comm mpi_communicator,
599 const AdditionalData &data = AdditionalData());
600
601 protected:
606
611 virtual void
612 set_solver_type(KSP &ksp) const override;
613 };
614
615
616
622 class SolverTFQMR : public SolverBase
623 {
624 public:
629 {};
630
639 const AdditionalData &data = AdditionalData());
640
649 const MPI_Comm mpi_communicator,
650 const AdditionalData &data = AdditionalData());
651
652 protected:
657
662 virtual void
663 set_solver_type(KSP &ksp) const override;
664 };
665
666
667
678 class SolverTCQMR : public SolverBase
679 {
680 public:
685 {};
686
695 const AdditionalData &data = AdditionalData());
696
705 const MPI_Comm mpi_communicator,
706 const AdditionalData &data = AdditionalData());
707
708 protected:
713
718 virtual void
719 set_solver_type(KSP &ksp) const override;
720 };
721
722
723
729 class SolverCR : public SolverBase
730 {
731 public:
736 {};
737
746
755 const MPI_Comm mpi_communicator,
756 const AdditionalData &data = AdditionalData());
757
758 protected:
763
768 virtual void
769 set_solver_type(KSP &ksp) const override;
770 };
771
772
773
780 class SolverLSQR : public SolverBase
781 {
782 public:
787 {};
788
797 const AdditionalData &data = AdditionalData());
798
807 const MPI_Comm mpi_communicator,
808 const AdditionalData &data = AdditionalData());
809
810 protected:
815
820 virtual void
821 set_solver_type(KSP &ksp) const override;
822 };
823
824
837 {
838 public:
843 {};
844
853 const AdditionalData &data = AdditionalData());
854
863 const MPI_Comm mpi_communicator,
864 const AdditionalData &data = AdditionalData());
865
866 protected:
871
876 virtual void
877 set_solver_type(KSP &ksp) const override;
878 };
879
904 {
905 public:
910 {};
911
916 const AdditionalData &data = AdditionalData());
917
926 const MPI_Comm mpi_communicator,
927 const AdditionalData &data = AdditionalData());
928
932 void
933 solve(const MatrixBase &A, VectorBase &x, const VectorBase &b);
934
958 void
959 set_symmetric_mode(const bool matrix_is_symmetric);
960
961 protected:
966
967 virtual void
968 set_solver_type(KSP &ksp) const override;
969
975 };
976} // namespace PETScWrappers
977
979
980#endif // DEAL_II_WITH_PETSC
981
982#endif
SolverControl & control() const
void initialize(const PreconditionBase &preconditioner)
void perhaps_set_convergence_test() const
static PetscErrorCode convergence_test(KSP ksp, const PetscInt iteration, const PetscReal residual_norm, KSPConvergedReason *reason, void *solver_control)
SmartPointer< SolverControl, SolverBase > solver_control
void set_prefix(const std::string &prefix)
void initialize_ksp_with_comm(const MPI_Comm comm)
void solve(const MatrixBase &A, VectorBase &x, const VectorBase &b, const PreconditionBase &preconditioner)
virtual void set_solver_type(KSP &ksp) const
SolverBiCG(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
virtual void set_solver_type(KSP &ksp) const override
virtual void set_solver_type(KSP &ksp) const override
const AdditionalData additional_data
SolverBicgstab(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
virtual void set_solver_type(KSP &ksp) const override
SolverCGS(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
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
SolverCR(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
virtual void set_solver_type(KSP &ksp) const override
SolverChebychev(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
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
SolverPreOnly(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
virtual void set_solver_type(KSP &ksp) const override
virtual void set_solver_type(KSP &ksp) const override
const AdditionalData additional_data
SolverRichardson(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
virtual void set_solver_type(KSP &ksp) const override
SolverTCQMR(SolverControl &cn, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
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
void set_symmetric_mode(const bool matrix_is_symmetric)
SparseDirectMUMPS(SolverControl &cn, const AdditionalData &data=AdditionalData())
#define DEAL_II_DEPRECATED
Definition config.h:207
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
*braid_SplitCommworld & comm
AdditionalData(const unsigned int restart_parameter=30, const bool right_preconditioning=false)