Storm 1.11.1.1
A Modern Probabilistic Model Checker
Loading...
Searching...
No Matches
EigenLinearEquationSolver.cpp
Go to the documentation of this file.
2
8
9namespace storm {
10namespace solver {
11
12template<typename ValueType>
16
17template<typename ValueType>
21
22template<typename ValueType>
26
27template<typename ValueType>
29 eigenA = storm::adapters::EigenAdapter::toEigenSparseMatrix<ValueType>(A);
30 this->clearCache();
31}
32
33template<typename ValueType>
35 // Take ownership of the matrix so it is destroyed after we have translated it to Eigen's format.
36 storm::storage::SparseMatrix<ValueType> localA(std::move(A));
37 this->setMatrix(localA);
38 this->clearCache();
39}
40
41template<typename ValueType>
42EigenLinearEquationSolverMethod EigenLinearEquationSolver<ValueType>::getMethod(Environment const& env, bool isExactMode) const {
43 // Adjust the method if none was specified and we are using rational numbers.
44 auto method = env.solver().eigen().getMethod();
45
46 if (isExactMode && method != EigenLinearEquationSolverMethod::SparseLU) {
47 if (env.solver().eigen().isMethodSetFromDefault()) {
48 STORM_LOG_INFO("Selecting 'SparseLU' as the solution technique to guarantee exact results.");
49 } else {
50 STORM_LOG_WARN("The selected solution method does not guarantee exact results. Falling back to SparseLU.");
51 }
52 method = EigenLinearEquationSolverMethod::SparseLU;
53 } else if (env.solver().isForceSoundness() && method != EigenLinearEquationSolverMethod::SparseLU) {
54 if (env.solver().eigen().isMethodSetFromDefault()) {
56 "Selecting 'SparseLU' as the solution technique to guarantee sound results. If you want to override this, please explicitly specify a "
57 "different method.");
58 method = EigenLinearEquationSolverMethod::SparseLU;
59 } else {
60 STORM_LOG_WARN("The selected solution method does not guarantee sound results.");
61 }
62 }
63 return method;
64}
65
66#ifdef STORM_HAVE_CARL
67// Specialization for storm::RationalNumber
68template<>
69bool EigenLinearEquationSolver<storm::RationalNumber>::internalSolveEquations(Environment const& env, std::vector<storm::RationalNumber>& x,
70 std::vector<storm::RationalNumber> const& b) const {
71 auto solutionMethod = getMethod(env, true);
72 STORM_LOG_WARN_COND(solutionMethod == EigenLinearEquationSolverMethod::SparseLU, "Switching method to SparseLU.");
73 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with with rational numbers using LU factorization (Eigen library).");
74
75 // Map the input vectors to Eigen's format.
76 auto eigenX = Eigen::Matrix<storm::RationalNumber, Eigen::Dynamic, 1>::Map(x.data(), x.size());
77 auto eigenB = Eigen::Matrix<storm::RationalNumber, Eigen::Dynamic, 1>::Map(b.data(), b.size());
78
79 Eigen::SparseLU<Eigen::SparseMatrix<storm::RationalNumber>, Eigen::COLAMDOrdering<int>> solver;
80 solver.compute(*eigenA);
81 solver._solve_impl(eigenB, eigenX);
82
83 return solver.info() == Eigen::ComputationInfo::Success;
84}
85
86// Specialization for storm::RationalFunction
87template<>
88bool EigenLinearEquationSolver<storm::RationalFunction>::internalSolveEquations(Environment const& env, std::vector<storm::RationalFunction>& x,
89 std::vector<storm::RationalFunction> const& b) const {
90 auto solutionMethod = getMethod(env, true);
91 STORM_LOG_WARN_COND(solutionMethod == EigenLinearEquationSolverMethod::SparseLU, "Switching method to SparseLU.");
92 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with rational functions using LU factorization (Eigen library).");
93
94 // Map the input vectors to Eigen's format.
95 auto eigenX = Eigen::Matrix<storm::RationalFunction, Eigen::Dynamic, 1>::Map(x.data(), x.size());
96 auto eigenB = Eigen::Matrix<storm::RationalFunction, Eigen::Dynamic, 1>::Map(b.data(), b.size());
97
98 Eigen::SparseLU<Eigen::SparseMatrix<storm::RationalFunction>, Eigen::COLAMDOrdering<int>> solver;
99 solver.compute(*eigenA);
100 solver._solve_impl(eigenB, eigenX);
101 return solver.info() == Eigen::ComputationInfo::Success;
102}
103#endif
104
105template<typename ValueType>
106bool EigenLinearEquationSolver<ValueType>::internalSolveEquations(Environment const& env, std::vector<ValueType>& x, std::vector<ValueType> const& b) const {
107 // Map the input vectors to Eigen's format.
108 auto eigenX = Eigen::Matrix<ValueType, Eigen::Dynamic, 1>::Map(x.data(), x.size());
109 auto eigenB = Eigen::Matrix<ValueType, Eigen::Dynamic, 1>::Map(b.data(), b.size());
110
111 auto solutionMethod = getMethod(env, env.solver().isForceExact());
112 if (solutionMethod == EigenLinearEquationSolverMethod::SparseLU) {
113 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with sparse LU factorization (Eigen library).");
114 Eigen::SparseLU<Eigen::SparseMatrix<ValueType>, Eigen::COLAMDOrdering<int>> solver;
115 solver.compute(*this->eigenA);
116 solver._solve_impl(eigenB, eigenX);
117 } else {
118 bool converged = false;
119 uint64_t numberOfIterations = 0;
120 Eigen::Index maxIter = std::numeric_limits<Eigen::Index>::max();
121 if (env.solver().eigen().getMaximalNumberOfIterations() < static_cast<uint64_t>(maxIter)) {
122 maxIter = env.solver().eigen().getMaximalNumberOfIterations();
123 }
124 uint64_t restartThreshold = env.solver().eigen().getRestartThreshold();
125 ValueType precision = storm::utility::convertNumber<ValueType>(env.solver().eigen().getPrecision());
126 EigenLinearEquationSolverPreconditioner preconditioner = env.solver().eigen().getPreconditioner();
127 if (solutionMethod == EigenLinearEquationSolverMethod::Bicgstab) {
128 if (preconditioner == EigenLinearEquationSolverPreconditioner::Ilu) {
129 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with BiCGSTAB with Ilu preconditioner (Eigen library).");
130
131 Eigen::BiCGSTAB<Eigen::SparseMatrix<ValueType>, Eigen::IncompleteLUT<ValueType>> solver;
132 solver.compute(*this->eigenA);
133 solver.setTolerance(precision);
134 solver.setMaxIterations(maxIter);
135 eigenX = solver.solveWithGuess(eigenB, eigenX);
136 converged = solver.info() == Eigen::ComputationInfo::Success;
137 numberOfIterations = solver.iterations();
138 } else if (preconditioner == EigenLinearEquationSolverPreconditioner::Diagonal) {
139 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with BiCGSTAB with Diagonal preconditioner (Eigen library).");
140
141 Eigen::BiCGSTAB<Eigen::SparseMatrix<ValueType>, Eigen::DiagonalPreconditioner<ValueType>> solver;
142 solver.setTolerance(precision);
143 solver.setMaxIterations(maxIter);
144 solver.compute(*this->eigenA);
145 eigenX = solver.solveWithGuess(eigenB, eigenX);
146 converged = solver.info() == Eigen::ComputationInfo::Success;
147 numberOfIterations = solver.iterations();
148 } else {
149 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with BiCGSTAB with identity preconditioner (Eigen library).");
150
151 Eigen::BiCGSTAB<Eigen::SparseMatrix<ValueType>, Eigen::IdentityPreconditioner> solver;
152 solver.setTolerance(precision);
153 solver.setMaxIterations(maxIter);
154 solver.compute(*this->eigenA);
155 eigenX = solver.solveWithGuess(eigenB, eigenX);
156 numberOfIterations = solver.iterations();
157 converged = solver.info() == Eigen::ComputationInfo::Success;
158 }
159 } else if (solutionMethod == EigenLinearEquationSolverMethod::DGmres) {
160 if (preconditioner == EigenLinearEquationSolverPreconditioner::Ilu) {
161 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with DGMRES with Ilu preconditioner (Eigen library).");
162 Eigen::DGMRES<Eigen::SparseMatrix<ValueType>, Eigen::IncompleteLUT<ValueType>> solver;
163 solver.setTolerance(precision);
164 solver.setMaxIterations(maxIter);
165 solver.set_restart(restartThreshold);
166 solver.compute(*this->eigenA);
167 eigenX = solver.solveWithGuess(eigenB, eigenX);
168 converged = solver.info() == Eigen::ComputationInfo::Success;
169 numberOfIterations = solver.iterations();
170 } else if (preconditioner == EigenLinearEquationSolverPreconditioner::Diagonal) {
171 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with DGMRES with Diagonal preconditioner (Eigen library).");
172
173 Eigen::DGMRES<Eigen::SparseMatrix<ValueType>, Eigen::DiagonalPreconditioner<ValueType>> solver;
174 solver.setTolerance(precision);
175 solver.setMaxIterations(maxIter);
176 solver.set_restart(restartThreshold);
177 solver.compute(*this->eigenA);
178 eigenX = solver.solveWithGuess(eigenB, eigenX);
179 converged = solver.info() == Eigen::ComputationInfo::Success;
180 numberOfIterations = solver.iterations();
181 } else {
182 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with DGMRES with identity preconditioner (Eigen library).");
183
184 Eigen::DGMRES<Eigen::SparseMatrix<ValueType>, Eigen::IdentityPreconditioner> solver;
185 solver.setTolerance(precision);
186 solver.setMaxIterations(maxIter);
187 solver.set_restart(restartThreshold);
188 solver.compute(*this->eigenA);
189 eigenX = solver.solveWithGuess(eigenB, eigenX);
190 converged = solver.info() == Eigen::ComputationInfo::Success;
191 numberOfIterations = solver.iterations();
192 }
193
194 } else if (solutionMethod == EigenLinearEquationSolverMethod::Gmres) {
195 if (preconditioner == EigenLinearEquationSolverPreconditioner::Ilu) {
196 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with GMRES with Ilu preconditioner (Eigen library).");
197
198 Eigen::GMRES<Eigen::SparseMatrix<ValueType>, Eigen::IncompleteLUT<ValueType>> solver;
199 solver.setTolerance(precision);
200 solver.setMaxIterations(maxIter);
201 solver.set_restart(restartThreshold);
202 solver.compute(*this->eigenA);
203 eigenX = solver.solveWithGuess(eigenB, eigenX);
204 converged = solver.info() == Eigen::ComputationInfo::Success;
205 numberOfIterations = solver.iterations();
206 } else if (preconditioner == EigenLinearEquationSolverPreconditioner::Diagonal) {
207 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with GMRES with Diagonal preconditioner (Eigen library).");
208
209 Eigen::GMRES<Eigen::SparseMatrix<ValueType>, Eigen::DiagonalPreconditioner<ValueType>> solver;
210 solver.setTolerance(precision);
211 solver.setMaxIterations(maxIter);
212 solver.set_restart(restartThreshold);
213 solver.compute(*this->eigenA);
214 eigenX = solver.solveWithGuess(eigenB, eigenX);
215 converged = solver.info() == Eigen::ComputationInfo::Success;
216 numberOfIterations = solver.iterations();
217 } else {
218 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with GMRES with identity preconditioner (Eigen library).");
219
220 Eigen::GMRES<Eigen::SparseMatrix<ValueType>, Eigen::IdentityPreconditioner> solver;
221 solver.setTolerance(precision);
222 solver.setMaxIterations(maxIter);
223 solver.set_restart(restartThreshold);
224 solver.compute(*this->eigenA);
225 eigenX = solver.solveWithGuess(eigenB, eigenX);
226 converged = solver.info() == Eigen::ComputationInfo::Success;
227 numberOfIterations = solver.iterations();
228 }
229 }
230
231 // Make sure that all results conform to the (global) bounds.
232 storm::utility::vector::clip(x, this->lowerBound, this->upperBound);
233
234 // Check if the solver converged and issue a warning otherwise.
235 if (converged) {
236 STORM_LOG_INFO("Iterative solver converged after " << numberOfIterations << " iterations.");
237 return true;
238 } else {
239 STORM_LOG_WARN("Iterative solver did not converge.");
240 return false;
241 }
242 }
243
244 return true;
245}
246
247template<typename ValueType>
251
252template<typename ValueType>
254 return eigenA->rows();
255}
256
257template<typename ValueType>
258uint64_t EigenLinearEquationSolver<ValueType>::getMatrixColumnCount() const {
259 return eigenA->cols();
260}
261
262template<typename ValueType>
263std::unique_ptr<storm::solver::LinearEquationSolver<ValueType>> EigenLinearEquationSolverFactory<ValueType>::create(Environment const&) const {
264 return std::make_unique<storm::solver::EigenLinearEquationSolver<ValueType>>();
265}
266
267template<typename ValueType>
268std::unique_ptr<LinearEquationSolverFactory<ValueType>> EigenLinearEquationSolverFactory<ValueType>::clone() const {
269 return std::make_unique<EigenLinearEquationSolverFactory<ValueType>>(*this);
270}
271
274
277
280} // namespace solver
281} // namespace storm
storm::solver::EigenLinearEquationSolverMethod const & getMethod() const
uint64_t const & getMaximalNumberOfIterations() const
uint64_t const & getRestartThreshold() const
storm::RationalNumber const & getPrecision() const
storm::solver::EigenLinearEquationSolverPreconditioner const & getPreconditioner() const
SolverEnvironment & solver()
EigenSolverEnvironment & eigen()
virtual std::unique_ptr< LinearEquationSolverFactory< ValueType > > clone() const override
Creates a copy of this factory.
virtual std::unique_ptr< storm::solver::LinearEquationSolver< ValueType > > create(Environment const &env) const override
Creates an equation solver with the current settings, but without a matrix.
A class that uses the Eigen library to implement the LinearEquationSolver interface.
virtual void setMatrix(storm::storage::SparseMatrix< ValueType > const &A) override
virtual bool internalSolveEquations(Environment const &env, std::vector< ValueType > &x, std::vector< ValueType > const &b) const override
virtual LinearEquationSolverProblemFormat getEquationProblemFormat(Environment const &env) const override
Retrieves the format in which this solver expects to solve equations.
A class that holds a possibly non-square matrix in the compressed row storage format.
#define STORM_LOG_INFO(message)
Definition logging.h:24
#define STORM_LOG_WARN(message)
Definition logging.h:25
#define STORM_LOG_WARN_COND(cond, message)
Definition macros.h:38
void clip(std::vector< ValueType > &x, boost::optional< ValueType > const &lowerBound, boost::optional< ValueType > const &upperBound)
Takes the input vector and ensures that all entries conform to the bounds.
Definition vector.h:888