20 #ifndef LIBMESH_PETSC_LINEAR_SOLVER_H 21 #define LIBMESH_PETSC_LINEAR_SOLVER_H 23 #include "libmesh/libmesh_config.h" 25 #ifdef LIBMESH_HAVE_PETSC 27 #include "libmesh/petsc_macro.h" 28 #include "libmesh/petsc_solver_exception.h" 29 #include "libmesh/wrapped_petsc.h" 33 # define LIBMESH_SAW_I 37 # undef I // Avoid complex.h contamination 41 #include "libmesh/linear_solver.h" 68 #if LIBMESH_ENABLE_DEPRECATED 88 template <
typename T>
class PetscMatrixBase;
115 virtual void clear ()
override;
122 virtual void init (
const char *
name =
nullptr)
override;
128 const char *
name =
nullptr);
154 virtual std::pair<unsigned int, Real>
158 const std::optional<double> tol = std::nullopt,
159 const std::optional<unsigned int> m_its = std::nullopt)
override 161 return this->
solve(matrix_in, matrix_in, solution_in, rhs_in, tol, m_its);
169 virtual std::pair<unsigned int, Real>
173 const std::optional<double> tol = std::nullopt,
174 const std::optional<unsigned int> m_its = std::nullopt)
override;
192 virtual std::pair<unsigned int, Real>
197 const std::optional<double> tol = std::nullopt,
198 const std::optional<unsigned int> m_its = std::nullopt)
override;
203 virtual std::pair<unsigned int, Real>
207 const std::optional<double> tol = std::nullopt,
208 const std::optional<unsigned int> m_its = std::nullopt)
override;
215 virtual std::pair<unsigned int, Real>
220 const std::optional<double> tol = std::nullopt,
221 const std::optional<unsigned int> m_its = std::nullopt)
override;
267 virtual std::pair<unsigned int, Real>
272 const double rel_tol,
273 const double abs_tol,
274 const unsigned int m_its,
281 virtual std::pair<unsigned int, Real>
286 const double rel_tol,
287 const double abs_tol,
288 const unsigned int m_its);
293 std::pair<unsigned int, Real>
299 const double rel_tol,
300 const double abs_tol,
301 const unsigned int m_its,
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.
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...
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)
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.
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...
Real get_initial_residual()
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.