libMesh
petsc_linear_solver.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef LIBMESH_PETSC_LINEAR_SOLVER_H
21 #define LIBMESH_PETSC_LINEAR_SOLVER_H
22 
23 #include "libmesh/libmesh_config.h"
24 
25 #ifdef LIBMESH_HAVE_PETSC
26 
27 #include "libmesh/petsc_macro.h"
28 #include "libmesh/petsc_solver_exception.h"
29 #include "libmesh/wrapped_petsc.h"
30 
31 // Petsc include files.
32 #ifdef I
33 # define LIBMESH_SAW_I
34 #endif
35 #include <petscksp.h>
36 #ifndef LIBMESH_SAW_I
37 # undef I // Avoid complex.h contamination
38 #endif
39 
40 // Local includes
41 #include "libmesh/linear_solver.h"
42 
43 // C++ includes
44 #include <cstddef>
45 #include <vector>
46 #include <memory>
47 
48 //--------------------------------------------------------------------
49 // Functions with C linkage to pass to PETSc. PETSc will call these
50 // methods as needed for preconditioning
51 //
52 // Since they must have C linkage they have no knowledge of a namespace.
53 // Give them an obscure name to avoid namespace pollution.
54 extern "C"
55 {
60  PetscErrorCode libmesh_petsc_preconditioner_setup (PC);
61 
66  PetscErrorCode libmesh_petsc_preconditioner_apply(PC, Vec x, Vec y);
67 
68 #if LIBMESH_ENABLE_DEPRECATED
69 
73  PetscErrorCode __libmesh_petsc_preconditioner_setup (PC);
74 
79  PetscErrorCode __libmesh_petsc_preconditioner_apply(PC, Vec x, Vec y);
80 #endif
81 } // end extern "C"
82 
83 
84 namespace libMesh
85 {
86 
87 // forward declarations
88 template <typename T> class PetscMatrixBase;
89 
98 template <typename T>
100 {
101 public:
106 
110  virtual ~PetscLinearSolver () = default;
111 
115  virtual void clear () override;
116 
122  virtual void init (const char * name = nullptr) override;
123 
127  void init (PetscMatrixBase<T> * matrix,
128  const char * name = nullptr);
129 
138  virtual void init_names (const System &) override;
139 
147  virtual void restrict_solve_to (const std::vector<unsigned int> * const dofs,
148  const SubsetSolveMode subset_solve_mode=SUBSET_ZERO) override;
149 
154  virtual std::pair<unsigned int, Real>
155  solve (SparseMatrix<T> & matrix_in,
156  NumericVector<T> & solution_in,
157  NumericVector<T> & rhs_in,
158  const std::optional<double> tol = std::nullopt,
159  const std::optional<unsigned int> m_its = std::nullopt) override
160  {
161  return this->solve(matrix_in, matrix_in, solution_in, rhs_in, tol, m_its);
162  }
163 
164 
169  virtual std::pair<unsigned int, Real>
170  adjoint_solve (SparseMatrix<T> & matrix_in,
171  NumericVector<T> & solution_in,
172  NumericVector<T> & rhs_in,
173  const std::optional<double> tol = std::nullopt,
174  const std::optional<unsigned int> m_its = std::nullopt) override;
175 
192  virtual std::pair<unsigned int, Real>
193  solve (SparseMatrix<T> & matrix,
194  SparseMatrix<T> & preconditioner,
195  NumericVector<T> & solution,
196  NumericVector<T> & rhs,
197  const std::optional<double> tol = std::nullopt,
198  const std::optional<unsigned int> m_its = std::nullopt) override;
199 
203  virtual std::pair<unsigned int, Real>
204  solve (const ShellMatrix<T> & shell_matrix,
205  NumericVector<T> & solution_in,
206  NumericVector<T> & rhs_in,
207  const std::optional<double> tol = std::nullopt,
208  const std::optional<unsigned int> m_its = std::nullopt) override;
209 
215  virtual std::pair<unsigned int, Real>
216  solve (const ShellMatrix<T> & shell_matrix,
217  const SparseMatrix<T> & precond_matrix,
218  NumericVector<T> & solution_in,
219  NumericVector<T> & rhs_in,
220  const std::optional<double> tol = std::nullopt,
221  const std::optional<unsigned int> m_its = std::nullopt) override;
222 
230  PC pc() { this->init(); return _pc; }
231 
238  KSP ksp();
239 
244  void get_residual_history(std::vector<double> & hist);
245 
254 
258  virtual LinearConvergenceReason get_converged_reason() const override;
259 
260 private:
261 
262  typedef PetscErrorCode (*ksp_solve_func_type)(KSP,Vec,Vec);
263 
264  /*
265  * Helper function to run solve() and adjoint_solve()
266  */
267  virtual std::pair<unsigned int, Real>
268  solve_common (SparseMatrix<T> & matrix_in,
269  SparseMatrix<T> & precond_in,
270  NumericVector<T> & solution_in,
271  NumericVector<T> & rhs_in,
272  const double rel_tol,
273  const double abs_tol,
274  const unsigned int m_its,
275  ksp_solve_func_type solve_func);
276 
277  /*
278  * Helper function to run ShellMatrix solve() with or without a
279  * separate preconditioner matrix
280  */
281  virtual std::pair<unsigned int, Real>
282  shell_solve_common (const ShellMatrix<T> & shell_matrix,
283  PetscMatrixBase<T> * precond_matrix,
284  NumericVector<T> & solution_in,
285  NumericVector<T> & rhs_in,
286  const double rel_tol,
287  const double abs_tol,
288  const unsigned int m_its);
289 
290  /*
291  * Helper function for the helper functions
292  */
293  std::pair<unsigned int, Real>
294  solve_base (SparseMatrix<T> * matrix,
295  PetscMatrixBase<T> * precond,
296  Mat mat,
297  NumericVector<T> & solution_in,
298  NumericVector<T> & rhs_in,
299  const double rel_tol,
300  const double abs_tol,
301  const unsigned int m_its,
302  ksp_solve_func_type solve_func);
303 
304 
309  void set_petsc_solver_type ();
310 
314  static PetscErrorCode _petsc_shell_matrix_mult(Mat mat, Vec arg, Vec dest);
315 
319  static PetscErrorCode _petsc_shell_matrix_mult_add(Mat mat, Vec arg, Vec add, Vec dest);
320 
324  static PetscErrorCode _petsc_shell_matrix_get_diagonal(Mat mat, Vec dest);
325 
329  PC _pc;
330 
335 
341 
348 
352  PetscInt restrict_solve_to_is_local_size() const;
353 
359  void create_complement_is (const NumericVector<T> & vec_in);
360 
366 };
367 
368 } // namespace libMesh
369 
370 
371 #endif // #ifdef LIBMESH_HAVE_PETSC
372 #endif // LIBMESH_PETSC_LINEAR_SOLVER_H
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
void get_residual_history(std::vector< double > &hist)
Fills the input vector with the sequence of residual norms from the latest iterative solve...
This class provides a nice interface to the PETSc C-based data structures for parallel, sparse matrices.
void create_complement_is(const NumericVector< T > &vec_in)
Creates _restrict_solve_to_is_complement to contain all indices that are local in vec_in...
WrappedPetsc< IS > _restrict_solve_to_is
PETSc index set containing the dofs on which to solve (nullptr means solve on all dofs)...
virtual void init(const char *name=nullptr) override
Initialize data structures if not done so already.
PetscErrorCode libmesh_petsc_preconditioner_apply(PC, Vec x, Vec y)
This function is called by PETSc to actually apply the preconditioner.
static PetscErrorCode _petsc_shell_matrix_mult_add(Mat mat, Vec arg, Vec add, Vec dest)
Internal function if shell matrix mode is used.
LinearConvergenceReason
Linear solver convergence flags (taken from the PETSc flags).
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
Definition: vector_fe_ex5.C:44
std::pair< unsigned int, Real > solve_base(SparseMatrix< T > *matrix, PetscMatrixBase< T > *precond, Mat mat, NumericVector< T > &solution_in, NumericVector< T > &rhs_in, const double rel_tol, const double abs_tol, const unsigned int m_its, ksp_solve_func_type solve_func)
PetscInt restrict_solve_to_is_local_size() const
virtual std::pair< unsigned int, Real > solve(SparseMatrix< T > &matrix_in, NumericVector< T > &solution_in, NumericVector< T > &rhs_in, const std::optional< double > tol=std::nullopt, const std::optional< unsigned int > m_its=std::nullopt) override
Call the Petsc solver.
The libMesh namespace provides an interface to certain functionality in the library.
virtual std::pair< unsigned int, Real > adjoint_solve(SparseMatrix< T > &matrix_in, NumericVector< T > &solution_in, NumericVector< T > &rhs_in, const std::optional< double > tol=std::nullopt, const std::optional< unsigned int > m_its=std::nullopt) override
Call the Petsc solver.
PetscErrorCode(* ksp_solve_func_type)(KSP, Vec, Vec)
Generic sparse matrix.
Definition: vector_fe_ex5.C:46
This base class can be inherited from to provide interfaces to linear solvers from different packages...
virtual std::pair< unsigned int, Real > solve_common(SparseMatrix< T > &matrix_in, SparseMatrix< T > &precond_in, NumericVector< T > &solution_in, NumericVector< T > &rhs_in, const double rel_tol, const double abs_tol, const unsigned int m_its, ksp_solve_func_type solve_func)
PetscLinearSolver(const libMesh::Parallel::Communicator &comm_in)
Constructor.
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:96
void set_petsc_solver_type()
Tells PETSC to use the user-specified solver stored in _solver_type.
virtual void restrict_solve_to(const std::vector< unsigned int > *const dofs, const SubsetSolveMode subset_solve_mode=SUBSET_ZERO) override
After calling this method, all successive solves will be restricted to the given set of dofs...
virtual ~PetscLinearSolver()=default
Destructor.
PC _pc
Preconditioner context.
PetscErrorCode __libmesh_petsc_preconditioner_apply(PC, Vec x, Vec y)
This function is called by PETSc to actually apply the preconditioner.
virtual std::pair< unsigned int, Real > shell_solve_common(const ShellMatrix< T > &shell_matrix, PetscMatrixBase< T > *precond_matrix, NumericVector< T > &solution_in, NumericVector< T > &rhs_in, const double rel_tol, const double abs_tol, const unsigned int m_its)
virtual void init_names(const System &) override
Apply names to the system to be solved.
static PetscErrorCode _petsc_shell_matrix_mult(Mat mat, Vec arg, Vec dest)
Internal function if shell matrix mode is used.
PetscErrorCode libmesh_petsc_preconditioner_setup(PC)
This function is called by PETSc to initialize the preconditioner.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
PetscErrorCode __libmesh_petsc_preconditioner_setup(PC)
This function is called by PETSc to initialize the preconditioner.
This class provides an interface to PETSc iterative solvers that is compatible with the libMesh Linea...
SubsetSolveMode _subset_solve_mode
If restrict-solve-to-subset mode is active, this member decides what happens with the dofs outside th...
virtual LinearConvergenceReason get_converged_reason() const override
WrappedPetsc< IS > _restrict_solve_to_is_complement
PETSc index set, complement to _restrict_solve_to_is.
Generic shell matrix, i.e.
SubsetSolveMode
defines an enum for the question what happens to the dofs outside the given subset when a system is s...
static PetscErrorCode _petsc_shell_matrix_get_diagonal(Mat mat, Vec dest)
Internal function if shell matrix mode is used.
WrappedPetsc< KSP > _ksp
Krylov subspace context.
virtual void clear() override
Release all memory and clear data structures.