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