Storm 1.11.1.1
A Modern Probabilistic Model Checker
Loading...
Searching...
No Matches
IterativeMinMaxLinearEquationSolver.cpp
Go to the documentation of this file.
2
3#include <functional>
4#include <type_traits>
5
29
30namespace storm {
31namespace solver {
32
33template<typename ValueType, typename SolutionType>
35 std::unique_ptr<LinearEquationSolverFactory<SolutionType>>&& linearEquationSolverFactory)
36 : linearEquationSolverFactory(std::move(linearEquationSolverFactory)) {
37 // Intentionally left empty
38}
39
40template<typename ValueType, typename SolutionType>
42 storm::storage::SparseMatrix<ValueType> const& A, std::unique_ptr<LinearEquationSolverFactory<SolutionType>>&& linearEquationSolverFactory)
43 : StandardMinMaxLinearEquationSolver<ValueType, SolutionType>(A), linearEquationSolverFactory(std::move(linearEquationSolverFactory)) {
44 // Intentionally left empty.
45}
46
47template<typename ValueType, typename SolutionType>
49 storm::storage::SparseMatrix<ValueType>&& A, std::unique_ptr<LinearEquationSolverFactory<SolutionType>>&& linearEquationSolverFactory)
50 : StandardMinMaxLinearEquationSolver<ValueType, SolutionType>(std::move(A)), linearEquationSolverFactory(std::move(linearEquationSolverFactory)) {
51 // Intentionally left empty.
52}
53
54template<typename ValueType, typename SolutionType>
55MinMaxMethod IterativeMinMaxLinearEquationSolver<ValueType, SolutionType>::getMethod(Environment const& env, bool isExactMode) const {
56 // Adjust the method if none was specified and we want exact or sound computations.
57 auto method = env.solver().minMax().getMethod();
58
59 if (isExactMode && method != MinMaxMethod::PolicyIteration && method != MinMaxMethod::RationalSearch && method != MinMaxMethod::ViToPi) {
60 if (env.solver().minMax().isMethodSetFromDefault()) {
62 "Selecting 'Policy iteration' as the solution technique to guarantee exact results. If you want to override this, please explicitly specify a "
63 "different method.");
64 method = MinMaxMethod::PolicyIteration;
65 } else {
66 STORM_LOG_WARN("The selected solution method " << toString(method) << " does not guarantee exact results.");
67 }
68 } else if (env.solver().isForceSoundness() && method != MinMaxMethod::SoundValueIteration && method != MinMaxMethod::IntervalIteration &&
69 method != MinMaxMethod::PolicyIteration && method != MinMaxMethod::RationalSearch && method != MinMaxMethod::OptimisticValueIteration &&
70 method != MinMaxMethod::GuessingValueIteration) {
71 if (env.solver().minMax().isMethodSetFromDefault()) {
72 method = MinMaxMethod::OptimisticValueIteration;
74 "Selecting '"
75 << toString(method)
76 << "' as the solution technique to guarantee sound results. If you want to override this, please explicitly specify a different method.");
77 } else {
78 STORM_LOG_WARN("The selected solution method does not guarantee sound results.");
79 }
80 }
81 STORM_LOG_THROW(method == MinMaxMethod::ValueIteration || method == MinMaxMethod::PolicyIteration || method == MinMaxMethod::RationalSearch ||
82 method == MinMaxMethod::SoundValueIteration || method == MinMaxMethod::IntervalIteration ||
83 method == MinMaxMethod::OptimisticValueIteration || method == MinMaxMethod::GuessingValueIteration || method == MinMaxMethod::ViToPi,
84 storm::exceptions::InvalidEnvironmentException, "This solver does not support the selected method '" << toString(method) << "'.");
85 return method;
86}
87
88template<typename ValueType, typename SolutionType>
90 std::vector<SolutionType>& x, std::vector<ValueType> const& b) const {
91 bool result = false;
92 switch (getMethod(env, storm::NumberTraits<ValueType>::IsExact || env.solver().isForceExact())) {
93 case MinMaxMethod::ValueIteration:
94 result = solveEquationsValueIteration(env, dir, x, b);
95 break;
96 case MinMaxMethod::OptimisticValueIteration:
97 result = solveEquationsOptimisticValueIteration(env, dir, x, b);
98 break;
99 case MinMaxMethod::GuessingValueIteration:
100 result = solveEquationsGuessingValueIteration(env, dir, x, b);
101 break;
102 case MinMaxMethod::PolicyIteration:
103 result = solveEquationsPolicyIteration(env, dir, x, b);
104 break;
105 case MinMaxMethod::RationalSearch:
106 result = solveEquationsRationalSearch(env, dir, x, b);
107 break;
108 case MinMaxMethod::IntervalIteration:
109 result = solveEquationsIntervalIteration(env, dir, x, b);
110 break;
111 case MinMaxMethod::SoundValueIteration:
112 result = solveEquationsSoundValueIteration(env, dir, x, b);
113 break;
114 case MinMaxMethod::ViToPi:
115 result = solveEquationsViToPi(env, dir, x, b);
116 break;
117 default:
118 STORM_LOG_THROW(false, storm::exceptions::InvalidEnvironmentException, "This solver does not implement the selected solution method");
119 }
120
121 return result;
122}
123
124template<typename ValueType, typename SolutionType>
126 if (!viOperatorTriv && !viOperatorNontriv) {
127 if (this->A->hasTrivialRowGrouping()) {
128 // The trivial row grouping minmax operator makes sense over intervals.
129 viOperatorTriv = std::make_shared<helper::ValueIterationOperator<ValueType, true, SolutionType>>();
130 viOperatorTriv->setMatrixBackwards(*this->A);
131 if constexpr (!std::is_same_v<ValueType, storm::Interval>) {
132 // It might be that someone is using a minmaxlinearequationsolver with an advanced VI algorithm
133 // but is just passing a DTMC over doubles. In this case we need to populate this VI operator.
134 // It behaves exactly the same as the trivial row grouping operator, but it is currently hardcoded
135 // to be used by, e.g., optimistic value iteration.
136 viOperatorNontriv = std::make_shared<helper::ValueIterationOperator<ValueType, false, SolutionType>>();
137 viOperatorNontriv->setMatrixBackwards(*this->A);
138 }
139 } else {
140 // The nontrivial row grouping minmax operator makes sense for MDPs.
141 viOperatorNontriv = std::make_shared<helper::ValueIterationOperator<ValueType, false, SolutionType>>();
142 viOperatorNontriv->setMatrixBackwards(*this->A);
143 }
144 }
145 if (this->choiceFixedForRowGroup) {
146 // Ignore those rows that are not selected
147 assert(this->initialScheduler);
148 auto callback = [&](uint64_t groupIndex, uint64_t localRowIndex) {
149 return this->choiceFixedForRowGroup->get(groupIndex) && this->initialScheduler->at(groupIndex) != localRowIndex;
150 };
151 if (viOperatorTriv) {
152 viOperatorTriv->setIgnoredRows(true, callback);
153 }
154 if (viOperatorNontriv) {
155 viOperatorNontriv->setIgnoredRows(true, callback);
156 }
157 }
158}
159
160template<typename ValueType, typename SolutionType>
161void IterativeMinMaxLinearEquationSolver<ValueType, SolutionType>::extractScheduler(std::vector<SolutionType>& x, std::vector<ValueType> const& b,
162 OptimizationDirection const& dir, bool updateX, bool robust) const {
163 if (std::is_same_v<ValueType, storm::Interval> && this->A->hasTrivialRowGrouping()) {
164 // Create robust scheduler index if it doesn't exist yet
165 if (!this->robustSchedulerIndex) {
166 this->robustSchedulerIndex = std::vector<uint64_t>(x.size(), 0);
167 } else {
168 this->robustSchedulerIndex->resize(x.size(), 0);
169 }
170 uint64_t numSchedulerChoices = 0;
171 for (uint64_t row = 0; row < this->A->getRowCount(); ++row) {
172 this->robustSchedulerIndex->at(row) = numSchedulerChoices;
173 numSchedulerChoices += this->A->getRow(row).getNumberOfEntries();
174 }
175 // Make sure that storage for scheduler choices is available
176 if (!this->schedulerChoices) {
177 this->schedulerChoices = std::vector<uint64_t>(numSchedulerChoices, 0);
178 } else {
179 this->schedulerChoices->resize(numSchedulerChoices, 0);
180 }
181 } else {
182 // Make sure that storage for scheduler choices is available
183 if (!this->schedulerChoices) {
184 this->schedulerChoices = std::vector<uint64_t>(x.size(), 0);
185 } else {
186 this->schedulerChoices->resize(x.size(), 0);
187 }
188 }
189
190 // Set the correct choices.
191 if (!viOperatorTriv && !viOperatorNontriv) {
192 STORM_LOG_WARN("Expected VI operator to be initialized for scheduler extraction. Initializing now, but this is inefficient.");
193 setUpViOperator();
194 }
195 if (viOperatorTriv) {
196 if constexpr (std::is_same<ValueType, storm::Interval>() && std::is_same<SolutionType, double>()) {
198 schedHelper.computeScheduler(x, b, dir, *this->schedulerChoices, robust, updateX ? &x : nullptr, this->robustSchedulerIndex);
199 } else {
200 STORM_LOG_ERROR("SchedulerTrackingHelper not implemented for this setting (trivial row grouping but not Interval->double).");
201 }
202 }
203 if (viOperatorNontriv) {
205 schedHelper.computeScheduler(x, b, dir, *this->schedulerChoices, robust, updateX ? &x : nullptr, this->robustSchedulerIndex);
206 }
207}
208
209template<typename ValueType, typename SolutionType>
210bool IterativeMinMaxLinearEquationSolver<ValueType, SolutionType>::solveInducedEquationSystem(
211 Environment const& env, std::unique_ptr<LinearEquationSolver<SolutionType>>& linearEquationSolver, std::vector<uint64_t> const& scheduler,
212 std::vector<SolutionType>& x, std::vector<ValueType>& subB, std::vector<ValueType> const& originalB, OptimizationDirection dir) const {
213 if constexpr (std::is_same_v<ValueType, storm::Interval>) {
214 if constexpr (std::is_same_v<SolutionType, storm::Interval>) {
215 STORM_LOG_THROW(false, storm::exceptions::NotImplementedException,
216 "We did not implement solving induced equation systems for interval-based models outside of robust VI.");
217 // Implementing this requires linear equation systems with different value types and solution types (or some appropriate casting)
218 return false;
219 }
220 STORM_LOG_ASSERT(subB.size() == x.size(), "Sizes of subB and x do not coincide.");
221 STORM_LOG_ASSERT(this->linearEquationSolverFactory != nullptr, "Wrong constructor was called.");
222 STORM_LOG_THROW(this->linearEquationSolverFactory->getEquationProblemFormat(env) == LinearEquationSolverProblemFormat::FixedPointSystem,
223 storm::exceptions::NotImplementedException, "Solving induced system of Interval Model not supported for the selected equation solver");
224
225 storm::storage::SparseMatrixBuilder<SolutionType> newMatrixBuilder(this->A->getRowCount(), this->A->getColumnCount(), this->A->getEntryCount());
226
227 // Robust VI scheduler is an order, compute the correct values for this order
228 auto schedulerIterator = scheduler.begin();
229 for (uint64_t rowIndex = 0; rowIndex < this->A->getRowCount(); rowIndex++) {
230 auto const& row = this->A->getRow(rowIndex);
231
232 static std::vector<SolutionType> tmp;
233 tmp.clear();
234
235 SolutionType probLeft = storm::utility::one<SolutionType>();
236
237 for (auto const& entry : row) {
238 tmp.push_back(entry.getValue().lower());
239 probLeft -= entry.getValue().lower();
240 }
241
242 auto const& rowIter = row.begin();
243 for (uint64_t i = 0; i < row.getNumberOfEntries(); i++, schedulerIterator++) {
244 if (!utility::isZero(probLeft)) {
245 auto const& entry = rowIter[*schedulerIterator];
246 auto const diameter = entry.getValue().upper() - entry.getValue().lower();
247 auto const value = utility::min(probLeft, diameter);
248 tmp[*schedulerIterator] += value;
249 probLeft -= value;
250 } else {
251 // Intentionally left empty: advance schedulerIterator to end of row
252 }
253 }
254
255 for (uint64_t i = 0; i < row.getNumberOfEntries(); i++) {
256 auto const& entry = rowIter[i];
257 newMatrixBuilder.addNextValue(rowIndex, entry.getColumn(), tmp[i]);
258 }
259 }
260 STORM_LOG_ASSERT(schedulerIterator == scheduler.end(), "Offset issue in scheduler?");
261
262 subB = originalB;
263
264 std::vector<SolutionType> b;
265 for (auto const& entry : subB) {
266 if (dir == OptimizationDirection::Maximize) {
267 b.push_back(entry.upper());
268 } else {
269 b.push_back(entry.lower());
270 }
271 }
272
273 auto const submatrix = newMatrixBuilder.build();
274
275 // Check whether the linear equation solver is already initialized
276 if (!linearEquationSolver) {
277 // Initialize the equation solver
278 linearEquationSolver = this->linearEquationSolverFactory->create(env, std::move(submatrix));
279 linearEquationSolver->setBoundsFromOtherSolver(*this);
280 linearEquationSolver->setCachingEnabled(true);
281 } else {
282 // If the equation solver is already initialized, it suffices to update the matrix
283 linearEquationSolver->setMatrix(std::move(submatrix));
284 }
285 // Solve the equation system for the 'DTMC' and return true upon success
286 return linearEquationSolver->solveEquations(env, x, b);
287 } else {
288 STORM_LOG_ASSERT(subB.size() == x.size(), "Sizes of subB and x do not coincide.");
289 STORM_LOG_ASSERT(this->linearEquationSolverFactory != nullptr, "Wrong constructor was called.");
290
291 // Resolve the nondeterminism according to the given scheduler.
292 bool convertToEquationSystem = this->linearEquationSolverFactory->getEquationProblemFormat(env) == LinearEquationSolverProblemFormat::EquationSystem;
294
295 submatrix = this->A->selectRowsFromRowGroups(scheduler, convertToEquationSystem);
296 if (convertToEquationSystem) {
297 submatrix.convertToEquationSystem();
298 }
299 storm::utility::vector::selectVectorValues<ValueType>(subB, scheduler, this->A->getRowGroupIndices(), originalB);
300
301 // Check whether the linear equation solver is already initialized
302 if (!linearEquationSolver) {
303 // Initialize the equation solver
304 linearEquationSolver = this->linearEquationSolverFactory->create(env, std::move(submatrix));
305 linearEquationSolver->setBoundsFromOtherSolver(*this);
306 linearEquationSolver->setCachingEnabled(true);
307 } else {
308 // If the equation solver is already initialized, it suffices to update the matrix
309 linearEquationSolver->setMatrix(std::move(submatrix));
310 }
311 // Solve the equation system for the 'DTMC' and return true upon success
312 return linearEquationSolver->solveEquations(env, x, subB);
313 }
314}
315
316template<typename ValueType, typename SolutionType>
317bool IterativeMinMaxLinearEquationSolver<ValueType, SolutionType>::solveEquationsPolicyIteration(Environment const& env, OptimizationDirection dir,
318 std::vector<SolutionType>& x,
319 std::vector<ValueType> const& b) const {
320 std::vector<storm::storage::sparse::state_type> scheduler =
321 this->hasInitialScheduler() ? this->getInitialScheduler() : std::vector<storm::storage::sparse::state_type>(this->A->getRowGroupCount());
322 return performPolicyIteration(env, dir, x, b, std::move(scheduler));
323}
324
325template<typename ValueType, typename SolutionType>
326bool IterativeMinMaxLinearEquationSolver<ValueType, SolutionType>::performPolicyIteration(
327 Environment const& env, OptimizationDirection dir, std::vector<SolutionType>& x, std::vector<ValueType> const& b,
328 std::vector<storm::storage::sparse::state_type>&& initialPolicy) const {
329 if constexpr (std::is_same_v<ValueType, storm::Interval>) {
330 STORM_LOG_THROW(false, storm::exceptions::NotImplementedException, "We did not implement policy iteration for interval-based models.");
331 return false;
332 } else {
333 std::vector<storm::storage::sparse::state_type> scheduler = std::move(initialPolicy);
334 // Get a vector for storing the right-hand side of the inner equation system.
335 if (!auxiliaryRowGroupVector) {
336 auxiliaryRowGroupVector = std::make_unique<std::vector<ValueType>>(this->A->getRowGroupCount());
337 }
338 std::vector<ValueType>& subB = *auxiliaryRowGroupVector;
339
340 // The solver that we will use throughout the procedure.
341 std::unique_ptr<storm::solver::LinearEquationSolver<ValueType>> solver;
342 // The linear equation solver should be at least as precise as this solver
343 std::unique_ptr<storm::Environment> environmentOfSolverStorage;
344 auto precOfSolver = env.solver().getPrecisionOfLinearEquationSolver(env.solver().getLinearEquationSolverType());
346 bool changePrecision = precOfSolver.first && precOfSolver.first.get() > env.solver().minMax().getPrecision();
347 bool changeRelative = precOfSolver.second && !precOfSolver.second.get() && env.solver().minMax().getRelativeTerminationCriterion();
348 if (changePrecision || changeRelative) {
349 environmentOfSolverStorage = std::make_unique<storm::Environment>(env);
350 boost::optional<storm::RationalNumber> newPrecision;
351 boost::optional<bool> newRelative;
352 if (changePrecision) {
353 newPrecision = env.solver().minMax().getPrecision();
354 }
355 if (changeRelative) {
356 newRelative = true;
357 }
358 environmentOfSolverStorage->solver().setLinearEquationSolverPrecision(newPrecision, newRelative);
359 }
360 }
361 storm::Environment const& environmentOfSolver = environmentOfSolverStorage ? *environmentOfSolverStorage : env;
362
364 uint64_t iterations = 0;
365 this->startMeasureProgress();
366 do {
367 // Solve the equation system for the 'DTMC'.
368 solveInducedEquationSystem(environmentOfSolver, solver, scheduler, x, subB, b, dir);
369
370 // Go through the multiplication result and see whether we can improve any of the choices.
371 bool schedulerImproved = false;
372 // Group refers to the state number
373 for (uint_fast64_t group = 0; group < this->A->getRowGroupCount(); ++group) {
374 if (!this->choiceFixedForRowGroup || !this->choiceFixedForRowGroup.get()[group]) {
375 // Only update when the choice is not fixed
376 uint_fast64_t currentChoice = scheduler[group];
377 for (uint_fast64_t choice = this->A->getRowGroupIndices()[group]; choice < this->A->getRowGroupIndices()[group + 1]; ++choice) {
378 // If the choice is the currently selected one, we can skip it.
379 if (choice - this->A->getRowGroupIndices()[group] == currentChoice) {
380 continue;
381 }
382
383 // Create the value of the choice.
384 ValueType choiceValue = storm::utility::zero<ValueType>();
385 for (auto const& entry : this->A->getRow(choice)) {
386 choiceValue += entry.getValue() * x[entry.getColumn()];
387 }
388 choiceValue += b[choice];
389
390 // If the value is strictly better than the solution of the inner system, we need to improve the scheduler.
391 // TODO: If the underlying solver is not precise, this might run forever (i.e. when a state has two choices where the (exact) values are
392 // equal). only changing the scheduler if the values are not equal (modulo precision) would make this unsound.
393 if (valueImproved(dir, x[group], choiceValue)) {
394 schedulerImproved = true;
395 scheduler[group] = choice - this->A->getRowGroupIndices()[group];
396 x[group] = std::move(choiceValue);
397 }
398 }
399 }
400 }
401
402 // If the scheduler did not improve, we are done.
403 if (!schedulerImproved) {
405 }
406
407 // Update environment variables.
408 ++iterations;
409 status =
410 this->updateStatus(status, x, dir == storm::OptimizationDirection::Minimize ? SolverGuarantee::GreaterOrEqual : SolverGuarantee::LessOrEqual,
411 iterations, env.solver().minMax().getMaximalNumberOfIterations());
412
413 // Potentially show progress.
414 this->showProgressIterative(iterations);
415 } while (status == SolverStatus::InProgress);
416
417 STORM_LOG_INFO("Number of iterations: " << iterations);
418 this->reportStatus(status, iterations);
419
420 // If requested, we store the scheduler for retrieval.
421 if (this->isTrackSchedulerSet()) {
422 this->schedulerChoices = std::move(scheduler);
423 }
424
425 if (!this->isCachingEnabled()) {
426 clearCache();
427 }
428
429 return status == SolverStatus::Converged || status == SolverStatus::TerminatedEarly;
430 }
431}
432
433template<typename ValueType, typename SolutionType>
434bool IterativeMinMaxLinearEquationSolver<ValueType, SolutionType>::valueImproved(OptimizationDirection dir, ValueType const& value1,
435 ValueType const& value2) const {
436 if (dir == OptimizationDirection::Minimize) {
437 return value2 < value1;
438 } else {
439 return value2 > value1;
440 }
441}
442
443template<typename ValueType, typename SolutionType>
445 Environment const& env, boost::optional<storm::solver::OptimizationDirection> const& direction, bool const& hasInitialScheduler) const {
446 auto method = getMethod(env, storm::NumberTraits<ValueType>::IsExact || env.solver().isForceExact());
447
448 // Check whether a linear equation solver is needed and potentially start with its requirements
449 bool needsLinEqSolver = false;
450 needsLinEqSolver |= method == MinMaxMethod::PolicyIteration;
451 needsLinEqSolver |= method == MinMaxMethod::ValueIteration && (this->hasInitialScheduler() || hasInitialScheduler);
452 needsLinEqSolver |= method == MinMaxMethod::ViToPi;
453
455 if constexpr (std::is_same_v<ValueType, storm::Interval>) {
456 STORM_LOG_ASSERT(!needsLinEqSolver, "Intervals should not require a linear equation solver.");
457 // nothing to be done;
458 } else if (needsLinEqSolver) {
459 requirements = MinMaxLinearEquationSolverRequirements(this->linearEquationSolverFactory->getRequirements(env));
460 } else {
461 // nothing to be done.
462 }
463
464 if (method == MinMaxMethod::ValueIteration) {
465 if (!this->hasUniqueSolution()) { // Traditional value iteration has no requirements if the solution is unique.
466 // Computing a scheduler is only possible if the solution is unique
467 if (env.solver().minMax().isForceRequireUnique() || this->isTrackSchedulerSet()) {
468 requirements.requireUniqueSolution();
469 } else {
470 // As we want the smallest (largest) solution for maximizing (minimizing) equation systems, we have to approach the solution from below (above).
471 if (!direction || direction.get() == OptimizationDirection::Maximize) {
472 requirements.requireLowerBounds();
473 }
474 if (!direction || direction.get() == OptimizationDirection::Minimize) {
475 requirements.requireUpperBounds();
476 }
477 }
478 }
479 } else if (method == MinMaxMethod::OptimisticValueIteration) {
480 // OptimisticValueIteration always requires lower bounds and a unique solution.
481 if (!this->hasUniqueSolution()) {
482 requirements.requireUniqueSolution();
483 }
484 requirements.requireLowerBounds();
485
486 } else if (method == MinMaxMethod::GuessingValueIteration) {
487 // Guessing value iteration requires a unique solution and lower+upper bounds
488 if (!this->hasUniqueSolution()) {
489 requirements.requireUniqueSolution();
490 }
491 requirements.requireBounds();
492 } else if (method == MinMaxMethod::IntervalIteration) {
493 // Interval iteration requires a unique solution and lower+upper bounds
494 if (!this->hasUniqueSolution()) {
495 requirements.requireUniqueSolution();
496 }
497 requirements.requireBounds();
498 } else if (method == MinMaxMethod::RationalSearch) {
499 // Rational search needs to approach the solution from below.
500 requirements.requireLowerBounds();
501 // The solution needs to be unique in case of minimizing or in cases where we want a scheduler.
502 if (!this->hasUniqueSolution()) {
503 // RationalSearch guesses and verifies a fixpoint and terminates once a fixpoint is found. To ensure that the guessed fixpoint is the
504 // correct one, we enforce uniqueness.
505 requirements.requireUniqueSolution();
506 }
507 } else if (method == MinMaxMethod::PolicyIteration) {
508 // The initial scheduler shall not select an end component
509 if (!this->hasUniqueSolution() && env.solver().minMax().isForceRequireUnique()) {
510 requirements.requireUniqueSolution();
511 }
512 if (!this->hasNoEndComponents() && !this->hasInitialScheduler()) {
513 requirements.requireValidInitialScheduler();
514 }
515 } else if (method == MinMaxMethod::SoundValueIteration) {
516 if (!this->hasUniqueSolution()) {
517 requirements.requireUniqueSolution();
518 }
519 requirements.requireBounds(false);
520 } else if (method == MinMaxMethod::ViToPi) {
521 // Since we want to use value iteration to extract an initial scheduler, the solution has to be unique.
522 if (!this->hasUniqueSolution()) {
523 requirements.requireUniqueSolution();
524 }
525 } else {
526 STORM_LOG_THROW(false, storm::exceptions::InvalidEnvironmentException, "Unsupported technique for iterative MinMax linear equation solver.");
527 }
528 return requirements;
529}
530
531template<typename ValueType, typename SolutionType>
532ValueType computeMaxAbsDiff(std::vector<ValueType> const& allValues, storm::storage::BitVector const& relevantValues, std::vector<ValueType> const& oldValues) {
533 ValueType result = storm::utility::zero<ValueType>();
534 auto oldValueIt = oldValues.begin();
535 for (auto value : relevantValues) {
536 result = storm::utility::max<ValueType>(result, storm::utility::abs<ValueType>(allValues[value] - *oldValueIt));
537 ++oldValueIt;
538 }
539 return result;
540}
541
542template<typename ValueType, typename SolutionType>
543ValueType computeMaxAbsDiff(std::vector<ValueType> const& allOldValues, std::vector<ValueType> const& allNewValues,
544 storm::storage::BitVector const& relevantValues) {
545 ValueType result = storm::utility::zero<ValueType>();
546 for (auto value : relevantValues) {
547 result = storm::utility::max<ValueType>(result, storm::utility::abs<ValueType>(allNewValues[value] - allOldValues[value]));
548 }
549 return result;
550}
551
552template<typename ValueType, typename SolutionType>
553bool IterativeMinMaxLinearEquationSolver<ValueType, SolutionType>::solveEquationsOptimisticValueIteration(Environment const& env, OptimizationDirection dir,
554 std::vector<SolutionType>& x,
555 std::vector<ValueType> const& b) const {
556 if constexpr (std::is_same_v<ValueType, storm::Interval>) {
557 STORM_LOG_THROW(false, storm::exceptions::NotImplementedException, "We did not implement optimistic value iteration for interval-based models.");
558 return false;
559 } else {
561 // If all entries are zero, OVI might run in an endless loop. However, the result is easy in this case.
562 x.assign(x.size(), storm::utility::zero<SolutionType>());
563 if (this->isTrackSchedulerSet()) {
564 this->schedulerChoices = std::vector<uint_fast64_t>(x.size(), 0);
565 }
566 return true;
567 }
568
569 setUpViOperator();
570
571 helper::OptimisticValueIterationHelper<ValueType, false> oviHelper(viOperatorNontriv);
572 auto prec = storm::utility::convertNumber<ValueType>(env.solver().minMax().getPrecision());
573 std::optional<ValueType> lowerBound, upperBound;
574 if (this->hasLowerBound()) {
575 lowerBound = this->getLowerBound(true);
576 }
577 if (this->hasUpperBound()) {
578 upperBound = this->getUpperBound(true);
579 }
580 uint64_t numIterations{0};
581 auto oviCallback = [&](SolverStatus const& current, std::vector<ValueType> const& v) {
582 this->showProgressIterative(numIterations);
583 return this->updateStatus(current, v, SolverGuarantee::LessOrEqual, numIterations, env.solver().minMax().getMaximalNumberOfIterations());
584 };
585 this->createLowerBoundsVector(x);
586 std::optional<ValueType> guessingFactor;
587 if (env.solver().ovi().getUpperBoundGuessingFactor()) {
588 guessingFactor = storm::utility::convertNumber<ValueType>(*env.solver().ovi().getUpperBoundGuessingFactor());
589 }
590 this->startMeasureProgress();
591 auto status = oviHelper.OVI(x, b, numIterations, env.solver().minMax().getRelativeTerminationCriterion(), prec, dir, guessingFactor, lowerBound,
592 upperBound, oviCallback);
593 this->reportStatus(status, numIterations);
594
595 // If requested, we store the scheduler for retrieval.
596 if (this->isTrackSchedulerSet()) {
597 this->extractScheduler(x, b, dir, this->isUncertaintyRobust());
598 }
599
600 if (!this->isCachingEnabled()) {
601 clearCache();
602 }
603
604 return status == SolverStatus::Converged || status == SolverStatus::TerminatedEarly;
605 }
606}
607
608template<typename ValueType, typename SolutionType>
609bool IterativeMinMaxLinearEquationSolver<ValueType, SolutionType>::solveEquationsGuessingValueIteration(Environment const& env, OptimizationDirection dir,
610 std::vector<SolutionType>& x,
611 std::vector<ValueType> const& b) const {
612 if constexpr (std::is_same_v<ValueType, storm::Interval>) {
613 STORM_LOG_THROW(false, storm::exceptions::NotImplementedException, "We did not implement guessing value iteration for interval-based models.");
614 return false;
615 } else {
616 setUpViOperator();
617
618 auto& lowerX = x;
619 auto upperX = std::make_unique<std::vector<SolutionType>>(x.size());
620
622
623 uint64_t numIterations{0};
624
625 auto gviCallback = [&](helper::GVIData<ValueType> const& data) {
626 this->showProgressIterative(numIterations);
627 bool terminateEarly = this->hasCustomTerminationCondition() && this->getTerminationCondition().terminateNow(data.x, SolverGuarantee::LessOrEqual) &&
628 this->getTerminationCondition().terminateNow(data.y, SolverGuarantee::GreaterOrEqual);
629 return this->updateStatus(data.status, terminateEarly, numIterations, env.solver().minMax().getMaximalNumberOfIterations());
630 };
631
632 this->createLowerBoundsVector(lowerX);
633 this->createUpperBoundsVector(*upperX);
634
635 this->startMeasureProgress();
636 auto statusIters = helper.solveEquations(lowerX, *upperX, b, numIterations,
637 storm::utility::convertNumber<ValueType>(env.solver().minMax().getPrecision()), dir, gviCallback);
638 auto two = storm::utility::convertNumber<ValueType>(2.0);
639 storm::utility::vector::applyPointwise<ValueType, ValueType, ValueType>(
640 lowerX, *upperX, x, [&two](ValueType const& a, ValueType const& b) -> ValueType { return (a + b) / two; });
641
642 this->reportStatus(statusIters, numIterations);
643
644 // If requested, we store the scheduler for retrieval.
645 if (this->isTrackSchedulerSet()) {
646 this->extractScheduler(x, b, dir, this->isUncertaintyRobust());
647 }
648
649 if (!this->isCachingEnabled()) {
650 clearCache();
651 }
652
653 return statusIters == SolverStatus::Converged || statusIters == SolverStatus::TerminatedEarly;
654 }
655}
656
657template<typename ValueType, typename SolutionType>
658bool IterativeMinMaxLinearEquationSolver<ValueType, SolutionType>::solveEquationsValueIteration(Environment const& env, OptimizationDirection dir,
659 std::vector<SolutionType>& x,
660 std::vector<ValueType> const& b) const {
661 setUpViOperator();
662 // By default, we can not provide any guarantee
664
665 if (this->hasInitialScheduler()) {
666 if (!auxiliaryRowGroupVector) {
667 auxiliaryRowGroupVector = std::make_unique<std::vector<ValueType>>(this->A->getRowGroupCount());
668 }
669 // Solve the equation system induced by the initial scheduler.
670 std::unique_ptr<storm::solver::LinearEquationSolver<SolutionType>> linEqSolver;
671 // The linear equation solver should be at least as precise as this solver
672 std::unique_ptr<storm::Environment> environmentOfSolverStorage;
673 auto precOfSolver = env.solver().getPrecisionOfLinearEquationSolver(env.solver().getLinearEquationSolverType());
675 bool changePrecision = precOfSolver.first && precOfSolver.first.get() > env.solver().minMax().getPrecision();
676 bool changeRelative = precOfSolver.second && !precOfSolver.second.get() && env.solver().minMax().getRelativeTerminationCriterion();
677 if (changePrecision || changeRelative) {
678 environmentOfSolverStorage = std::make_unique<storm::Environment>(env);
679 boost::optional<storm::RationalNumber> newPrecision;
680 boost::optional<bool> newRelative;
681 if (changePrecision) {
682 newPrecision = env.solver().minMax().getPrecision();
683 }
684 if (changeRelative) {
685 newRelative = true;
686 }
687 environmentOfSolverStorage->solver().setLinearEquationSolverPrecision(newPrecision, newRelative);
688 }
689 }
690 storm::Environment const& environmentOfSolver = environmentOfSolverStorage ? *environmentOfSolverStorage : env;
691
692 bool success = solveInducedEquationSystem(environmentOfSolver, linEqSolver, this->getInitialScheduler(), x, *auxiliaryRowGroupVector, b, dir);
693 if (success) {
694 // If we were given an initial scheduler and are maximizing (minimizing), our current solution becomes
695 // always less-or-equal (greater-or-equal) than the actual solution.
697 } else {
698 guarantee = SolverGuarantee::None;
699 }
700 } else if (!this->hasUniqueSolution()) {
701 if (maximize(dir)) {
702 this->createLowerBoundsVector(x);
704 } else {
705 this->createUpperBoundsVector(x);
707 }
708 } else if (this->hasCustomTerminationCondition()) {
709 if (this->getTerminationCondition().requiresGuarantee(SolverGuarantee::LessOrEqual) && this->hasLowerBound()) {
710 this->createLowerBoundsVector(x);
712 } else if (this->getTerminationCondition().requiresGuarantee(SolverGuarantee::GreaterOrEqual) && this->hasUpperBound()) {
713 this->createUpperBoundsVector(x);
715 }
716 }
717
718 uint64_t numIterations{0};
719 auto viCallback = [&](SolverStatus const& current) {
720 this->showProgressIterative(numIterations);
721 return this->updateStatus(current, x, guarantee, numIterations, env.solver().minMax().getMaximalNumberOfIterations());
722 };
723 this->startMeasureProgress();
724 // This code duplication is necessary because the helper class is different for the two cases.
725 if (this->A->hasTrivialRowGrouping()) {
727
728 auto status = viHelper.VI(x, b, numIterations, env.solver().minMax().getRelativeTerminationCriterion(),
729 storm::utility::convertNumber<SolutionType>(env.solver().minMax().getPrecision()), dir, viCallback,
730 env.solver().minMax().getMultiplicationStyle(), this->isUncertaintyRobust());
731 this->reportStatus(status, numIterations);
732
733 // If requested, we store the scheduler for retrieval.
734 if (this->isTrackSchedulerSet()) {
735 this->extractScheduler(x, b, dir, this->isUncertaintyRobust());
736 }
737
738 if (!this->isCachingEnabled()) {
739 clearCache();
740 }
741
742 return status == SolverStatus::Converged || status == SolverStatus::TerminatedEarly;
743 } else {
745
746 auto status = viHelper.VI(x, b, numIterations, env.solver().minMax().getRelativeTerminationCriterion(),
747 storm::utility::convertNumber<SolutionType>(env.solver().minMax().getPrecision()), dir, viCallback,
748 env.solver().minMax().getMultiplicationStyle(), this->isUncertaintyRobust());
749 this->reportStatus(status, numIterations);
750
751 // If requested, we store the scheduler for retrieval.
752 if (this->isTrackSchedulerSet()) {
753 this->extractScheduler(x, b, dir, this->isUncertaintyRobust());
754 }
755
756 if (!this->isCachingEnabled()) {
757 clearCache();
758 }
759
760 return status == SolverStatus::Converged || status == SolverStatus::TerminatedEarly;
761 }
762}
763
764template<typename ValueType, typename SolutionType>
765void preserveOldRelevantValues(std::vector<ValueType> const& allValues, storm::storage::BitVector const& relevantValues, std::vector<ValueType>& oldValues) {
766 storm::utility::vector::selectVectorValues(oldValues, relevantValues, allValues);
767}
768
775template<typename ValueType, typename SolutionType>
776bool IterativeMinMaxLinearEquationSolver<ValueType, SolutionType>::solveEquationsIntervalIteration(Environment const& env, OptimizationDirection dir,
777 std::vector<SolutionType>& x,
778 std::vector<ValueType> const& b) const {
779 if constexpr (std::is_same_v<ValueType, storm::Interval>) {
780 STORM_LOG_THROW(false, storm::exceptions::NotImplementedException, "We did not implement intervaliteration for interval-based models");
781 return false;
782 } else {
783 setUpViOperator();
784 helper::IntervalIterationHelper<ValueType, false> iiHelper(viOperatorNontriv);
785 auto prec = storm::utility::convertNumber<ValueType>(env.solver().minMax().getPrecision());
786 auto lowerBoundsCallback = [&](std::vector<SolutionType>& vector) { this->createLowerBoundsVector(vector); };
787 auto upperBoundsCallback = [&](std::vector<SolutionType>& vector) { this->createUpperBoundsVector(vector); };
788
789 uint64_t numIterations{0};
790 auto iiCallback = [&](helper::IIData<ValueType> const& data) {
791 this->showProgressIterative(numIterations);
792 bool terminateEarly = this->hasCustomTerminationCondition() && this->getTerminationCondition().terminateNow(data.x, SolverGuarantee::LessOrEqual) &&
793 this->getTerminationCondition().terminateNow(data.y, SolverGuarantee::GreaterOrEqual);
794 return this->updateStatus(data.status, terminateEarly, numIterations, env.solver().minMax().getMaximalNumberOfIterations());
795 };
796 std::optional<storm::storage::BitVector> optionalRelevantValues;
797 if (this->hasRelevantValues()) {
798 optionalRelevantValues = this->getRelevantValues();
799 }
800 this->startMeasureProgress();
801 auto status = iiHelper.II(x, b, numIterations, env.solver().minMax().getRelativeTerminationCriterion(), prec, lowerBoundsCallback, upperBoundsCallback,
802 dir, iiCallback, optionalRelevantValues);
803 this->reportStatus(status, numIterations);
804
805 // If requested, we store the scheduler for retrieval.
806 if (this->isTrackSchedulerSet()) {
807 this->extractScheduler(x, b, dir, this->isUncertaintyRobust());
808 }
809
810 if (!this->isCachingEnabled()) {
811 clearCache();
812 }
813
814 return status == SolverStatus::Converged || status == SolverStatus::TerminatedEarly;
815 }
816}
817
818template<typename ValueType, typename SolutionType>
819bool IterativeMinMaxLinearEquationSolver<ValueType, SolutionType>::solveEquationsSoundValueIteration(Environment const& env, OptimizationDirection dir,
820 std::vector<SolutionType>& x,
821 std::vector<ValueType> const& b) const {
822 if constexpr (std::is_same_v<ValueType, storm::Interval>) {
823 STORM_LOG_THROW(false, storm::exceptions::NotImplementedException, "SoundVI does not handle interval-based models");
824 return false;
825 } else {
826 // Prepare the solution vectors and the helper.
827 assert(x.size() == this->A->getRowGroupCount());
828
829 std::optional<ValueType> lowerBound, upperBound;
830 if (this->hasLowerBound()) {
831 lowerBound = this->getLowerBound(true);
832 }
833 if (this->hasUpperBound()) {
834 upperBound = this->getUpperBound(true);
835 }
836
837 setUpViOperator();
838
839 auto precision = storm::utility::convertNumber<ValueType>(env.solver().minMax().getPrecision());
840 uint64_t numIterations{0};
841 auto sviCallback = [&](typename helper::SoundValueIterationHelper<ValueType, false>::SVIData const& current) {
842 this->showProgressIterative(numIterations);
843 return this->updateStatus(current.status,
844 this->hasCustomTerminationCondition() && current.checkCustomTerminationCondition(this->getTerminationCondition()),
845 numIterations, env.solver().minMax().getMaximalNumberOfIterations());
846 };
847 this->startMeasureProgress();
848 helper::SoundValueIterationHelper<ValueType, false> sviHelper(viOperatorNontriv);
849 std::optional<storm::storage::BitVector> optionalRelevantValues;
850 if (this->hasRelevantValues()) {
851 optionalRelevantValues = this->getRelevantValues();
852 }
853 auto status = sviHelper.SVI(x, b, numIterations, env.solver().minMax().getRelativeTerminationCriterion(), precision, dir, lowerBound, upperBound,
854 sviCallback, optionalRelevantValues);
855
856 // If requested, we store the scheduler for retrieval.
857 if (this->isTrackSchedulerSet()) {
858 this->extractScheduler(x, b, dir, this->isUncertaintyRobust());
859 }
860
861 this->reportStatus(status, numIterations);
862
863 if (!this->isCachingEnabled()) {
864 clearCache();
865 }
866
867 return status == SolverStatus::Converged || status == SolverStatus::TerminatedEarly;
868 }
869}
870
871template<typename ValueType, typename SolutionType>
872bool IterativeMinMaxLinearEquationSolver<ValueType, SolutionType>::solveEquationsViToPi(Environment const& env, OptimizationDirection dir,
873 std::vector<SolutionType>& x, std::vector<ValueType> const& b) const {
874 if constexpr (std::is_same_v<ValueType, storm::Interval>) {
875 STORM_LOG_THROW(false, storm::exceptions::NotImplementedException, "ViToPi does not handle interval-based models");
876 return false;
877 }
878 // First create an (inprecise) vi solver to get a good initial strategy for the (potentially precise) policy iteration solver.
879 std::vector<storm::storage::sparse::state_type> initialSched;
880 {
881 Environment viEnv = env;
882 viEnv.solver().minMax().setMethod(MinMaxMethod::ValueIteration);
883 viEnv.solver().setForceExact(false);
884 viEnv.solver().setForceSoundness(false);
885 auto impreciseSolver = GeneralMinMaxLinearEquationSolverFactory<double>().create(viEnv, this->A->template toValueType<double>());
886 impreciseSolver->setHasUniqueSolution(this->hasUniqueSolution());
887 impreciseSolver->setTrackScheduler(true);
888 if (this->hasInitialScheduler()) {
889 auto initSched = this->getInitialScheduler();
890 impreciseSolver->setInitialScheduler(std::move(initSched));
891 }
892 auto impreciseSolverReq = impreciseSolver->getRequirements(viEnv, dir);
893 STORM_LOG_THROW(!impreciseSolverReq.hasEnabledCriticalRequirement(), storm::exceptions::UnmetRequirementException,
894 "The value-iteration based solver has an unmet requirement: " << impreciseSolverReq.getEnabledRequirementsAsString());
895 impreciseSolver->setRequirementsChecked(true);
896 auto xVi = storm::utility::vector::convertNumericVector<double>(x);
897 auto bVi = storm::utility::vector::convertNumericVector<double>(b);
898 impreciseSolver->solveEquations(viEnv, dir, xVi, bVi);
899 initialSched = impreciseSolver->getSchedulerChoices();
900 }
901 STORM_LOG_INFO("Found initial policy using Value Iteration. Starting Policy iteration now.");
902 return performPolicyIteration(env, dir, x, b, std::move(initialSched));
903}
904
905template<typename ValueType, typename SolutionType>
906bool IterativeMinMaxLinearEquationSolver<ValueType, SolutionType>::solveEquationsRationalSearch(Environment const& env, OptimizationDirection dir,
907 std::vector<SolutionType>& x,
908 std::vector<ValueType> const& b) const {
909 if constexpr (std::is_same_v<ValueType, storm::Interval>) {
910 STORM_LOG_THROW(false, storm::exceptions::NotImplementedException, "Rational search does not handle interval-based models");
911 return false;
912 } else {
913 // Set up two value iteration operators. One for exact and one for imprecise computations
914 setUpViOperator();
915 std::shared_ptr<helper::ValueIterationOperator<storm::RationalNumber, false>> exactOp;
916 std::shared_ptr<helper::ValueIterationOperator<double, false>> impreciseOp;
917 std::function<bool(uint64_t, uint64_t)> fixedChoicesCallback;
918 if (this->choiceFixedForRowGroup) {
919 // Ignore those rows that are not selected
920 assert(this->initialScheduler);
921 fixedChoicesCallback = [&](uint64_t groupIndex, uint64_t localRowIndex) {
922 return this->choiceFixedForRowGroup->get(groupIndex) && this->initialScheduler->at(groupIndex) != localRowIndex;
923 };
924 }
925
926 if constexpr (std::is_same_v<ValueType, storm::RationalNumber>) {
927 exactOp = viOperatorNontriv;
928 impreciseOp = std::make_shared<helper::ValueIterationOperator<double, false>>();
929 impreciseOp->setMatrixBackwards(this->A->template toValueType<double>(), &this->A->getRowGroupIndices());
930 if (this->choiceFixedForRowGroup) {
931 impreciseOp->setIgnoredRows(true, fixedChoicesCallback);
932 }
933 } else if constexpr (std::is_same_v<ValueType, double>) {
934 impreciseOp = viOperatorNontriv;
935 exactOp = std::make_shared<helper::ValueIterationOperator<storm::RationalNumber, false>>();
936 exactOp->setMatrixBackwards(this->A->template toValueType<storm::RationalNumber>(), &this->A->getRowGroupIndices());
937 if (this->choiceFixedForRowGroup) {
938 exactOp->setIgnoredRows(true, fixedChoicesCallback);
939 }
940 }
941
943 uint64_t numIterations{0};
944 auto rsCallback = [&](SolverStatus const& current) {
945 this->showProgressIterative(numIterations);
946 return this->updateStatus(current, x, SolverGuarantee::None, numIterations, env.solver().minMax().getMaximalNumberOfIterations());
947 };
948 this->startMeasureProgress();
949 auto status = rsHelper.RS(x, b, numIterations, storm::utility::convertNumber<ValueType>(env.solver().minMax().getPrecision()), dir, rsCallback);
950
951 this->reportStatus(status, numIterations);
952
953 // If requested, we store the scheduler for retrieval.
954 if (this->isTrackSchedulerSet()) {
955 this->extractScheduler(x, b, dir, this->isUncertaintyRobust());
956 }
957
958 if (!this->isCachingEnabled()) {
959 clearCache();
960 }
961
962 return status == SolverStatus::Converged || status == SolverStatus::TerminatedEarly;
963 }
964}
965
966template<typename ValueType, typename SolutionType>
968 auxiliaryRowGroupVector.reset();
969 if (viOperatorTriv) {
970 viOperatorTriv.reset();
971 }
972 if (viOperatorNontriv) {
973 viOperatorNontriv.reset();
974 }
976}
977
981
982} // namespace solver
983} // namespace storm
SolverEnvironment & solver()
uint64_t const & getMaximalNumberOfIterations() const
storm::RationalNumber const & getPrecision() const
storm::solver::MinMaxMethod const & getMethod() const
bool const & getRelativeTerminationCriterion() const
std::optional< storm::RationalNumber > const & getUpperBoundGuessingFactor() const
OviSolverEnvironment const & ovi() const
MinMaxSolverEnvironment & minMax()
virtual void clearCache() const override
Clears the currently cached data that has been stored during previous calls of the solver.
IterativeMinMaxLinearEquationSolver(std::unique_ptr< LinearEquationSolverFactory< SolutionType > > &&linearEquationSolverFactory)
virtual bool internalSolveEquations(Environment const &env, OptimizationDirection dir, std::vector< SolutionType > &x, std::vector< ValueType > const &b) const override
virtual MinMaxLinearEquationSolverRequirements getRequirements(Environment const &env, boost::optional< storm::solver::OptimizationDirection > const &direction=boost::none, bool const &hasInitialScheduler=false) const override
Retrieves the requirements of this solver for solving equations with the current settings.
virtual void clearCache() const
Clears the currently cached data that has been stored during previous calls of the solver.
MinMaxLinearEquationSolverRequirements & requireBounds(bool critical=true)
MinMaxLinearEquationSolverRequirements & requireUniqueSolution(bool critical=true)
MinMaxLinearEquationSolverRequirements & requireLowerBounds(bool critical=true)
MinMaxLinearEquationSolverRequirements & requireValidInitialScheduler(bool critical=true)
MinMaxLinearEquationSolverRequirements & requireUpperBounds(bool critical=true)
Helper class to extract optimal scheduler choices from a MinMax equation system solution.
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.
A class that holds a possibly non-square matrix in the compressed row storage format.
void convertToEquationSystem()
Transforms the matrix into an equation system.
SparseMatrix selectRowsFromRowGroups(std::vector< index_type > const &rowGroupToRowIndexMapping, bool insertDiagonalEntries=true) const
Selects exactly one row from each row group of this matrix and returns the resulting matrix.
#define STORM_LOG_INFO(message)
Definition logging.h:24
#define STORM_LOG_WARN(message)
Definition logging.h:25
#define STORM_LOG_ERROR(message)
Definition logging.h:26
#define STORM_LOG_ASSERT(cond, message)
Definition macros.h:11
#define STORM_LOG_THROW(cond, exception, message)
Definition macros.h:30
SFTBDDChecker::ValueType ValueType
bool constexpr maximize(OptimizationDirection d)
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)
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
bool hasNonZeroEntry(std::vector< T > const &v)
Definition vector.h:1099
ValueType min(ValueType const &first, ValueType const &second)
bool isZero(ValueType const &a)
Definition constants.cpp:39