Storm
A Modern Probabilistic Model Checker
Loading...
Searching...
No Matches
NativeLinearEquationSolver.cpp
Go to the documentation of this file.
2
3#include <limits>
4
7
21
22namespace storm {
23namespace solver {
24
25template<typename ValueType>
27 // Intentionally left empty.
28}
29
30template<typename ValueType>
34
35template<typename ValueType>
39
40template<typename ValueType>
42 localA.reset();
43 this->A = &A;
44 clearCache();
45}
46
47template<typename ValueType>
49 localA = std::make_unique<storm::storage::SparseMatrix<ValueType>>(std::move(A));
50 this->A = localA.get();
51 clearCache();
52}
53
54template<typename ValueType>
56 if (!viOperator) {
57 viOperator = std::make_shared<helper::ValueIterationOperator<ValueType, true>>();
58 viOperator->setMatrixBackwards(*this->A);
59 }
60}
61
62template<typename ValueType>
63bool NativeLinearEquationSolver<ValueType>::solveEquationsSOR(Environment const& env, std::vector<ValueType>& x, std::vector<ValueType> const& b,
64 ValueType const& omega) const {
65 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with NativeLinearEquationSolver (Gauss-Seidel, SOR omega = " << omega << ")");
66
67 if (!this->cachedRowVector) {
68 this->cachedRowVector = std::make_unique<std::vector<ValueType>>(getMatrixRowCount());
69 }
70
71 ValueType precision = storm::utility::convertNumber<ValueType>(env.solver().native().getPrecision());
72 uint64_t maxIter = env.solver().native().getMaximalNumberOfIterations();
73 bool relative = env.solver().native().getRelativeTerminationCriterion();
74
75 // Set up additional environment variables.
76 uint_fast64_t iterations = 0;
78
79 this->startMeasureProgress();
80 while (status == SolverStatus::InProgress && iterations < maxIter) {
81 A->performSuccessiveOverRelaxationStep(omega, x, b);
82
83 // Now check if the process already converged within our precision.
84 if (storm::utility::vector::equalModuloPrecision<ValueType>(*this->cachedRowVector, x, precision, relative)) {
86 }
87 // If we did not yet converge, we need to backup the contents of x.
88 if (status != SolverStatus::Converged) {
89 *this->cachedRowVector = x;
90 }
91
92 // Potentially show progress.
93 this->showProgressIterative(iterations);
94
95 // Increase iteration count so we can abort if convergence is too slow.
96 ++iterations;
97
98 status = this->updateStatus(status, x, SolverGuarantee::None, iterations, maxIter);
99 }
100
101 if (!this->isCachingEnabled()) {
102 clearCache();
103 }
104
105 this->reportStatus(status, iterations);
106
107 return status == SolverStatus::Converged;
108}
109
110template<typename ValueType>
111NativeLinearEquationSolver<ValueType>::JacobiDecomposition::JacobiDecomposition(Environment const& env, storm::storage::SparseMatrix<ValueType> const& A) {
112 auto decomposition = A.getJacobiDecomposition();
113 this->LUMatrix = std::move(decomposition.first);
114 this->DVector = std::move(decomposition.second);
115 this->multiplier = storm::solver::MultiplierFactory<ValueType>().create(env, this->LUMatrix);
116}
117
118template<typename ValueType>
119bool NativeLinearEquationSolver<ValueType>::solveEquationsJacobi(Environment const& env, std::vector<ValueType>& x, std::vector<ValueType> const& b) const {
120 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with NativeLinearEquationSolver (Jacobi)");
121
122 if (!this->cachedRowVector) {
123 this->cachedRowVector = std::make_unique<std::vector<ValueType>>(getMatrixRowCount());
124 }
125
126 // Get a Jacobi decomposition of the matrix A.
127 if (!jacobiDecomposition) {
128 jacobiDecomposition = std::make_unique<JacobiDecomposition>(env, *A);
129 }
130
131 ValueType precision = storm::utility::convertNumber<ValueType>(env.solver().native().getPrecision());
132 uint64_t maxIter = env.solver().native().getMaximalNumberOfIterations();
133 bool relative = env.solver().native().getRelativeTerminationCriterion();
134
135 std::vector<ValueType>* currentX = &x;
136 std::vector<ValueType>* nextX = this->cachedRowVector.get();
137
138 // Set up additional environment variables.
139 uint_fast64_t iterations = 0;
141
142 this->startMeasureProgress();
143 while (status == SolverStatus::InProgress && iterations < maxIter) {
144 // Compute D^-1 * (b - LU * x) and store result in nextX.
145 jacobiDecomposition->multiplier->multiply(env, *currentX, nullptr, *nextX);
147 storm::utility::vector::multiplyVectorsPointwise(jacobiDecomposition->DVector, *nextX, *nextX);
148
149 // Now check if the process already converged within our precision.
150 if (storm::utility::vector::equalModuloPrecision<ValueType>(*currentX, *nextX, precision, relative)) {
152 }
153 // Swap the two pointers as a preparation for the next iteration.
154 std::swap(nextX, currentX);
155
156 // Potentially show progress.
157 this->showProgressIterative(iterations);
158
159 // Increase iteration count so we can abort if convergence is too slow.
160 ++iterations;
161
162 status = this->updateStatus(status, *currentX, SolverGuarantee::None, iterations, maxIter);
163 }
164
165 // If the last iteration did not write to the original x we have to swap the contents, because the
166 // output has to be written to the input parameter x.
167 if (currentX == this->cachedRowVector.get()) {
168 std::swap(x, *currentX);
169 }
170
171 if (!this->isCachingEnabled()) {
172 clearCache();
173 }
174
175 this->reportStatus(status, iterations);
176
177 return status == SolverStatus::Converged;
178}
179
180template<typename ValueType>
181NativeLinearEquationSolver<ValueType>::WalkerChaeData::WalkerChaeData(Environment const& env, storm::storage::SparseMatrix<ValueType> const& originalMatrix,
182 std::vector<ValueType> const& originalB)
183 : t(storm::utility::convertNumber<ValueType>(1000.0)) {
184 computeWalkerChaeMatrix(originalMatrix);
185 computeNewB(originalB);
186 precomputeAuxiliaryData();
187 multiplier = storm::solver::MultiplierFactory<ValueType>().create(env, this->matrix);
188}
189
190template<typename ValueType>
191void NativeLinearEquationSolver<ValueType>::WalkerChaeData::computeWalkerChaeMatrix(storm::storage::SparseMatrix<ValueType> const& originalMatrix) {
192 storm::storage::BitVector columnsWithNegativeEntries(originalMatrix.getColumnCount());
193 ValueType zero = storm::utility::zero<ValueType>();
194 for (auto const& e : originalMatrix) {
195 if (e.getValue() < zero) {
196 columnsWithNegativeEntries.set(e.getColumn());
197 }
198 }
199 std::vector<uint64_t> columnsWithNegativeEntriesBefore = columnsWithNegativeEntries.getNumberOfSetBitsBeforeIndices();
200
201 // We now build an extended equation system matrix that only has non-negative coefficients.
203
204 uint64_t row = 0;
205 for (; row < originalMatrix.getRowCount(); ++row) {
206 for (auto const& entry : originalMatrix.getRow(row)) {
207 if (entry.getValue() < zero) {
208 builder.addNextValue(row, originalMatrix.getRowCount() + columnsWithNegativeEntriesBefore[entry.getColumn()], -entry.getValue());
209 } else {
210 builder.addNextValue(row, entry.getColumn(), entry.getValue());
211 }
212 }
213 }
214 ValueType one = storm::utility::one<ValueType>();
215 for (auto column : columnsWithNegativeEntries) {
216 builder.addNextValue(row, column, one);
217 builder.addNextValue(row, originalMatrix.getRowCount() + columnsWithNegativeEntriesBefore[column], one);
218 ++row;
219 }
220
221 matrix = builder.build();
222}
223
224template<typename ValueType>
225void NativeLinearEquationSolver<ValueType>::WalkerChaeData::computeNewB(std::vector<ValueType> const& originalB) {
226 b = std::vector<ValueType>(originalB);
227 b.resize(matrix.getRowCount());
228}
229
230template<typename ValueType>
231void NativeLinearEquationSolver<ValueType>::WalkerChaeData::precomputeAuxiliaryData() {
232 columnSums = std::vector<ValueType>(matrix.getColumnCount());
233 for (auto const& e : matrix) {
234 columnSums[e.getColumn()] += e.getValue();
235 }
236
237 newX.resize(matrix.getRowCount());
238}
239
240template<typename ValueType>
241bool NativeLinearEquationSolver<ValueType>::solveEquationsWalkerChae(Environment const& env, std::vector<ValueType>& x, std::vector<ValueType> const& b) const {
242 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with NativeLinearEquationSolver (WalkerChae)");
243
244 // (1) Compute an equivalent equation system that has only non-negative coefficients.
245 if (!walkerChaeData) {
246 walkerChaeData = std::make_unique<WalkerChaeData>(env, *this->A, b);
247 }
248
249 // (2) Enlarge the vectors x and b to account for additional variables.
250 x.resize(walkerChaeData->matrix.getRowCount());
251
252 // Square the error bound, so we can use it to check for convergence. We take the squared error, because we
253 // do not want to compute the root in the 2-norm computation.
254 ValueType squaredErrorBound = storm::utility::pow(storm::utility::convertNumber<ValueType>(env.solver().native().getPrecision()), 2);
255
256 uint64_t maxIter = env.solver().native().getMaximalNumberOfIterations();
257
258 // Set up references to the x-vectors used in the iteration loop.
259 std::vector<ValueType>* currentX = &x;
260 std::vector<ValueType>* nextX = &walkerChaeData->newX;
261
262 std::vector<ValueType> tmp = walkerChaeData->matrix.getRowSumVector();
263 storm::utility::vector::applyPointwise(tmp, walkerChaeData->b, walkerChaeData->b,
264 [this](ValueType const& first, ValueType const& second) -> ValueType { return walkerChaeData->t * first + second; });
265
266 // Add t to all entries of x.
267 storm::utility::vector::applyPointwise(x, x, [this](ValueType const& value) -> ValueType { return value + walkerChaeData->t; });
268
269 // Create a vector that always holds Ax.
270 std::vector<ValueType> currentAx(x.size());
271 walkerChaeData->multiplier->multiply(env, *currentX, nullptr, currentAx);
272
273 // (3) Perform iterations until convergence.
275 uint64_t iterations = 0;
276 this->startMeasureProgress();
277 while (status == SolverStatus::InProgress && iterations < maxIter) {
278 // Perform one Walker-Chae step.
279 walkerChaeData->matrix.performWalkerChaeStep(*currentX, walkerChaeData->columnSums, walkerChaeData->b, currentAx, *nextX);
280
281 // Compute new Ax.
282 walkerChaeData->multiplier->multiply(env, *nextX, nullptr, currentAx);
283
284 // Check for convergence.
285 if (storm::utility::vector::computeSquaredNorm2Difference(currentAx, walkerChaeData->b) <= squaredErrorBound) {
287 }
288
289 // Swap the x vectors for the next iteration.
290 std::swap(currentX, nextX);
291
292 // Potentially show progress.
293 this->showProgressIterative(iterations);
294
295 // Increase iteration count so we can abort if convergence is too slow.
296 ++iterations;
297
299 status = SolverStatus::Aborted;
300 }
301 }
302
303 // If the last iteration did not write to the original x we have to swap the contents, because the
304 // output has to be written to the input parameter x.
305 if (currentX == &walkerChaeData->newX) {
306 std::swap(x, *currentX);
307 }
308
309 // Resize the solution to the right size.
310 x.resize(this->A->getRowCount());
311
312 // Finalize solution vector.
313 storm::utility::vector::applyPointwise(x, x, [this](ValueType const& value) -> ValueType { return value - walkerChaeData->t; });
314
315 if (!this->isCachingEnabled()) {
316 clearCache();
317 }
318
319 this->reportStatus(status, iterations);
320
321 return status == SolverStatus::Converged;
322}
323
324template<typename ValueType>
325typename NativeLinearEquationSolver<ValueType>::PowerIterationResult NativeLinearEquationSolver<ValueType>::performPowerIteration(
326 Environment const& env, std::vector<ValueType>*& currentX, std::vector<ValueType>*& newX, std::vector<ValueType> const& b, ValueType const& precision,
327 bool relative, SolverGuarantee const& guarantee, uint64_t currentIterations, uint64_t maxIterations,
328 storm::solver::MultiplicationStyle const& multiplicationStyle) const {
329 bool useGaussSeidelMultiplication = multiplicationStyle == storm::solver::MultiplicationStyle::GaussSeidel;
330
331 uint64_t iterations = currentIterations;
332 SolverStatus status = this->terminateNow(*currentX, guarantee) ? SolverStatus::TerminatedEarly : SolverStatus::InProgress;
333 while (status == SolverStatus::InProgress && iterations < maxIterations) {
334 if (useGaussSeidelMultiplication) {
335 *newX = *currentX;
336 this->multiplier->multiplyGaussSeidel(env, *newX, &b);
337 } else {
338 this->multiplier->multiply(env, *currentX, &b, *newX);
339 }
340
341 // Check for convergence.
342 if (storm::utility::vector::equalModuloPrecision<ValueType>(*currentX, *newX, precision, relative)) {
344 }
345
346 // Check for termination.
347 std::swap(currentX, newX);
348 ++iterations;
349
350 status = this->updateStatus(status, *currentX, guarantee, iterations, maxIterations);
351
352 // Potentially show progress.
353 this->showProgressIterative(iterations);
354 }
355
356 return PowerIterationResult(iterations - currentIterations, status);
357}
358
359template<typename ValueType>
360bool NativeLinearEquationSolver<ValueType>::solveEquationsPower(Environment const& env, std::vector<ValueType>& x, std::vector<ValueType> const& b) const {
361 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with NativeLinearEquationSolver (Power)");
362 // Prepare the solution vectors.
363 setUpViOperator();
364
366 if (this->hasCustomTerminationCondition()) {
367 if (this->getTerminationCondition().requiresGuarantee(SolverGuarantee::LessOrEqual) && this->hasLowerBound()) {
368 this->createLowerBoundsVector(x);
370 } else if (this->getTerminationCondition().requiresGuarantee(SolverGuarantee::GreaterOrEqual) && this->hasUpperBound()) {
371 this->createUpperBoundsVector(x);
373 }
374 }
375
377 uint64_t numIterations{0};
378 auto viCallback = [&](SolverStatus const& current) {
379 this->showProgressIterative(numIterations);
380 return this->updateStatus(current, x, guarantee, numIterations, env.solver().native().getMaximalNumberOfIterations());
381 };
382 this->startMeasureProgress();
383 auto status = viHelper.VI(x, b, numIterations, env.solver().native().getRelativeTerminationCriterion(),
384 storm::utility::convertNumber<ValueType>(env.solver().native().getPrecision()), {}, viCallback,
385 env.solver().native().getPowerMethodMultiplicationStyle());
386
387 this->reportStatus(status, numIterations);
388
389 if (!this->isCachingEnabled()) {
390 clearCache();
391 }
392
393 return status == SolverStatus::Converged || status == SolverStatus::TerminatedEarly;
394}
395
396template<typename ValueType>
397void preserveOldRelevantValues(std::vector<ValueType> const& allValues, storm::storage::BitVector const& relevantValues, std::vector<ValueType>& oldValues) {
398 storm::utility::vector::selectVectorValues(oldValues, relevantValues, allValues);
399}
400
401template<typename ValueType>
402ValueType computeMaxAbsDiff(std::vector<ValueType> const& allValues, storm::storage::BitVector const& relevantValues, std::vector<ValueType> const& oldValues) {
403 ValueType result = storm::utility::zero<ValueType>();
404 auto oldValueIt = oldValues.begin();
405 for (auto value : relevantValues) {
406 result = storm::utility::max<ValueType>(result, storm::utility::abs<ValueType>(allValues[value] - *oldValueIt));
407 }
408 return result;
409}
410
411template<typename ValueType>
412ValueType computeMaxAbsDiff(std::vector<ValueType> const& allOldValues, std::vector<ValueType> const& allNewValues,
413 storm::storage::BitVector const& relevantValues) {
414 ValueType result = storm::utility::zero<ValueType>();
415 for (auto value : relevantValues) {
416 result = storm::utility::max<ValueType>(result, storm::utility::abs<ValueType>(allNewValues[value] - allOldValues[value]));
417 }
418 return result;
419}
420
421template<typename ValueType>
422bool NativeLinearEquationSolver<ValueType>::solveEquationsIntervalIteration(Environment const& env, std::vector<ValueType>& x,
423 std::vector<ValueType> const& b) const {
424 STORM_LOG_THROW(this->hasLowerBound(), storm::exceptions::UnmetRequirementException, "Solver requires lower bound, but none was given.");
425 STORM_LOG_THROW(this->hasUpperBound(), storm::exceptions::UnmetRequirementException, "Solver requires upper bound, but none was given.");
426 STORM_LOG_INFO("Solving linear equation system (" << x.size() << " rows) with NativeLinearEquationSolver (IntervalIteration)");
427 setUpViOperator();
428 helper::IntervalIterationHelper<ValueType, true> iiHelper(viOperator);
429 auto prec = storm::utility::convertNumber<ValueType>(env.solver().native().getPrecision());
430 auto lowerBoundsCallback = [&](std::vector<ValueType>& vector) { this->createLowerBoundsVector(vector); };
431 auto upperBoundsCallback = [&](std::vector<ValueType>& vector) { this->createUpperBoundsVector(vector); };
432
433 uint64_t numIterations{0};
434 auto iiCallback = [&](helper::IIData<ValueType> const& data) {
435 this->showProgressIterative(numIterations);
436 bool terminateEarly = this->hasCustomTerminationCondition() && this->getTerminationCondition().terminateNow(data.x, SolverGuarantee::LessOrEqual) &&
437 this->getTerminationCondition().terminateNow(data.y, SolverGuarantee::GreaterOrEqual);
438 return this->updateStatus(data.status, terminateEarly, numIterations, env.solver().native().getMaximalNumberOfIterations());
439 };
440 std::optional<storm::storage::BitVector> optionalRelevantValues;
441 if (this->hasRelevantValues()) {
442 optionalRelevantValues = this->getRelevantValues();
443 }
444 this->startMeasureProgress();
445 auto status = iiHelper.II(x, b, numIterations, env.solver().native().getRelativeTerminationCriterion(), prec, lowerBoundsCallback, upperBoundsCallback, {},
446 iiCallback, optionalRelevantValues);
447 this->reportStatus(status, numIterations);
448
449 if (!this->isCachingEnabled()) {
450 clearCache();
451 }
452
453 return status == SolverStatus::Converged || status == SolverStatus::TerminatedEarly;
454}
455
456template<typename ValueType>
457bool NativeLinearEquationSolver<ValueType>::solveEquationsSoundValueIteration(Environment const& env, std::vector<ValueType>& x,
458 std::vector<ValueType> const& b) const {
459 // Prepare the solution vectors and the helper.
460 assert(x.size() == this->A->getRowGroupCount());
461
462 std::optional<ValueType> lowerBound, upperBound;
463 if (this->hasLowerBound()) {
464 lowerBound = this->getLowerBound(true);
465 }
466 if (this->hasUpperBound()) {
467 upperBound = this->getUpperBound(true);
468 }
469
470 setUpViOperator();
471
472 auto precision = storm::utility::convertNumber<ValueType>(env.solver().native().getPrecision());
473 uint64_t numIterations{0};
474 auto sviCallback = [&](typename helper::SoundValueIterationHelper<ValueType, true>::SVIData const& current) {
475 this->showProgressIterative(numIterations);
476 return this->updateStatus(current.status,
477 this->hasCustomTerminationCondition() && current.checkCustomTerminationCondition(this->getTerminationCondition()),
478 numIterations, env.solver().native().getMaximalNumberOfIterations());
479 };
480 std::optional<storm::storage::BitVector> optionalRelevantValues;
481 if (this->hasRelevantValues()) {
482 optionalRelevantValues = this->getRelevantValues();
483 }
484 this->startMeasureProgress();
485 helper::SoundValueIterationHelper<ValueType, true> sviHelper(viOperator);
486 auto status = sviHelper.SVI(x, b, numIterations, env.solver().native().getRelativeTerminationCriterion(), precision, {}, lowerBound, upperBound,
487 sviCallback, optionalRelevantValues);
488
489 this->reportStatus(status, numIterations);
490
491 if (!this->isCachingEnabled()) {
492 clearCache();
493 }
494
495 return status == SolverStatus::Converged || status == SolverStatus::TerminatedEarly;
496}
497
498template<typename ValueType>
499bool NativeLinearEquationSolver<ValueType>::solveEquationsOptimisticValueIteration(Environment const& env, std::vector<ValueType>& x,
500 std::vector<ValueType> const& b) const {
502 // If all entries are zero, OVI might run in an endless loop. However, the result is easy in this case.
503 x.assign(x.size(), storm::utility::zero<ValueType>());
504 return true;
505 }
506
507 setUpViOperator();
508
509 helper::OptimisticValueIterationHelper<ValueType, true> oviHelper(viOperator);
510 auto prec = storm::utility::convertNumber<ValueType>(env.solver().native().getPrecision());
511 std::optional<ValueType> lowerBound, upperBound;
512 if (this->hasLowerBound()) {
513 lowerBound = this->getLowerBound(true);
514 }
515 if (this->hasUpperBound()) {
516 upperBound = this->getUpperBound(true);
517 }
518 uint64_t numIterations{0};
519 auto oviCallback = [&](SolverStatus const& current, std::vector<ValueType> const& v) {
520 this->showProgressIterative(numIterations);
521 return this->updateStatus(current, v, SolverGuarantee::LessOrEqual, numIterations, env.solver().native().getMaximalNumberOfIterations());
522 };
523 this->createLowerBoundsVector(x);
524 std::optional<ValueType> guessingFactor;
525 if (env.solver().ovi().getUpperBoundGuessingFactor()) {
526 guessingFactor = storm::utility::convertNumber<ValueType>(*env.solver().ovi().getUpperBoundGuessingFactor());
527 }
528 this->startMeasureProgress();
529 auto status = oviHelper.OVI(x, b, env.solver().native().getRelativeTerminationCriterion(), prec, {}, guessingFactor, lowerBound, upperBound, oviCallback);
530 this->reportStatus(status, numIterations);
531
532 if (!this->isCachingEnabled()) {
533 clearCache();
534 }
535
536 return status == SolverStatus::Converged || status == SolverStatus::TerminatedEarly;
537}
538
539template<typename ValueType>
540bool NativeLinearEquationSolver<ValueType>::solveEquationsRationalSearch(Environment const& env, std::vector<ValueType>& x,
541 std::vector<ValueType> const& b) const {
542 // Set up two value iteration operators. One for exact and one for imprecise computations
543 setUpViOperator();
544 std::shared_ptr<helper::ValueIterationOperator<storm::RationalNumber, true>> exactOp;
545 std::shared_ptr<helper::ValueIterationOperator<double, true>> impreciseOp;
546
547 if constexpr (std::is_same_v<ValueType, storm::RationalNumber>) {
548 exactOp = viOperator;
549 impreciseOp = std::make_shared<helper::ValueIterationOperator<double, true>>();
550 impreciseOp->setMatrixBackwards(this->A->template toValueType<double>());
551 } else {
552 impreciseOp = viOperator;
553 exactOp = std::make_shared<helper::ValueIterationOperator<storm::RationalNumber, true>>();
554 exactOp->setMatrixBackwards(this->A->template toValueType<storm::RationalNumber>());
555 }
556
558 uint64_t numIterations{0};
559 auto rsCallback = [&](SolverStatus const& current) {
560 this->showProgressIterative(numIterations);
561 return this->updateStatus(current, x, SolverGuarantee::None, numIterations, env.solver().native().getMaximalNumberOfIterations());
562 };
563 this->startMeasureProgress();
564 auto status = rsHelper.RS(x, b, numIterations, storm::utility::convertNumber<ValueType>(env.solver().native().getPrecision()), {}, rsCallback);
565
566 this->reportStatus(status, numIterations);
567
568 if (!this->isCachingEnabled()) {
569 clearCache();
570 }
571
572 return status == SolverStatus::Converged || status == SolverStatus::TerminatedEarly;
573}
574
575template<typename ValueType>
576NativeLinearEquationSolverMethod NativeLinearEquationSolver<ValueType>::getMethod(Environment const& env, bool isExactMode) const {
577 // Adjust the method if none was specified and we want exact or sound computations
578 auto method = env.solver().native().getMethod();
579
580 if (isExactMode && method != NativeLinearEquationSolverMethod::RationalSearch) {
581 if (env.solver().native().isMethodSetFromDefault()) {
582 method = NativeLinearEquationSolverMethod::RationalSearch;
584 "Selecting '" + toString(method) +
585 "' as the solution technique to guarantee exact results. If you want to override this, please explicitly specify a different method.");
586 } else {
587 STORM_LOG_WARN("The selected solution method does not guarantee exact results.");
588 }
589 } else if (env.solver().isForceSoundness() && method != NativeLinearEquationSolverMethod::SoundValueIteration &&
590 method != NativeLinearEquationSolverMethod::OptimisticValueIteration && method != NativeLinearEquationSolverMethod::IntervalIteration &&
591 method != NativeLinearEquationSolverMethod::RationalSearch) {
592 if (env.solver().native().isMethodSetFromDefault()) {
593 method = NativeLinearEquationSolverMethod::OptimisticValueIteration;
595 "Selecting '" + toString(method) +
596 "' as the solution technique to guarantee sound results. If you want to override this, please explicitly specify a different method.");
597 } else {
598 STORM_LOG_WARN("The selected solution method does not guarantee sound results.");
599 }
600 }
601 return method;
602}
603
604template<typename ValueType>
605bool NativeLinearEquationSolver<ValueType>::internalSolveEquations(Environment const& env, std::vector<ValueType>& x, std::vector<ValueType> const& b) const {
606 switch (getMethod(env, storm::NumberTraits<ValueType>::IsExact || env.solver().isForceExact())) {
607 case NativeLinearEquationSolverMethod::SOR:
608 return this->solveEquationsSOR(env, x, b, storm::utility::convertNumber<ValueType>(env.solver().native().getSorOmega()));
609 case NativeLinearEquationSolverMethod::GaussSeidel:
610 return this->solveEquationsSOR(env, x, b, storm::utility::one<ValueType>());
611 case NativeLinearEquationSolverMethod::Jacobi:
612 return this->solveEquationsJacobi(env, x, b);
613 case NativeLinearEquationSolverMethod::WalkerChae:
614 return this->solveEquationsWalkerChae(env, x, b);
615 case NativeLinearEquationSolverMethod::Power:
616 return this->solveEquationsPower(env, x, b);
617 case NativeLinearEquationSolverMethod::SoundValueIteration:
618 return this->solveEquationsSoundValueIteration(env, x, b);
619 case NativeLinearEquationSolverMethod::OptimisticValueIteration:
620 return this->solveEquationsOptimisticValueIteration(env, x, b);
621 case NativeLinearEquationSolverMethod::IntervalIteration:
622 return this->solveEquationsIntervalIteration(env, x, b);
623 case NativeLinearEquationSolverMethod::RationalSearch:
624 return this->solveEquationsRationalSearch(env, x, b);
625 }
626 STORM_LOG_THROW(false, storm::exceptions::InvalidEnvironmentException, "Unknown solving technique.");
627 return false;
628}
629
630template<typename ValueType>
632 auto method = getMethod(env, storm::NumberTraits<ValueType>::IsExact || env.solver().isForceExact());
633 if (method == NativeLinearEquationSolverMethod::Power || method == NativeLinearEquationSolverMethod::SoundValueIteration ||
634 method == NativeLinearEquationSolverMethod::OptimisticValueIteration || method == NativeLinearEquationSolverMethod::RationalSearch ||
635 method == NativeLinearEquationSolverMethod::IntervalIteration) {
637 } else {
639 }
640}
641
642template<typename ValueType>
645 auto method = getMethod(env, storm::NumberTraits<ValueType>::IsExact || env.solver().isForceExact());
646 if (method == NativeLinearEquationSolverMethod::IntervalIteration) {
647 requirements.requireBounds();
648 } else if (method == NativeLinearEquationSolverMethod::RationalSearch || method == NativeLinearEquationSolverMethod::OptimisticValueIteration) {
649 requirements.requireLowerBounds();
650 } else if (method == NativeLinearEquationSolverMethod::SoundValueIteration) {
651 requirements.requireBounds(false);
652 }
653 return requirements;
654}
655
656template<typename ValueType>
658 jacobiDecomposition.reset();
659 walkerChaeData.reset();
660 multiplier.reset();
661 viOperator.reset();
663}
664
665template<typename ValueType>
667 return this->A->getRowCount();
668}
669
670template<typename ValueType>
671uint64_t NativeLinearEquationSolver<ValueType>::getMatrixColumnCount() const {
672 return this->A->getColumnCount();
673}
674
675template<typename ValueType>
676std::unique_ptr<storm::solver::LinearEquationSolver<ValueType>> NativeLinearEquationSolverFactory<ValueType>::create(Environment const&) const {
677 return std::make_unique<storm::solver::NativeLinearEquationSolver<ValueType>>();
678}
679
680template<typename ValueType>
681std::unique_ptr<LinearEquationSolverFactory<ValueType>> NativeLinearEquationSolverFactory<ValueType>::clone() const {
682 return std::make_unique<NativeLinearEquationSolverFactory<ValueType>>(*this);
683}
684
685// Explicitly instantiate the linear equation solver.
688
689#ifdef STORM_HAVE_CARL
692
693#endif
694} // namespace solver
695} // namespace storm
SolverEnvironment & solver()
storm::RationalNumber const & getPrecision() const
uint64_t const & getMaximalNumberOfIterations() const
bool const & getRelativeTerminationCriterion() const
storm::RationalNumber const & getSorOmega() const
NativeSolverEnvironment & native()
LinearEquationSolverRequirements & requireLowerBounds(bool critical=true)
LinearEquationSolverRequirements & requireBounds(bool critical=true)
std::unique_ptr< Multiplier< ValueType > > create(Environment const &env, storm::storage::SparseMatrix< ValueType > const &matrix)
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.
virtual std::unique_ptr< LinearEquationSolverFactory< ValueType > > clone() const override
Creates a copy of this factory.
A class that uses storm's native matrix operations to implement the LinearEquationSolver interface.
virtual LinearEquationSolverProblemFormat getEquationProblemFormat(storm::Environment const &env) const override
Retrieves the format in which this solver expects to solve equations.
virtual bool internalSolveEquations(storm::Environment const &env, std::vector< ValueType > &x, std::vector< ValueType > const &b) const override
virtual void setMatrix(storm::storage::SparseMatrix< ValueType > const &A) override
virtual LinearEquationSolverRequirements getRequirements(Environment const &env) const override
Retrieves the requirements of the solver under the current settings.
A bit vector that is internally represented as a vector of 64-bit values.
Definition BitVector.h:18
A class that can be used to build a sparse matrix by adding value by value.
void addNextValue(index_type row, index_type column, value_type const &value)
Sets the matrix entry at the given row and column to the given value.
SparseMatrix< value_type > build(index_type overriddenRowCount=0, index_type overriddenColumnCount=0, index_type overriddenRowGroupCount=0)
A class that holds a possibly non-square matrix in the compressed row storage format.
index_type getColumnCount() const
Returns the number of columns of the matrix.
std::pair< storm::storage::SparseMatrix< value_type >, std::vector< value_type > > getJacobiDecomposition() const
Calculates the Jacobi decomposition of this sparse matrix.
index_type getRowCount() const
Returns the number of rows of the matrix.
#define STORM_LOG_INFO(message)
Definition logging.h:29
#define STORM_LOG_WARN(message)
Definition logging.h:30
#define STORM_LOG_THROW(cond, exception, message)
Definition macros.h:30
SFTBDDChecker::ValueType ValueType
void preserveOldRelevantValues(std::vector< ValueType > const &allValues, storm::storage::BitVector const &relevantValues, std::vector< ValueType > &oldValues)
std::string toString(GurobiSolverMethod const &method)
Yields a string representation of the GurobiSolverMethod.
ValueType computeMaxAbsDiff(std::vector< ValueType > const &allValues, storm::storage::BitVector const &relevantValues, std::vector< ValueType > const &oldValues)
bool isTerminate()
Check whether the program should terminate (due to some abort signal).
void selectVectorValues(std::vector< T > &vector, storm::storage::BitVector const &positions, std::vector< T > const &values)
Selects the elements from a vector at the specified positions and writes them consecutively into anot...
Definition vector.h:189
T computeSquaredNorm2Difference(std::vector< T > const &b1, std::vector< T > const &b2)
Definition vector.h:1021
void multiplyVectorsPointwise(std::vector< InValueType1 > const &firstOperand, std::vector< InValueType2 > const &secondOperand, std::vector< OutValueType > &target)
Multiplies the two given vectors (pointwise) and writes the result to the target vector.
Definition vector.h:467
void applyPointwise(std::vector< InValueType1 > const &firstOperand, std::vector< InValueType2 > const &secondOperand, std::vector< OutValueType > &target, Operation f=Operation())
Applies the given operation pointwise on the two given vectors and writes the result to the third vec...
Definition vector.h:398
void subtractVectors(std::vector< InValueType1 > const &firstOperand, std::vector< InValueType2 > const &secondOperand, std::vector< OutValueType > &target)
Subtracts the two given vectors and writes the result to the target vector.
Definition vector.h:455
bool hasNonZeroEntry(std::vector< T > const &v)
Definition vector.h:1252
ValueType zero()
Definition constants.cpp:26
ValueType pow(ValueType const &value, int_fast64_t exponent)
ValueType one()
Definition constants.cpp:21
TargetType convertNumber(SourceType const &number)
Definition constants.cpp:98
LabParser.cpp.
Definition cli.cpp:18