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