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