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// Specialization for storm::RationalNumber
67template<>
69 std::vector<storm::RationalNumber> const& b) const {
70 auto solutionMethod = getMethod(env, true);
71 STORM_LOG_WARN_COND(solutionMethod == EigenLinearEquationSolverMethod::SparseLU, "Switching method to SparseLU.");
72 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with with rational numbers using LU factorization (Eigen library).");
73
74 // Map the input vectors to Eigen's format.
75 auto eigenX = Eigen::Matrix<storm::RationalNumber, Eigen::Dynamic, 1>::Map(x.data(), x.size());
76 auto eigenB = Eigen::Matrix<storm::RationalNumber, Eigen::Dynamic, 1>::Map(b.data(), b.size());
77
78 Eigen::SparseLU<Eigen::SparseMatrix<storm::RationalNumber>, Eigen::COLAMDOrdering<int>> solver;
79 solver.compute(*eigenA);
80 solver._solve_impl(eigenB, eigenX);
81
82 return solver.info() == Eigen::ComputationInfo::Success;
83}
84
85// Specialization for storm::RationalFunction
86template<>
87bool EigenLinearEquationSolver<storm::RationalFunction>::internalSolveEquations(Environment const& env, std::vector<storm::RationalFunction>& x,
88 std::vector<storm::RationalFunction> const& b) const {
89 auto solutionMethod = getMethod(env, true);
90 STORM_LOG_WARN_COND(solutionMethod == EigenLinearEquationSolverMethod::SparseLU, "Switching method to SparseLU.");
91 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with rational functions using LU factorization (Eigen library).");
92
93 // Map the input vectors to Eigen's format.
94 auto eigenX = Eigen::Matrix<storm::RationalFunction, Eigen::Dynamic, 1>::Map(x.data(), x.size());
95 auto eigenB = Eigen::Matrix<storm::RationalFunction, Eigen::Dynamic, 1>::Map(b.data(), b.size());
96
97 Eigen::SparseLU<Eigen::SparseMatrix<storm::RationalFunction>, Eigen::COLAMDOrdering<int>> solver;
98 solver.compute(*eigenA);
99 solver._solve_impl(eigenB, eigenX);
100 return solver.info() == Eigen::ComputationInfo::Success;
101}
102
103template<typename ValueType>
104bool EigenLinearEquationSolver<ValueType>::internalSolveEquations(Environment const& env, std::vector<ValueType>& x, std::vector<ValueType> const& b) const {
105 // Map the input vectors to Eigen's format.
106 auto eigenX = Eigen::Matrix<ValueType, Eigen::Dynamic, 1>::Map(x.data(), x.size());
107 auto eigenB = Eigen::Matrix<ValueType, Eigen::Dynamic, 1>::Map(b.data(), b.size());
108
109 auto solutionMethod = getMethod(env, env.solver().isForceExact());
110 if (solutionMethod == EigenLinearEquationSolverMethod::SparseLU) {
111 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with sparse LU factorization (Eigen library).");
112 Eigen::SparseLU<Eigen::SparseMatrix<ValueType>, Eigen::COLAMDOrdering<int>> solver;
113 solver.compute(*this->eigenA);
114 solver._solve_impl(eigenB, eigenX);
115 } else {
116 bool converged = false;
117 uint64_t numberOfIterations = 0;
118 Eigen::Index maxIter = std::numeric_limits<Eigen::Index>::max();
119 if (env.solver().eigen().getMaximalNumberOfIterations() < static_cast<uint64_t>(maxIter)) {
120 maxIter = env.solver().eigen().getMaximalNumberOfIterations();
121 }
122 uint64_t restartThreshold = env.solver().eigen().getRestartThreshold();
123 ValueType precision = storm::utility::convertNumber<ValueType>(env.solver().eigen().getPrecision());
124 EigenLinearEquationSolverPreconditioner preconditioner = env.solver().eigen().getPreconditioner();
125 if (solutionMethod == EigenLinearEquationSolverMethod::Bicgstab) {
126 if (preconditioner == EigenLinearEquationSolverPreconditioner::Ilu) {
127 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with BiCGSTAB with Ilu preconditioner (Eigen library).");
128
129 Eigen::BiCGSTAB<Eigen::SparseMatrix<ValueType>, Eigen::IncompleteLUT<ValueType>> solver;
130 solver.compute(*this->eigenA);
131 solver.setTolerance(precision);
132 solver.setMaxIterations(maxIter);
133 eigenX = solver.solveWithGuess(eigenB, eigenX);
134 converged = solver.info() == Eigen::ComputationInfo::Success;
135 numberOfIterations = solver.iterations();
136 } else if (preconditioner == EigenLinearEquationSolverPreconditioner::Diagonal) {
137 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with BiCGSTAB with Diagonal preconditioner (Eigen library).");
138
139 Eigen::BiCGSTAB<Eigen::SparseMatrix<ValueType>, Eigen::DiagonalPreconditioner<ValueType>> solver;
140 solver.setTolerance(precision);
141 solver.setMaxIterations(maxIter);
142 solver.compute(*this->eigenA);
143 eigenX = solver.solveWithGuess(eigenB, eigenX);
144 converged = solver.info() == Eigen::ComputationInfo::Success;
145 numberOfIterations = solver.iterations();
146 } else {
147 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with BiCGSTAB with identity preconditioner (Eigen library).");
148
149 Eigen::BiCGSTAB<Eigen::SparseMatrix<ValueType>, Eigen::IdentityPreconditioner> solver;
150 solver.setTolerance(precision);
151 solver.setMaxIterations(maxIter);
152 solver.compute(*this->eigenA);
153 eigenX = solver.solveWithGuess(eigenB, eigenX);
154 numberOfIterations = solver.iterations();
155 converged = solver.info() == Eigen::ComputationInfo::Success;
156 }
157 } else if (solutionMethod == EigenLinearEquationSolverMethod::DGmres) {
158 if (preconditioner == EigenLinearEquationSolverPreconditioner::Ilu) {
159 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with DGMRES with Ilu preconditioner (Eigen library).");
160 Eigen::DGMRES<Eigen::SparseMatrix<ValueType>, Eigen::IncompleteLUT<ValueType>> solver;
161 solver.setTolerance(precision);
162 solver.setMaxIterations(maxIter);
163 solver.set_restart(restartThreshold);
164 solver.compute(*this->eigenA);
165 eigenX = solver.solveWithGuess(eigenB, eigenX);
166 converged = solver.info() == Eigen::ComputationInfo::Success;
167 numberOfIterations = solver.iterations();
168 } else if (preconditioner == EigenLinearEquationSolverPreconditioner::Diagonal) {
169 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with DGMRES with Diagonal preconditioner (Eigen library).");
170
171 Eigen::DGMRES<Eigen::SparseMatrix<ValueType>, Eigen::DiagonalPreconditioner<ValueType>> solver;
172 solver.setTolerance(precision);
173 solver.setMaxIterations(maxIter);
174 solver.set_restart(restartThreshold);
175 solver.compute(*this->eigenA);
176 eigenX = solver.solveWithGuess(eigenB, eigenX);
177 converged = solver.info() == Eigen::ComputationInfo::Success;
178 numberOfIterations = solver.iterations();
179 } else {
180 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with DGMRES with identity preconditioner (Eigen library).");
181
182 Eigen::DGMRES<Eigen::SparseMatrix<ValueType>, Eigen::IdentityPreconditioner> solver;
183 solver.setTolerance(precision);
184 solver.setMaxIterations(maxIter);
185 solver.set_restart(restartThreshold);
186 solver.compute(*this->eigenA);
187 eigenX = solver.solveWithGuess(eigenB, eigenX);
188 converged = solver.info() == Eigen::ComputationInfo::Success;
189 numberOfIterations = solver.iterations();
190 }
191
192 } else if (solutionMethod == EigenLinearEquationSolverMethod::Gmres) {
193 if (preconditioner == EigenLinearEquationSolverPreconditioner::Ilu) {
194 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with GMRES with Ilu preconditioner (Eigen library).");
195
196 Eigen::GMRES<Eigen::SparseMatrix<ValueType>, Eigen::IncompleteLUT<ValueType>> solver;
197 solver.setTolerance(precision);
198 solver.setMaxIterations(maxIter);
199 solver.set_restart(restartThreshold);
200 solver.compute(*this->eigenA);
201 eigenX = solver.solveWithGuess(eigenB, eigenX);
202 converged = solver.info() == Eigen::ComputationInfo::Success;
203 numberOfIterations = solver.iterations();
204 } else if (preconditioner == EigenLinearEquationSolverPreconditioner::Diagonal) {
205 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with GMRES with Diagonal preconditioner (Eigen library).");
206
207 Eigen::GMRES<Eigen::SparseMatrix<ValueType>, Eigen::DiagonalPreconditioner<ValueType>> solver;
208 solver.setTolerance(precision);
209 solver.setMaxIterations(maxIter);
210 solver.set_restart(restartThreshold);
211 solver.compute(*this->eigenA);
212 eigenX = solver.solveWithGuess(eigenB, eigenX);
213 converged = solver.info() == Eigen::ComputationInfo::Success;
214 numberOfIterations = solver.iterations();
215 } else {
216 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with GMRES with identity preconditioner (Eigen library).");
217
218 Eigen::GMRES<Eigen::SparseMatrix<ValueType>, Eigen::IdentityPreconditioner> solver;
219 solver.setTolerance(precision);
220 solver.setMaxIterations(maxIter);
221 solver.set_restart(restartThreshold);
222 solver.compute(*this->eigenA);
223 eigenX = solver.solveWithGuess(eigenB, eigenX);
224 converged = solver.info() == Eigen::ComputationInfo::Success;
225 numberOfIterations = solver.iterations();
226 }
227 }
228
229 // Make sure that all results conform to the (global) bounds.
230 storm::utility::vector::clip(x, this->lowerBound, this->upperBound);
231
232 // Check if the solver converged and issue a warning otherwise.
233 if (converged) {
234 STORM_LOG_INFO("Iterative solver converged after " << numberOfIterations << " iterations.");
235 return true;
236 } else {
237 STORM_LOG_WARN("Iterative solver did not converge.");
238 return false;
239 }
240 }
241
242 return true;
243}
244
245template<typename ValueType>
249
250template<typename ValueType>
252 return eigenA->rows();
253}
254
255template<typename ValueType>
256uint64_t EigenLinearEquationSolver<ValueType>::getMatrixColumnCount() const {
257 return eigenA->cols();
258}
259
260template<typename ValueType>
261std::unique_ptr<storm::solver::LinearEquationSolver<ValueType>> EigenLinearEquationSolverFactory<ValueType>::create(Environment const&) const {
262 return std::make_unique<storm::solver::EigenLinearEquationSolver<ValueType>>();
263}
264
265template<typename ValueType>
266std::unique_ptr<LinearEquationSolverFactory<ValueType>> EigenLinearEquationSolverFactory<ValueType>::clone() const {
267 return std::make_unique<EigenLinearEquationSolverFactory<ValueType>>(*this);
268}
269
272
275
278} // namespace solver
279} // 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