libMesh
petsc_linear_solver.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 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 
30 // Petsc include files.
31 #ifdef I
32 # define LIBMESH_SAW_I
33 #endif
34 #include <petscksp.h>
35 #ifndef LIBMESH_SAW_I
36 # undef I // Avoid complex.h contamination
37 #endif
38 
39 // Local includes
40 #include "libmesh/linear_solver.h"
41 
42 // C++ includes
43 #include <cstddef>
44 #include <vector>
45 
46 //--------------------------------------------------------------------
47 // Functions with C linkage to pass to PETSc. PETSc will call these
48 // methods as needed for preconditioning
49 //
50 // Since they must have C linkage they have no knowledge of a namespace.
51 // Give them an obscure name to avoid namespace pollution.
52 extern "C"
53 {
54 #if PETSC_RELEASE_LESS_THAN(3,0,1)
55 
59  PetscErrorCode libmesh_petsc_preconditioner_setup (void * ctx);
60 
65  PetscErrorCode libmesh_petsc_preconditioner_apply(void * ctx, Vec x, Vec y);
66 #else
67  PetscErrorCode libmesh_petsc_preconditioner_setup (PC);
68  PetscErrorCode libmesh_petsc_preconditioner_apply(PC, Vec x, Vec y);
69 #endif
70 
71 #if LIBMESH_ENABLE_DEPRECATED
72 #if PETSC_RELEASE_LESS_THAN(3,0,1)
73 
77  PetscErrorCode __libmesh_petsc_preconditioner_setup (void * ctx);
78 
83  PetscErrorCode __libmesh_petsc_preconditioner_apply(void * ctx, Vec x, Vec y);
84 #else
85  PetscErrorCode __libmesh_petsc_preconditioner_setup (PC);
86  PetscErrorCode __libmesh_petsc_preconditioner_apply(PC, Vec x, Vec y);
87 #endif
88 #endif
89 } // end extern "C"
90 
91 
92 namespace libMesh
93 {
94 
95 // forward declarations
96 template <typename T> class PetscMatrix;
97 
106 template <typename T>
108 {
109 public:
113  PetscLinearSolver (const libMesh::Parallel::Communicator & comm_in);
114 
119 
123  virtual void clear () override;
124 
130  virtual void init (const char * name = nullptr) override;
131 
135  void init (PetscMatrix<T> * matrix,
136  const char * name = nullptr);
137 
146  virtual void init_names (const System &) override;
147 
155  virtual void restrict_solve_to (const std::vector<unsigned int> * const dofs,
156  const SubsetSolveMode subset_solve_mode=SUBSET_ZERO) override;
157 
162  virtual std::pair<unsigned int, Real>
163  solve (SparseMatrix<T> & matrix_in,
164  NumericVector<T> & solution_in,
165  NumericVector<T> & rhs_in,
166  const double tol,
167  const unsigned int m_its) override
168  {
169  return this->solve(matrix_in, matrix_in, solution_in, rhs_in, tol, m_its);
170  }
171 
172 
177  virtual std::pair<unsigned int, Real>
178  adjoint_solve (SparseMatrix<T> & matrix_in,
179  NumericVector<T> & solution_in,
180  NumericVector<T> & rhs_in,
181  const double tol,
182  const unsigned int m_its) override;
183 
200  virtual std::pair<unsigned int, Real>
201  solve (SparseMatrix<T> & matrix,
202  SparseMatrix<T> & preconditioner,
203  NumericVector<T> & solution,
204  NumericVector<T> & rhs,
205  const double tol,
206  const unsigned int m_its) override;
207 
211  virtual std::pair<unsigned int, Real>
212  solve (const ShellMatrix<T> & shell_matrix,
213  NumericVector<T> & solution_in,
214  NumericVector<T> & rhs_in,
215  const double tol,
216  const unsigned int m_its) override;
217 
223  virtual std::pair<unsigned int, Real>
224  solve (const ShellMatrix<T> & shell_matrix,
225  const SparseMatrix<T> & precond_matrix,
226  NumericVector<T> & solution_in,
227  NumericVector<T> & rhs_in,
228  const double tol,
229  const unsigned int m_its) override;
230 
238  PC pc() { this->init(); return _pc; }
239 
246  KSP ksp() { this->init(); return _ksp; }
247 
252  void get_residual_history(std::vector<double> & hist);
253 
262 
266  virtual LinearConvergenceReason get_converged_reason() const override;
267 
268 private:
269 
274  void set_petsc_solver_type ();
275 
279  static PetscErrorCode _petsc_shell_matrix_mult(Mat mat, Vec arg, Vec dest);
280 
284  static PetscErrorCode _petsc_shell_matrix_mult_add(Mat mat, Vec arg, Vec add, Vec dest);
285 
289  static PetscErrorCode _petsc_shell_matrix_get_diagonal(Mat mat, Vec dest);
290 
294  PC _pc;
295 
299  KSP _ksp;
300 
306 
313 
317  PetscInt _restrict_solve_to_is_local_size() const;
318 
324  void _create_complement_is (const NumericVector<T> & vec_in);
325 
331 };
332 
333 
334 /*----------------------- functions ----------------------------------*/
335 template <typename T>
336 inline
338 {
339  this->clear ();
340 }
341 
342 
343 
344 template <typename T>
345 inline PetscInt
347 {
348  libmesh_assert(_restrict_solve_to_is);
349 
350  PetscInt s;
351  int ierr = ISGetLocalSize(_restrict_solve_to_is, &s);
352  LIBMESH_CHKERR(ierr);
353 
354  return s;
355 }
356 
357 
358 
359 template <typename T>
360 void
362 #if PETSC_VERSION_LESS_THAN(3,0,0)
363  // unnamed to avoid compiler "unused parameter" warning
364 #else
365  vec_in
366 #endif
367  )
368 {
369  libmesh_assert(_restrict_solve_to_is);
370 #if PETSC_VERSION_LESS_THAN(3,0,0)
371  // No ISComplement in PETSc 2.3.3
372  libmesh_not_implemented();
373 #else
374  if (_restrict_solve_to_is_complement==nullptr)
375  {
376  int ierr = ISComplement(_restrict_solve_to_is,
377  vec_in.first_local_index(),
378  vec_in.last_local_index(),
379  &_restrict_solve_to_is_complement);
380  LIBMESH_CHKERR(ierr);
381  }
382 #endif
383 }
384 
385 
386 
387 } // namespace libMesh
388 
389 
390 #endif // #ifdef LIBMESH_HAVE_PETSC
391 #endif // LIBMESH_PETSC_LINEAR_SOLVER_H
libMesh::ShellMatrix
Generic shell matrix, i.e.
Definition: eigen_preconditioner.h:36
libMesh::System
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:100
libMesh::PetscLinearSolver::_restrict_solve_to_is_local_size
PetscInt _restrict_solve_to_is_local_size() const
Definition: petsc_linear_solver.h:346
libMesh::PetscLinearSolver::get_initial_residual
Real get_initial_residual()
Definition: petsc_linear_solver.C:1571
libMesh::PetscLinearSolver::_petsc_shell_matrix_get_diagonal
static PetscErrorCode _petsc_shell_matrix_get_diagonal(Mat mat, Vec dest)
Internal function if shell matrix mode is used.
Definition: petsc_linear_solver.C:1788
libMesh::PetscLinearSolver::~PetscLinearSolver
~PetscLinearSolver()
Destructor.
Definition: petsc_linear_solver.h:337
libMesh::PetscLinearSolver::_petsc_shell_matrix_mult_add
static PetscErrorCode _petsc_shell_matrix_mult_add(Mat mat, Vec arg, Vec add, Vec dest)
Internal function if shell matrix mode is used.
Definition: petsc_linear_solver.C:1758
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::PetscLinearSolver::get_residual_history
void get_residual_history(std::vector< double > &hist)
Fills the input vector with the sequence of residual norms from the latest iterative solve.
Definition: petsc_linear_solver.C:1538
libMesh::PetscLinearSolver::ksp
KSP ksp()
Definition: petsc_linear_solver.h:246
libMesh::PetscLinearSolver::pc
PC pc()
Definition: petsc_linear_solver.h:238
libMesh::ierr
ierr
Definition: petsc_dm_wrapper.C:72
libMesh::SparseMatrix
Generic sparse matrix.
Definition: vector_fe_ex5.C:45
libMesh::SUBSET_ZERO
Definition: enum_subset_solve_mode.h:37
libmesh_petsc_preconditioner_setup
PetscErrorCode libmesh_petsc_preconditioner_setup(void *ctx)
This function is called by PETSc to initialize the preconditioner.
Definition: petsc_linear_solver.C:49
libMesh::NumericVector
Provides a uniform interface to vector storage schemes for different linear algebra libraries.
Definition: vector_fe_ex5.C:43
libMesh::PetscLinearSolver::_ksp
KSP _ksp
Krylov subspace context.
Definition: petsc_linear_solver.h:299
libMesh::PetscLinearSolver::_restrict_solve_to_is
IS _restrict_solve_to_is
PETSc index set containing the dofs on which to solve (nullptr means solve on all dofs).
Definition: petsc_linear_solver.h:305
libMesh::PetscLinearSolver::_petsc_shell_matrix_mult
static PetscErrorCode _petsc_shell_matrix_mult(Mat mat, Vec arg, Vec dest)
Internal function if shell matrix mode is used.
Definition: petsc_linear_solver.C:1734
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::LinearConvergenceReason
LinearConvergenceReason
Linear solver convergence flags (taken from the PETSc flags).
Definition: enum_convergence_flags.h:33
libMesh::PetscLinearSolver::set_petsc_solver_type
void set_petsc_solver_type()
Tells PETSC to use the user-specified solver stored in _solver_type.
Definition: petsc_linear_solver.C:1601
libMesh::PetscLinearSolver::clear
virtual void clear() override
Release all memory and clear data structures.
Definition: petsc_linear_solver.C:150
__libmesh_petsc_preconditioner_apply
PetscErrorCode __libmesh_petsc_preconditioner_apply(void *ctx, Vec x, Vec y)
This function is called by PETSc to actually apply the preconditioner.
Definition: petsc_linear_solver.C:111
libMesh::PetscLinearSolver
This class provides an interface to PETSc iterative solvers that is compatible with the libMesh Linea...
Definition: petsc_linear_solver.h:107
libMesh::PetscLinearSolver::PetscLinearSolver
PetscLinearSolver(const libMesh::Parallel::Communicator &comm_in)
Constructor.
Definition: petsc_linear_solver.C:135
libMesh::PetscLinearSolver::_subset_solve_mode
SubsetSolveMode _subset_solve_mode
If restrict-solve-to-subset mode is active, this member decides what happens with the dofs outside th...
Definition: petsc_linear_solver.h:330
libMesh::PetscLinearSolver::solve
virtual std::pair< unsigned int, Real > solve(SparseMatrix< T > &matrix_in, NumericVector< T > &solution_in, NumericVector< T > &rhs_in, const double tol, const unsigned int m_its) override
Call the Petsc solver.
Definition: petsc_linear_solver.h:163
libmesh_petsc_preconditioner_apply
PetscErrorCode libmesh_petsc_preconditioner_apply(void *ctx, Vec x, Vec y)
This function is called by PETSc to actually apply the preconditioner.
Definition: petsc_linear_solver.C:62
libMesh::PetscLinearSolver::init
virtual void init(const char *name=nullptr) override
Initialize data structures if not done so already.
Definition: petsc_linear_solver.C:193
libMesh::PetscLinearSolver::get_converged_reason
virtual LinearConvergenceReason get_converged_reason() const override
Definition: petsc_linear_solver.C:1685
libMesh::PetscMatrix
This class provides a nice interface to the PETSc C-based data structures for parallel,...
Definition: petsc_matrix.h:87
libMesh::LinearSolver
This base class can be inherited from to provide interfaces to linear solvers from different packages...
Definition: linear_solver.h:69
libMesh::PetscLinearSolver::restrict_solve_to
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,...
Definition: petsc_linear_solver.C:400
libMesh::PetscLinearSolver::_pc
PC _pc
Preconditioner context.
Definition: petsc_linear_solver.h:294
libMesh::SubsetSolveMode
SubsetSolveMode
Definition: enum_subset_solve_mode.h:35
__libmesh_petsc_preconditioner_setup
PetscErrorCode __libmesh_petsc_preconditioner_setup(void *ctx)
This function is called by PETSc to initialize the preconditioner.
Definition: petsc_linear_solver.C:105
libMesh::ctx
void * ctx
Definition: petsc_dm_wrapper.C:71
libMesh::PetscLinearSolver::_create_complement_is
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,...
Definition: petsc_linear_solver.h:361
libMesh::PetscLinearSolver::_restrict_solve_to_is_complement
IS _restrict_solve_to_is_complement
PETSc index set, complement to _restrict_solve_to_is.
Definition: petsc_linear_solver.h:312
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::PetscLinearSolver::init_names
virtual void init_names(const System &) override
Apply names to the system to be solved.
Definition: petsc_linear_solver.C:391
libMesh::PetscLinearSolver::adjoint_solve
virtual std::pair< unsigned int, Real > adjoint_solve(SparseMatrix< T > &matrix_in, NumericVector< T > &solution_in, NumericVector< T > &rhs_in, const double tol, const unsigned int m_its) override
Call the Petsc solver.
Definition: petsc_linear_solver.C:709
libMesh::Quality::name
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42