20 #ifndef LIBMESH_EIGEN_SPARSE_LINEAR_SOLVER_H 21 #define LIBMESH_EIGEN_SPARSE_LINEAR_SOLVER_H 23 #include "libmesh/libmesh_common.h" 25 #ifdef LIBMESH_HAVE_EIGEN 30 #include "libmesh/linear_solver.h" 31 #include "libmesh/eigen_sparse_vector.h" 32 #include "libmesh/eigen_sparse_matrix.h" 33 #include "libmesh/enum_convergence_flags.h" 48 class EigenSparseLinearSolver :
public LinearSolver<T>
64 virtual void clear ()
override;
69 virtual void init (
const char *
name=
nullptr)
override;
74 virtual std::pair<unsigned int, Real>
75 solve (SparseMatrix<T> & matrix,
76 NumericVector<T> & solution,
77 NumericVector<T> & rhs,
78 const std::optional<double> tol = std::nullopt,
79 const std::optional<unsigned int> m_its = std::nullopt)
override;
84 virtual std::pair<unsigned int, Real>
86 NumericVector<T> & solution,
87 NumericVector<T> & rhs,
88 const std::optional<double> tol = std::nullopt,
89 const std::optional<unsigned int> m_its = std::nullopt)
override;
94 virtual std::pair<unsigned int, Real>
95 solve (SparseMatrix<T> & matrix,
97 NumericVector<T> & solution,
98 NumericVector<T> & rhs,
99 const std::optional<double> tol = std::nullopt,
100 const std::optional<unsigned int> m_its = std::nullopt)
override;
105 virtual std::pair<unsigned int, Real>
106 solve (
const ShellMatrix<T> & shell_matrix,
107 NumericVector<T> & solution_in,
108 NumericVector<T> & rhs_in,
109 const std::optional<double> tol = std::nullopt,
110 const std::optional<unsigned int> m_its = std::nullopt)
override;
117 virtual std::pair<unsigned int, Real>
118 solve (
const ShellMatrix<T> & shell_matrix,
119 const SparseMatrix<T> & precond_matrix,
120 NumericVector<T> & solution_in,
121 NumericVector<T> & rhs_in,
122 const std::optional<double> tol = std::nullopt,
123 const std::optional<unsigned int> m_its = std::nullopt)
override;
152 static std::map<Eigen::ComputationInfo, LinearConvergenceReason>
build_map()
154 std::map<Eigen::ComputationInfo, LinearConvergenceReason> ret;
168 template <
typename T>
169 std::map<Eigen::ComputationInfo, LinearConvergenceReason>
174 template <
typename T>
183 template <
typename T>
185 std::pair<unsigned int, Real>
190 const std::optional<double>,
191 const std::optional<unsigned int>)
193 libmesh_error_msg(
"ERROR: Eigen does not support a user-supplied preconditioner!");
195 return std::pair<unsigned int, Real>();
200 #endif // #ifdef LIBMESH_HAVE_EIGEN 201 #endif // LIBMESH_EIGEN_SPARSE_LINEAR_SOLVER_H std::string name(const ElemQuality q)
This function returns a string containing some name for q.
~EigenSparseLinearSolver()
Destructor.
Eigen::ComputationInfo _comp_info
Store the result of the last solve.
EigenSparseLinearSolver(const libMesh::Parallel::Communicator &comm_in)
Constructor.
LinearConvergenceReason
Linear solver convergence flags (taken from the PETSc flags).
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
virtual void init(const char *name=nullptr) override
Initialize data structures if not done so already.
The libMesh namespace provides an interface to certain functionality in the library.
static std::map< Eigen::ComputationInfo, LinearConvergenceReason > build_map()
Static function used to initialize _convergence_reasons map.
virtual LinearConvergenceReason get_converged_reason() const override
static std::map< Eigen::ComputationInfo, LinearConvergenceReason > _convergence_reasons
Static map between Eigen ComputationInfo enumerations and libMesh LinearConvergenceReason enumeration...
virtual std::pair< unsigned int, Real > solve(SparseMatrix< T > &matrix, NumericVector< T > &solution, NumericVector< T > &rhs, const std::optional< double > tol=std::nullopt, const std::optional< unsigned int > m_its=std::nullopt) override
Call the Eigen solver.
virtual void clear() override
Release all memory and clear data structures.
void set_eigen_preconditioner_type()
Tells Eigen to use the user-specified preconditioner stored in _preconditioner_type.
virtual std::pair< unsigned int, Real > adjoint_solve(SparseMatrix< T > &matrix, NumericVector< T > &solution, NumericVector< T > &rhs, const std::optional< double > tol=std::nullopt, const std::optional< unsigned int > m_its=std::nullopt) override
Call the Eigen solver to solve A^T x = b.