Storm
A Modern Probabilistic Model Checker
Loading...
Searching...
No Matches
HybridMdpPrctlHelper.cpp
Go to the documentation of this file.
2
4
10
12#include "storm/utility/graph.h"
13
15
22
25
27
30
31namespace storm {
32namespace modelchecker {
33namespace helper {
34
35template<typename ValueType>
37 boost::optional<SparseMdpEndComponentInformation<ValueType>> ecInformation;
38 boost::optional<std::vector<uint64_t>> initialScheduler;
40 boost::optional<std::vector<ValueType>> oneStepTargetProbabilities;
41};
42
43template<typename ValueType>
45 std::vector<ValueType> const& b) {
46 uint64_t numberOfMaybeStates = transitionMatrix.getRowGroupCount();
47 std::vector<uint64_t> result(numberOfMaybeStates);
48 storm::storage::BitVector targetStates(numberOfMaybeStates);
49
50 for (uint64_t state = 0; state < numberOfMaybeStates; ++state) {
51 // Record all states with non-zero probability of moving directly to the target states.
52 for (uint64_t row = transitionMatrix.getRowGroupIndices()[state]; row < transitionMatrix.getRowGroupIndices()[state + 1]; ++row) {
53 if (!storm::utility::isZero(b[row])) {
54 targetStates.set(state);
55 result[state] = row - transitionMatrix.getRowGroupIndices()[state];
56 }
57 }
58 }
59
60 if (!targetStates.full()) {
61 storm::storage::Scheduler<ValueType> validScheduler(numberOfMaybeStates);
62 storm::storage::SparseMatrix<ValueType> backwardTransitions = transitionMatrix.transpose(true);
63 storm::utility::graph::computeSchedulerProbGreater0E(transitionMatrix, backwardTransitions, storm::storage::BitVector(numberOfMaybeStates, true),
64 targetStates, validScheduler, boost::none);
65
66 for (uint64_t state = 0; state < numberOfMaybeStates; ++state) {
67 if (!targetStates.get(state)) {
68 result[state] = validScheduler.getChoice(state).getDeterministicChoice();
69 }
70 }
71 }
72
73 return result;
74}
75
76template<typename ValueType>
77void eliminateExtendedStatesFromExplicitRepresentation(std::pair<storm::storage::SparseMatrix<ValueType>, std::vector<ValueType>>& explicitRepresentation,
78 boost::optional<std::vector<uint64_t>>& scheduler, storm::storage::BitVector const& properMaybeStates) {
79 if (scheduler) {
80 // Eliminate superfluous entries from the scheduler.
81 uint64_t position = 0;
82 for (auto state : properMaybeStates) {
83 scheduler.get()[position] = scheduler.get()[state];
84 position++;
85 }
86 scheduler.get().resize(properMaybeStates.getNumberOfSetBits());
87 scheduler.get().shrink_to_fit();
88 }
89
90 // Treat the matrix.
91 explicitRepresentation.first = explicitRepresentation.first.getSubmatrix(true, properMaybeStates, properMaybeStates);
92}
93
94template<typename ValueType>
96 std::pair<storm::storage::SparseMatrix<ValueType>, std::vector<ValueType>>& explicitRepresentation,
97 SolverRequirementsData<ValueType>& solverRequirementsData, storm::storage::BitVector const& targetStates) {
98 // Get easier handles to the data.
99 auto& transitionMatrix = explicitRepresentation.first;
100 auto& oneStepProbabilities = explicitRepresentation.second;
101
102 bool doDecomposition = !solverRequirementsData.properMaybeStates.empty();
103
105 if (doDecomposition) {
106 auto backwardTransitions = transitionMatrix.transpose(true);
107 // Get the set of states that (under some scheduler) can stay in the set of maybestates forever
108 storm::storage::BitVector candidateStates =
109 storm::utility::graph::performProb0E(transitionMatrix, transitionMatrix.getRowGroupIndices(), backwardTransitions,
110 solverRequirementsData.properMaybeStates, ~solverRequirementsData.properMaybeStates);
111
112 doDecomposition = !candidateStates.empty();
113
114 if (doDecomposition) {
115 // If there are candidates, compute the states that are in MECs with zero reward.
116 endComponentDecomposition = storm::storage::MaximalEndComponentDecomposition<ValueType>(transitionMatrix, backwardTransitions, candidateStates);
117 }
118 }
119
120 // Only do more work if there are actually end-components.
121 if (doDecomposition && !endComponentDecomposition.empty()) {
122 STORM_LOG_DEBUG("Eliminating " << endComponentDecomposition.size() << " EC(s).");
123 std::vector<ValueType> subvector;
124 SparseMdpEndComponentInformation<ValueType>::eliminateEndComponents(endComponentDecomposition, transitionMatrix,
125 solverRequirementsData.properMaybeStates, &targetStates, nullptr, nullptr,
126 explicitRepresentation.first, &subvector, nullptr);
127 oneStepProbabilities = std::move(subvector);
128 } else {
129 STORM_LOG_DEBUG("Not eliminating ECs as there are none.");
130 oneStepProbabilities = explicitRepresentation.first.getConstrainedRowGroupSumVector(solverRequirementsData.properMaybeStates, targetStates);
131 eliminateExtendedStatesFromExplicitRepresentation(explicitRepresentation, solverRequirementsData.initialScheduler,
132 solverRequirementsData.properMaybeStates);
133 }
134}
135
136template<storm::dd::DdType DdType, typename ValueType>
139 storm::dd::Add<DdType, ValueType> const& transitionMatrix, storm::dd::Bdd<DdType> const& phiStates, storm::dd::Bdd<DdType> const& psiStates,
140 bool qualitative) {
141 // We need to identify the states which have to be taken out of the matrix, i.e. all states that have
142 // probability 0 and 1 of satisfying the until-formula.
143 storm::dd::Bdd<DdType> transitionMatrixBdd = transitionMatrix.notZero();
144 std::pair<storm::dd::Bdd<DdType>, storm::dd::Bdd<DdType>> statesWithProbability01;
145 if (dir == OptimizationDirection::Minimize) {
146 statesWithProbability01 = storm::utility::graph::performProb01Min(model, transitionMatrixBdd, phiStates, psiStates);
147 } else {
148 statesWithProbability01 = storm::utility::graph::performProb01Max(model, transitionMatrixBdd, phiStates, psiStates);
149 }
150 storm::dd::Bdd<DdType> maybeStates = !statesWithProbability01.first && !statesWithProbability01.second && model.getReachableStates();
151
152 STORM_LOG_INFO("Preprocessing: " << statesWithProbability01.first.getNonZeroCount() << " states with probability 1, "
153 << statesWithProbability01.second.getNonZeroCount() << " with probability 0 (" << maybeStates.getNonZeroCount()
154 << " states remaining).");
155
156 // Check whether we need to compute exact probabilities for some states.
157 if (qualitative) {
158 // Set the values for all maybe-states to 0.5 to indicate that their probability values are neither 0 nor 1.
160 model.getReachableStates(),
161 statesWithProbability01.second.template toAdd<ValueType>() +
162 maybeStates.template toAdd<ValueType>() * model.getManager().getConstant(storm::utility::convertNumber<ValueType>(0.5))));
163 } else {
164 // If there are maybe states, we need to solve an equation system.
165 if (!maybeStates.isZero()) {
166 // If we minimize, we know that the solution to the equation system has no end components
167 bool hasNoEndComponents = dir == storm::solver::OptimizationDirection::Minimize;
168 // Check for requirements of the solver early so we can adjust the maybe state computation accordingly.
171 linearEquationSolverFactory.getRequirements(env, hasNoEndComponents, hasNoEndComponents, dir);
172 storm::solver::MinMaxLinearEquationSolverRequirements clearedRequirements = requirements;
173 SolverRequirementsData<ValueType> solverRequirementsData;
174 bool extendMaybeStates = false;
175
176 if (clearedRequirements.hasEnabledRequirement()) {
177 if (clearedRequirements.uniqueSolution()) {
178 STORM_LOG_DEBUG("Scheduling EC elimination, because the solver requires a unique solution.");
179 extendMaybeStates = true;
180 clearedRequirements.clearUniqueSolution();
181 hasNoEndComponents = true;
182 }
183 if (clearedRequirements.validInitialScheduler() && !hasNoEndComponents) {
184 STORM_LOG_DEBUG("Scheduling valid scheduler computation, because the solver requires it.");
185 clearedRequirements.clearValidInitialScheduler();
186 }
187 clearedRequirements.clearBounds();
188 STORM_LOG_THROW(!clearedRequirements.hasEnabledCriticalRequirement(), storm::exceptions::UncheckedRequirementException,
189 "Solver requirements " + clearedRequirements.getEnabledRequirementsAsString() + " not checked.");
190 }
191
192 storm::dd::Bdd<DdType> extendedMaybeStates = maybeStates;
193 if (extendMaybeStates) {
194 // Extend the maybe states by all non-maybe states that can be reached from a maybe state within one step (they
195 // either are states with probability 0 or 1).
196 extendedMaybeStates |= maybeStates.relationalProduct(transitionMatrixBdd.existsAbstract(model.getNondeterminismVariables()),
197 model.getRowVariables(), model.getColumnVariables());
198 }
199
200 storm::utility::Stopwatch conversionWatch;
201
202 // Create the ODD for the translation between symbolic and explicit storage.
203 conversionWatch.start();
204 storm::dd::Odd odd = extendedMaybeStates.createOdd();
205 conversionWatch.stop();
206
207 // Convert the maybe states BDD to an ADD.
208 storm::dd::Add<DdType, ValueType> maybeStatesAdd = maybeStates.template toAdd<ValueType>();
209
210 // Start by cutting away all rows that do not belong to maybe states. Note that this leaves columns targeting
211 // non-maybe states in the matrix.
212 storm::dd::Add<DdType, ValueType> submatrix = transitionMatrix * maybeStatesAdd;
213
214 // If the maybe states were extended, we generate the explicit representation slightly differently.
215 std::pair<storm::storage::SparseMatrix<ValueType>, std::vector<ValueType>> explicitRepresentation;
216 if (extendMaybeStates) {
217 // Eliminate all transitions to non-extended-maybe states.
218 submatrix *= extendedMaybeStates.template toAdd<ValueType>().swapVariables(model.getRowColumnMetaVariablePairs());
219
220 // Only translate the matrix for now.
221 conversionWatch.start();
222 explicitRepresentation.first = submatrix.toMatrix(model.getNondeterminismVariables(), odd, odd);
223
224 // Get all original maybe states in the extended matrix.
225 solverRequirementsData.properMaybeStates = maybeStates.toVector(odd);
226
227 // Compute the target states within the set of extended maybe states.
228 storm::storage::BitVector targetStates = (extendedMaybeStates && statesWithProbability01.second).toVector(odd);
229 conversionWatch.stop();
230
231 // Eliminate the end components and remove the states that are not interesting (target or non-filter).
232 eliminateEndComponentsAndExtendedStatesUntilProbabilities(explicitRepresentation, solverRequirementsData, targetStates);
233
234 } else {
235 // Then compute the vector that contains the one-step probabilities to a state with probability 1 for all
236 // maybe states.
237 storm::dd::Add<DdType, ValueType> prob1StatesAsColumn = statesWithProbability01.second.template toAdd<ValueType>();
238 prob1StatesAsColumn = prob1StatesAsColumn.swapVariables(model.getRowColumnMetaVariablePairs());
239 storm::dd::Add<DdType, ValueType> subvector = submatrix * prob1StatesAsColumn;
240 subvector = subvector.sumAbstract(model.getColumnVariables());
241
242 // Finally cut away all columns targeting non-maybe states.
243 submatrix *= maybeStatesAdd.swapVariables(model.getRowColumnMetaVariablePairs());
244
245 // Translate the symbolic matrix/vector to their explicit representations and solve the equation system.
246 conversionWatch.start();
247 explicitRepresentation = submatrix.toMatrixVector(subvector, model.getNondeterminismVariables(), odd, odd);
248 conversionWatch.stop();
249
250 if (requirements.validInitialScheduler()) {
251 solverRequirementsData.initialScheduler =
252 computeValidInitialSchedulerForUntilProbabilities<ValueType>(explicitRepresentation.first, explicitRepresentation.second);
253 }
254 }
255
256 STORM_LOG_INFO("Converting symbolic matrix/vector to explicit representation done in " << conversionWatch.getTimeInMilliseconds() << "ms.");
257
258 // Create the solution vector.
259 std::vector<ValueType> x(explicitRepresentation.first.getRowGroupCount(), storm::utility::zero<ValueType>());
260
261 std::unique_ptr<storm::solver::MinMaxLinearEquationSolver<ValueType>> solver =
262 linearEquationSolverFactory.create(env, std::move(explicitRepresentation.first));
263
264 // Set whether the equation system will have a unique solution / no end components
265 solver->setHasUniqueSolution(hasNoEndComponents);
266 solver->setHasNoEndComponents(hasNoEndComponents);
267
268 if (solverRequirementsData.initialScheduler) {
269 solver->setInitialScheduler(std::move(solverRequirementsData.initialScheduler.get()));
270 }
271 solver->setBounds(storm::utility::zero<ValueType>(), storm::utility::one<ValueType>());
272 solver->setRequirementsChecked();
273 solver->solveEquations(env, dir, x, explicitRepresentation.second);
274
275 // If we included some target and non-filter states in the ODD, we need to expand the result from the solver.
276 if (requirements.uniqueSolution() && solverRequirementsData.ecInformation) {
277 std::vector<ValueType> extendedVector(solverRequirementsData.properMaybeStates.getNumberOfSetBits());
278 solverRequirementsData.ecInformation.get().setValues(extendedVector, solverRequirementsData.properMaybeStates, x);
279 x = std::move(extendedVector);
280 }
281
282 // If we extended the maybe states, we create a new ODD containing only the propery maybe states.
283 if (extendMaybeStates) {
284 odd = maybeStates.createOdd();
285 }
286
287 // Return a hybrid check result that stores the numerical values explicitly.
288 return std::unique_ptr<CheckResult>(new storm::modelchecker::HybridQuantitativeCheckResult<DdType, ValueType>(
289 model.getReachableStates(), model.getReachableStates() && !maybeStates, statesWithProbability01.second.template toAdd<ValueType>(), maybeStates,
290 odd, x));
291 } else {
293 model.getReachableStates(), statesWithProbability01.second.template toAdd<ValueType>()));
294 }
295 }
296}
297
298template<storm::dd::DdType DdType, typename ValueType>
301 storm::dd::Add<DdType, ValueType> const& transitionMatrix, storm::dd::Bdd<DdType> const& psiStates, bool qualitative) {
302 std::unique_ptr<CheckResult> result =
303 computeUntilProbabilities(env, dir == OptimizationDirection::Minimize ? OptimizationDirection::Maximize : OptimizationDirection::Maximize, model,
304 transitionMatrix, model.getReachableStates(), !psiStates && model.getReachableStates(), qualitative);
305 result->asQuantitativeCheckResult<ValueType>().oneMinus();
306 return result;
307}
308
309template<storm::dd::DdType DdType, typename ValueType>
315
316template<storm::dd::DdType DdType, typename ValueType>
319 storm::dd::Add<DdType, ValueType> const& transitionMatrix, storm::dd::Bdd<DdType> const& phiStates, storm::dd::Bdd<DdType> const& psiStates,
320 uint_fast64_t stepBound) {
321 // We need to identify the states which have to be taken out of the matrix, i.e. all states that have
322 // probability 0 or 1 of satisfying the until-formula.
323 storm::dd::Bdd<DdType> statesWithProbabilityGreater0;
324 if (dir == OptimizationDirection::Minimize) {
325 statesWithProbabilityGreater0 = storm::utility::graph::performProbGreater0A(model, transitionMatrix.notZero(), phiStates, psiStates);
326 } else {
327 statesWithProbabilityGreater0 = storm::utility::graph::performProbGreater0E(model, transitionMatrix.notZero(), phiStates, psiStates);
328 }
329 storm::dd::Bdd<DdType> maybeStates = statesWithProbabilityGreater0 && !psiStates && model.getReachableStates();
330
331 STORM_LOG_INFO("Preprocessing: " << statesWithProbabilityGreater0.getNonZeroCount() << " states with probability greater 0.");
332
333 // If there are maybe states, we need to perform matrix-vector multiplications.
334 if (!maybeStates.isZero()) {
335 storm::utility::Stopwatch conversionWatch;
336
337 // Create the ODD for the translation between symbolic and explicit storage.
338 conversionWatch.start();
339 storm::dd::Odd odd = maybeStates.createOdd();
340 conversionWatch.stop();
341
342 // Create the matrix and the vector for the equation system.
343 storm::dd::Add<DdType, ValueType> maybeStatesAdd = maybeStates.template toAdd<ValueType>();
344
345 // Start by cutting away all rows that do not belong to maybe states. Note that this leaves columns targeting
346 // non-maybe states in the matrix.
347 storm::dd::Add<DdType, ValueType> submatrix = transitionMatrix * maybeStatesAdd;
348
349 // Then compute the vector that contains the one-step probabilities to a state with probability 1 for all
350 // maybe states.
351 storm::dd::Add<DdType, ValueType> prob1StatesAsColumn = psiStates.template toAdd<ValueType>().swapVariables(model.getRowColumnMetaVariablePairs());
352 storm::dd::Add<DdType, ValueType> subvector = (submatrix * prob1StatesAsColumn).sumAbstract(model.getColumnVariables());
353
354 // Finally cut away all columns targeting non-maybe states.
355 submatrix *= maybeStatesAdd.swapVariables(model.getRowColumnMetaVariablePairs());
356
357 // Create the solution vector.
358 std::vector<ValueType> x(maybeStates.getNonZeroCount(), storm::utility::zero<ValueType>());
359
360 // Translate the symbolic matrix/vector to their explicit representations.
361 conversionWatch.start();
362 std::pair<storm::storage::SparseMatrix<ValueType>, std::vector<ValueType>> explicitRepresentation =
363 submatrix.toMatrixVector(subvector, model.getNondeterminismVariables(), odd, odd);
364 conversionWatch.stop();
365 STORM_LOG_INFO("Converting symbolic matrix/vector to explicit representation done in " << conversionWatch.getTimeInMilliseconds() << "ms.");
366
367 auto multiplier = storm::solver::MultiplierFactory<ValueType>().create(env, explicitRepresentation.first);
368 multiplier->repeatedMultiplyAndReduce(env, dir, x, &explicitRepresentation.second, stepBound);
369
370 // Return a hybrid check result that stores the numerical values explicitly.
371 return std::unique_ptr<CheckResult>(new storm::modelchecker::HybridQuantitativeCheckResult<DdType, ValueType>(
372 model.getReachableStates(), model.getReachableStates() && !maybeStates, psiStates.template toAdd<ValueType>(), maybeStates, odd, x));
373 } else {
374 return std::unique_ptr<CheckResult>(
375 new storm::modelchecker::SymbolicQuantitativeCheckResult<DdType, ValueType>(model.getReachableStates(), psiStates.template toAdd<ValueType>()));
376 }
377}
378
379template<storm::dd::DdType DdType, typename ValueType>
382 storm::dd::Add<DdType, ValueType> const& transitionMatrix, RewardModelType const& rewardModel, uint_fast64_t stepBound) {
383 // Only compute the result if the model has at least one reward this->getModel().
384 STORM_LOG_THROW(rewardModel.hasStateRewards(), storm::exceptions::InvalidPropertyException,
385 "Computing instantaneous rewards for a reward model that does not define any state-rewards. The result is trivially 0.");
386
387 storm::utility::Stopwatch conversionWatch;
388
389 // Create the ODD for the translation between symbolic and explicit storage.
391
392 // Translate the symbolic matrix to its explicit representations.
393 storm::storage::SparseMatrix<ValueType> explicitMatrix = transitionMatrix.toMatrix(model.getNondeterminismVariables(), odd, odd);
394
395 // Create the solution vector (and initialize it to the state rewards of the model).
396 std::vector<ValueType> x = rewardModel.getStateRewardVector().toVector(odd);
397 conversionWatch.stop();
398 STORM_LOG_INFO("Converting symbolic matrix/vector to explicit representation done in " << conversionWatch.getTimeInMilliseconds() << "ms.");
399
400 // Perform the matrix-vector multiplication.
401 auto multiplier = storm::solver::MultiplierFactory<ValueType>().create(env, explicitMatrix);
402 multiplier->repeatedMultiplyAndReduce(env, dir, x, nullptr, stepBound);
403
404 // Return a hybrid check result that stores the numerical values explicitly.
405 return std::unique_ptr<CheckResult>(new HybridQuantitativeCheckResult<DdType, ValueType>(
406 model.getReachableStates(), model.getManager().getBddZero(), model.getManager().template getAddZero<ValueType>(), model.getReachableStates(), odd, x));
407}
408
409template<storm::dd::DdType DdType, typename ValueType>
412 storm::dd::Add<DdType, ValueType> const& transitionMatrix, RewardModelType const& rewardModel, uint_fast64_t stepBound) {
413 // Only compute the result if the model has at least one reward this->getModel().
414 STORM_LOG_THROW(!rewardModel.empty(), storm::exceptions::InvalidPropertyException, "Missing reward model for formula. Skipping formula.");
415
416 // Compute the reward vector to add in each step based on the available reward models.
417 storm::dd::Add<DdType, ValueType> totalRewardVector = rewardModel.getTotalRewardVector(transitionMatrix, model.getColumnVariables());
418
419 // Create the solution vector.
420 std::vector<ValueType> x(model.getNumberOfStates(), storm::utility::zero<ValueType>());
421
422 storm::utility::Stopwatch conversionWatch(true);
423
424 // Create the ODD for the translation between symbolic and explicit storage.
426
427 // Translate the symbolic matrix/vector to their explicit representations.
428 std::pair<storm::storage::SparseMatrix<ValueType>, std::vector<ValueType>> explicitRepresentation =
429 transitionMatrix.toMatrixVector(totalRewardVector, model.getNondeterminismVariables(), odd, odd);
430 conversionWatch.stop();
431 STORM_LOG_INFO("Converting symbolic matrix/vector to explicit representation done in " << conversionWatch.getTimeInMilliseconds() << "ms.");
432
433 // Perform the matrix-vector multiplication.
434 auto multiplier = storm::solver::MultiplierFactory<ValueType>().create(env, explicitRepresentation.first);
435 multiplier->repeatedMultiplyAndReduce(env, dir, x, &explicitRepresentation.second, stepBound);
436
437 // Return a hybrid check result that stores the numerical values explicitly.
438 return std::unique_ptr<CheckResult>(new HybridQuantitativeCheckResult<DdType, ValueType>(
439 model.getReachableStates(), model.getManager().getBddZero(), model.getManager().template getAddZero<ValueType>(), model.getReachableStates(), odd, x));
440}
441
442template<typename ValueType>
444 storm::storage::BitVector targetStates(transitionMatrix.getRowGroupCount());
445
446 // A state is a target state iff its row group is empty.
447 for (uint64_t rowGroup = 0; rowGroup < transitionMatrix.getRowGroupCount(); ++rowGroup) {
448 if (transitionMatrix.getRowGroupIndices()[rowGroup] == transitionMatrix.getRowGroupIndices()[rowGroup + 1]) {
449 targetStates.set(rowGroup);
450 }
451 }
452
453 return targetStates;
454}
455
456template<typename ValueType>
458 storm::storage::BitVector const& properMaybeStates,
459 storm::storage::BitVector const& targetStates) {
460 uint64_t numberOfMaybeStates = transitionMatrix.getRowGroupCount();
461 std::vector<uint64_t> result(numberOfMaybeStates);
462
463 storm::storage::Scheduler<ValueType> validScheduler(numberOfMaybeStates);
464 storm::storage::SparseMatrix<ValueType> backwardTransitions = transitionMatrix.transpose(true);
465 storm::utility::graph::computeSchedulerProb1E(storm::storage::BitVector(numberOfMaybeStates, true), transitionMatrix, backwardTransitions,
466 properMaybeStates, targetStates, validScheduler);
467
468 for (uint64_t state = 0; state < numberOfMaybeStates; ++state) {
469 if (!targetStates.get(state)) {
470 result[state] = validScheduler.getChoice(state).getDeterministicChoice();
471 }
472 }
473
474 return result;
475}
476
477template<typename ValueType>
479 storm::storage::BitVector const& properMaybeStates,
480 storm::storage::BitVector const& targetStates) {
481 return extendedMatrix.getConstrainedRowGroupSumVector(properMaybeStates, targetStates);
482}
483
484template<typename ValueType>
486 std::pair<storm::storage::SparseMatrix<ValueType>, std::vector<ValueType>>& explicitRepresentation,
487 SolverRequirementsData<ValueType>& solverRequirementsData, storm::storage::BitVector const& targetStates, bool computeOneStepTargetProbabilities) {
488 // Get easier handles to the data.
489 auto& transitionMatrix = explicitRepresentation.first;
490 auto& rewardVector = explicitRepresentation.second;
491
492 // Start by computing the choices with reward 0, as we only want ECs within this fragment.
493 storm::storage::BitVector zeroRewardChoices(transitionMatrix.getRowCount());
494
495 uint64_t index = 0;
496 for (auto const& e : rewardVector) {
497 if (storm::utility::isZero(e)) {
498 zeroRewardChoices.set(index);
499 }
500 ++index;
501 }
502
503 // Compute the states that have some zero reward choice.
504 storm::storage::BitVector candidateStates(solverRequirementsData.properMaybeStates);
505 for (auto state : solverRequirementsData.properMaybeStates) {
506 bool keepState = false;
507
508 for (auto row = transitionMatrix.getRowGroupIndices()[state], rowEnd = transitionMatrix.getRowGroupIndices()[state + 1]; row < rowEnd; ++row) {
509 if (zeroRewardChoices.get(row)) {
510 keepState = true;
511 break;
512 }
513 }
514
515 if (!keepState) {
516 candidateStates.set(state, false);
517 }
518 }
519
520 bool doDecomposition = !candidateStates.empty();
521
523 if (doDecomposition) {
524 auto backwardTransitions = transitionMatrix.transpose(true);
525
526 // Only keep the candidate states that (under some scheduler) can stay in the set of candidates forever
527 candidateStates = storm::utility::graph::performProb0E(transitionMatrix, transitionMatrix.getRowGroupIndices(), backwardTransitions, candidateStates,
528 ~candidateStates);
529
530 doDecomposition = !candidateStates.empty();
531
532 if (doDecomposition) {
533 // If there are candidates, compute the states that are in MECs with zero reward.
534 endComponentDecomposition =
535 storm::storage::MaximalEndComponentDecomposition<ValueType>(transitionMatrix, backwardTransitions, candidateStates, zeroRewardChoices);
536 }
537 }
538
539 // Only do more work if there are actually end-components.
540 if (doDecomposition && !endComponentDecomposition.empty()) {
541 STORM_LOG_DEBUG("Eliminating " << endComponentDecomposition.size() << " EC(s).");
542 if (computeOneStepTargetProbabilities) {
543 solverRequirementsData.oneStepTargetProbabilities =
544 std::vector<ValueType>(solverRequirementsData.properMaybeStates.getNumberOfSetBits(), storm::utility::zero<ValueType>());
545 }
546
547 std::vector<ValueType> subvector(solverRequirementsData.properMaybeStates.getNumberOfSetBits(), storm::utility::zero<ValueType>());
549 endComponentDecomposition, transitionMatrix, solverRequirementsData.properMaybeStates, computeOneStepTargetProbabilities ? &targetStates : nullptr,
550 nullptr, &rewardVector, transitionMatrix, computeOneStepTargetProbabilities ? &solverRequirementsData.oneStepTargetProbabilities.get() : nullptr,
551 &subvector);
552 rewardVector = std::move(subvector);
553 } else {
554 STORM_LOG_DEBUG("Not eliminating ECs as there are none.");
555 if (computeOneStepTargetProbabilities) {
557 explicitRepresentation.first, solverRequirementsData.properMaybeStates, targetStates);
558 }
559 eliminateExtendedStatesFromExplicitRepresentation(explicitRepresentation, solverRequirementsData.initialScheduler,
560 solverRequirementsData.properMaybeStates);
561 }
562}
563
564template<typename ValueType>
566 storm::storage::SparseMatrix<ValueType> const& submatrix, std::vector<ValueType> const& choiceRewards,
567 std::vector<ValueType> const& oneStepTargetProbabilities) {
568 // For the min-case, we use DS-MPI, for the max-case variant 2 of the Baier et al. paper (CAV'17).
569 if (direction == storm::OptimizationDirection::Minimize) {
570 DsMpiMdpUpperRewardBoundsComputer<ValueType> dsmpi(submatrix, choiceRewards, oneStepTargetProbabilities);
571 solver.setUpperBounds(dsmpi.computeUpperBounds());
572 } else {
573 BaierUpperRewardBoundsComputer<ValueType> baier(submatrix, choiceRewards, oneStepTargetProbabilities);
574 solver.setUpperBound(baier.computeUpperBound());
575 }
576}
577
578template<storm::dd::DdType DdType, typename ValueType>
581 storm::dd::Add<DdType, ValueType> const& transitionMatrix, RewardModelType const& rewardModel, storm::dd::Bdd<DdType> const& targetStates,
582 bool qualitative) {
583 // Only compute the result if there is at least one reward model.
584 STORM_LOG_THROW(!rewardModel.empty(), storm::exceptions::InvalidPropertyException, "Missing reward model for formula. Skipping formula.");
585
586 // Determine which states have a reward of infinity by definition.
587 storm::dd::Bdd<DdType> infinityStates;
588 storm::dd::Bdd<DdType> transitionMatrixBdd = transitionMatrix.notZero();
589 if (dir == OptimizationDirection::Minimize) {
591 model, transitionMatrixBdd, model.getReachableStates(), targetStates,
592 storm::utility::graph::performProbGreater0E(model, transitionMatrixBdd, model.getReachableStates(), targetStates));
593 } else {
595 model, transitionMatrixBdd, targetStates,
596 storm::utility::graph::performProbGreater0A(model, transitionMatrixBdd, model.getReachableStates(), targetStates));
597 }
598 infinityStates = !infinityStates && model.getReachableStates();
599 storm::dd::Bdd<DdType> maybeStatesWithTargetStates = !infinityStates && model.getReachableStates();
600 storm::dd::Bdd<DdType> maybeStates = !targetStates && maybeStatesWithTargetStates;
601
602 STORM_LOG_INFO("Preprocessing: " << infinityStates.getNonZeroCount() << " states with reward infinity, " << targetStates.getNonZeroCount()
603 << " target states (" << maybeStates.getNonZeroCount() << " states remaining).");
604
605 // Check whether we need to compute exact rewards for some states.
606 if (qualitative) {
607 // Set the values for all maybe-states to 1 to indicate that their reward values
608 // are neither 0 nor infinity.
609 return std::unique_ptr<CheckResult>(new SymbolicQuantitativeCheckResult<DdType, ValueType>(
610 model.getReachableStates(),
611 infinityStates.ite(model.getManager().getConstant(storm::utility::infinity<ValueType>()), model.getManager().template getAddZero<ValueType>()) +
612 maybeStates.template toAdd<ValueType>() * model.getManager().getConstant(storm::utility::one<ValueType>())));
613 } else {
614 // If there are maybe states, we need to solve an equation system.
615 if (!maybeStates.isZero()) {
616 // If we maximize, we know that the solution to the equation system is unique.
617 bool hasNoEndComponents = dir == storm::solver::OptimizationDirection::Maximize;
618 bool hasUniqueSolution = hasNoEndComponents;
619 // Check for requirements of the solver this early so we can adapt the maybe states accordingly.
622 linearEquationSolverFactory.getRequirements(env, hasUniqueSolution, hasNoEndComponents, dir);
623 storm::solver::MinMaxLinearEquationSolverRequirements clearedRequirements = requirements;
624 bool extendMaybeStates = false;
625 if (clearedRequirements.hasEnabledRequirement()) {
626 if (clearedRequirements.uniqueSolution()) {
627 STORM_LOG_DEBUG("Scheduling EC elimination, because the solver requires it.");
628 extendMaybeStates = true;
629 clearedRequirements.clearUniqueSolution();
630 hasUniqueSolution = true;
631 // There might still be end components in which reward is collected.
632 }
633 if (clearedRequirements.validInitialScheduler()) {
634 STORM_LOG_DEBUG("Computing valid scheduler, because the solver requires it.");
635 extendMaybeStates = true;
636 clearedRequirements.clearValidInitialScheduler();
637 }
638 clearedRequirements.clearLowerBounds();
639 if (clearedRequirements.upperBounds()) {
640 STORM_LOG_DEBUG("Computing upper bounds, because the solver requires it.");
641 extendMaybeStates = true;
642 clearedRequirements.clearUpperBounds();
643 }
644 STORM_LOG_THROW(!clearedRequirements.hasEnabledCriticalRequirement(), storm::exceptions::UncheckedRequirementException,
645 "Solver requirements " + clearedRequirements.getEnabledRequirementsAsString() + " not checked.");
646 }
647
648 // Compute the set of maybe states that we are required to keep in the translation to explicit.
649 storm::dd::Bdd<DdType> requiredMaybeStates = extendMaybeStates ? maybeStatesWithTargetStates : maybeStates;
650
651 storm::utility::Stopwatch conversionWatch;
652
653 // Create the ODD for the translation between symbolic and explicit storage.
654 conversionWatch.start();
655 storm::dd::Odd odd = requiredMaybeStates.createOdd();
656 conversionWatch.stop();
657
658 // Create the matrix and the vector for the equation system.
659 storm::dd::Add<DdType, ValueType> maybeStatesAdd = maybeStates.template toAdd<ValueType>();
660
661 // Start by getting rid of
662 // (a) transitions from non-maybe states, and
663 // (b) the choices in the transition matrix that lead to a state that is neither a maybe state
664 // nor a target state ('infinity choices').
665 storm::dd::Add<DdType, ValueType> choiceFilterAdd =
666 maybeStatesAdd * (transitionMatrixBdd && maybeStatesWithTargetStates.renameVariables(model.getRowVariables(), model.getColumnVariables()))
667 .existsAbstract(model.getColumnVariables())
668 .template toAdd<ValueType>();
669 storm::dd::Add<DdType, ValueType> submatrix = transitionMatrix * choiceFilterAdd;
670
671 // Then compute the reward vector to use in the computation.
673 rewardModel.getTotalRewardVector(maybeStatesAdd, choiceFilterAdd, submatrix, model.getColumnVariables());
674
675 conversionWatch.start();
676 std::vector<uint_fast64_t> rowGroupSizes = (submatrix.notZero().existsAbstract(model.getColumnVariables()) || subvector.notZero())
677 .template toAdd<uint_fast64_t>()
678 .sumAbstract(model.getNondeterminismVariables())
679 .toVector(odd);
680 conversionWatch.stop();
681
682 // Finally cut away all columns targeting non-maybe states (or non-(maybe or target) states, respectively).
683 submatrix *= extendMaybeStates ? maybeStatesWithTargetStates.swapVariables(model.getRowColumnMetaVariablePairs()).template toAdd<ValueType>()
684 : maybeStatesAdd.swapVariables(model.getRowColumnMetaVariablePairs());
685
686 // Translate the symbolic matrix/vector to their explicit representations.
687 conversionWatch.start();
688 std::pair<storm::storage::SparseMatrix<ValueType>, std::vector<ValueType>> explicitRepresentation = submatrix.toMatrixVector(
689 std::move(rowGroupSizes), subvector, model.getRowVariables(), model.getColumnVariables(), model.getNondeterminismVariables(), odd, odd);
690 conversionWatch.stop();
691 STORM_LOG_INFO("Converting symbolic matrix/vector to explicit representation done in " << conversionWatch.getTimeInMilliseconds() << "ms.");
692
693 // Fulfill the solver's requirements.
694 SolverRequirementsData<ValueType> solverRequirementsData;
695 if (extendMaybeStates) {
697 solverRequirementsData.properMaybeStates = ~targetStates;
698
699 if (requirements.uniqueSolution()) {
700 STORM_LOG_THROW(!requirements.validInitialScheduler(), storm::exceptions::UncheckedRequirementException,
701 "The underlying solver requires a unique solution and an initial valid scheduler. This is currently not supported for "
702 "expected reward properties.");
703 // eliminate the end components with reward 0.
704 // Note that this may also compute the oneStepTargetProbabilities if upper bounds are required.
705 eliminateEndComponentsAndTargetStatesReachabilityRewards(explicitRepresentation, solverRequirementsData, targetStates,
706 requirements.upperBounds());
707 // The solution becomes unique after end components have been eliminated.
708 hasUniqueSolution = true;
709 } else {
710 if (requirements.validInitialScheduler()) {
711 // Compute a valid initial scheduler.
712 solverRequirementsData.initialScheduler = computeValidInitialSchedulerForReachabilityRewards<ValueType>(
713 explicitRepresentation.first, solverRequirementsData.properMaybeStates, targetStates);
714 }
715
716 if (requirements.upperBounds()) {
718 explicitRepresentation.first, solverRequirementsData.properMaybeStates, targetStates);
719 }
720
721 // Since we needed the transitions to target states to be translated as well for the computation
722 // of the scheduler, we have to get rid of them now.
723 eliminateExtendedStatesFromExplicitRepresentation(explicitRepresentation, solverRequirementsData.initialScheduler,
724 solverRequirementsData.properMaybeStates);
725 }
726 }
727
728 // Create the solution vector.
729 std::vector<ValueType> x(explicitRepresentation.first.getRowGroupCount(), storm::utility::zero<ValueType>());
730
731 // Now solve the resulting equation system.
732 std::unique_ptr<storm::solver::MinMaxLinearEquationSolver<ValueType>> solver = linearEquationSolverFactory.create(env);
733
734 // Set whether the equation system will have a unique solution / no end components
735 solver->setHasUniqueSolution(hasUniqueSolution);
736 solver->setHasNoEndComponents(hasNoEndComponents);
737
738 // If the solver requires upper bounds, compute them now.
739 if (requirements.upperBounds()) {
740 setUpperRewardBounds(*solver, dir, explicitRepresentation.first, explicitRepresentation.second,
741 solverRequirementsData.oneStepTargetProbabilities.get());
742 }
743
744 solver->setMatrix(std::move(explicitRepresentation.first));
745
746 // Move the scheduler to the solver.
747 if (solverRequirementsData.initialScheduler) {
748 solver->setInitialScheduler(std::move(solverRequirementsData.initialScheduler.get()));
749 }
750
751 solver->setLowerBound(storm::utility::zero<ValueType>());
752 solver->setRequirementsChecked();
753 solver->solveEquations(env, dir, x, explicitRepresentation.second);
754
755 // If we eliminated end components, we need to extend the solution vector.
756 if (requirements.uniqueSolution() && solverRequirementsData.ecInformation) {
757 std::vector<ValueType> extendedVector(solverRequirementsData.properMaybeStates.getNumberOfSetBits());
758 solverRequirementsData.ecInformation.get().setValues(extendedVector, solverRequirementsData.properMaybeStates, x);
759 x = std::move(extendedVector);
760 }
761
762 // If we extended the maybe states, we create a new ODD that only contains proper maybe states.
763 if (extendMaybeStates) {
764 odd = maybeStates.createOdd();
765 }
766
767 // Return a hybrid check result that stores the numerical values explicitly.
768 return std::unique_ptr<CheckResult>(new storm::modelchecker::HybridQuantitativeCheckResult<DdType, ValueType>(
769 model.getReachableStates(), model.getReachableStates() && !maybeStates,
770 infinityStates.ite(model.getManager().getConstant(storm::utility::infinity<ValueType>()), model.getManager().template getAddZero<ValueType>()),
771 maybeStates, odd, x));
772 } else {
774 model.getReachableStates(), infinityStates.ite(model.getManager().getConstant(storm::utility::infinity<ValueType>()),
775 model.getManager().template getAddZero<ValueType>())));
776 }
777 }
778}
779
780template<storm::dd::DdType DdType, typename ValueType>
783 storm::dd::Add<DdType, ValueType> const& transitionMatrix, storm::dd::Bdd<DdType> const& targetStates, bool qualitative) {
784 RewardModelType rewardModel(model.getManager().getConstant(storm::utility::one<ValueType>()), boost::none, boost::none);
785 return computeReachabilityRewards(env, dir, model, transitionMatrix, rewardModel, targetStates, qualitative);
786}
787
790
792
793} // namespace helper
794} // namespace modelchecker
795} // namespace storm
Add< LibraryType, ValueType > swapVariables(std::vector< std::pair< storm::expressions::Variable, storm::expressions::Variable > > const &metaVariablePairs) const
Swaps the given pairs of meta variables in the ADD.
Definition Add.cpp:292
storm::storage::SparseMatrix< ValueType > toMatrix() const
Converts the ADD to a (sparse) matrix.
Definition Add.cpp:634
std::pair< storm::storage::SparseMatrix< ValueType >, std::vector< ValueType > > toMatrixVector(storm::dd::Add< LibraryType, ValueType > const &vector, std::set< storm::expressions::Variable > const &groupMetaVariables, storm::dd::Odd const &rowOdd, storm::dd::Odd const &columnOdd) const
Converts the ADD to a row-grouped (sparse) matrix and the given vector to a row-grouped vector.
Definition Add.cpp:916
Add< LibraryType, ValueType > sumAbstract(std::set< storm::expressions::Variable > const &metaVariables) const
Sum-abstracts from the given meta variables.
Definition Add.cpp:178
Bdd< LibraryType > notZero() const
Computes a BDD that represents the function in which all assignments with a function value unequal to...
Definition Add.cpp:431
Bdd< LibraryType > existsAbstract(std::set< storm::expressions::Variable > const &metaVariables) const
Existentially abstracts from the given meta variables.
Definition Bdd.cpp:178
bool isZero() const
Retrieves whether this DD represents the constant zero function.
Definition Bdd.cpp:547
virtual uint_fast64_t getNonZeroCount() const override
Retrieves the number of encodings that are mapped to a non-zero value.
Definition Bdd.cpp:513
Bdd< LibraryType > swapVariables(std::vector< std::pair< storm::expressions::Variable, storm::expressions::Variable > > const &metaVariablePairs) const
Swaps the given pairs of meta variables in the BDD.
Definition Bdd.cpp:302
Odd createOdd() const
Creates an ODD based on the current BDD.
Definition Bdd.cpp:571
storm::storage::BitVector toVector(storm::dd::Odd const &rowOdd) const
Converts the BDD to a bit vector.
Definition Bdd.cpp:497
Bdd< LibraryType > ite(Bdd< LibraryType > const &thenBdd, Bdd< LibraryType > const &elseBdd) const
Performs an if-then-else with the given operands, i.e.
Definition Bdd.cpp:113
Bdd< LibraryType > renameVariables(std::set< storm::expressions::Variable > const &from, std::set< storm::expressions::Variable > const &to) const
Renames the given meta variables in the BDD.
Definition Bdd.cpp:347
Bdd< LibraryType > relationalProduct(Bdd< LibraryType > const &relation, std::set< storm::expressions::Variable > const &rowMetaVariables, std::set< storm::expressions::Variable > const &columnMetaVariables) const
Computes the relational product of the current BDD and the given BDD representing a relation.
Definition Bdd.cpp:220
ValueType computeUpperBound()
Computes an upper bound on the expected rewards.
std::vector< ValueType > computeUpperBounds()
Computes upper bounds on the expected rewards.
static std::unique_ptr< CheckResult > computeGloballyProbabilities(Environment const &env, OptimizationDirection dir, storm::models::symbolic::NondeterministicModel< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, storm::dd::Bdd< DdType > const &psiStates, bool qualitative)
static std::unique_ptr< CheckResult > computeNextProbabilities(Environment const &env, OptimizationDirection dir, storm::models::symbolic::NondeterministicModel< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, storm::dd::Bdd< DdType > const &nextStates)
static std::unique_ptr< CheckResult > computeInstantaneousRewards(Environment const &env, OptimizationDirection dir, storm::models::symbolic::NondeterministicModel< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, RewardModelType const &rewardModel, uint_fast64_t stepBound)
storm::models::symbolic::NondeterministicModel< DdType, ValueType >::RewardModelType RewardModelType
static std::unique_ptr< CheckResult > computeReachabilityTimes(Environment const &env, OptimizationDirection dir, storm::models::symbolic::NondeterministicModel< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, storm::dd::Bdd< DdType > const &targetStates, bool qualitative)
static std::unique_ptr< CheckResult > computeCumulativeRewards(Environment const &env, OptimizationDirection dir, storm::models::symbolic::NondeterministicModel< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, RewardModelType const &rewardModel, uint_fast64_t stepBound)
static std::unique_ptr< CheckResult > computeBoundedUntilProbabilities(Environment const &env, OptimizationDirection dir, storm::models::symbolic::NondeterministicModel< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, storm::dd::Bdd< DdType > const &phiStates, storm::dd::Bdd< DdType > const &psiStates, uint_fast64_t stepBound)
static std::unique_ptr< CheckResult > computeUntilProbabilities(Environment const &env, OptimizationDirection dir, storm::models::symbolic::NondeterministicModel< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, storm::dd::Bdd< DdType > const &phiStates, storm::dd::Bdd< DdType > const &psiStates, bool qualitative)
static std::unique_ptr< CheckResult > computeReachabilityRewards(Environment const &env, OptimizationDirection dir, storm::models::symbolic::NondeterministicModel< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, RewardModelType const &rewardModel, storm::dd::Bdd< DdType > const &targetStates, bool qualitative)
static SparseMdpEndComponentInformation< ValueType > eliminateEndComponents(storm::storage::MaximalEndComponentDecomposition< ValueType > const &endComponentDecomposition, storm::storage::SparseMatrix< ValueType > const &transitionMatrix, storm::storage::BitVector const &maybeStates, storm::storage::BitVector const *sumColumns, storm::storage::BitVector const *selectedChoices, std::vector< ValueType > const *summand, storm::storage::SparseMatrix< ValueType > &submatrix, std::vector< ValueType > *columnSumVector, std::vector< ValueType > *summandResultVector, bool gatherExitChoices=false)
static std::unique_ptr< CheckResult > computeNextProbabilities(Environment const &env, OptimizationDirection dir, storm::models::symbolic::NondeterministicModel< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, storm::dd::Bdd< DdType > const &nextStates)
storm::dd::DdManager< Type > & getManager() const
Retrieves the manager responsible for the DDs that represent this model.
Definition Model.cpp:92
std::set< storm::expressions::Variable > const & getColumnVariables() const
Retrieves the meta variables used to encode the columns of the transition matrix and the vector indic...
Definition Model.cpp:197
std::vector< std::pair< storm::expressions::Variable, storm::expressions::Variable > > const & getRowColumnMetaVariablePairs() const
Retrieves the pairs of row and column meta variables.
Definition Model.cpp:223
std::set< storm::expressions::Variable > const & getRowVariables() const
Retrieves the meta variables used to encode the rows of the transition matrix and the vector indices.
Definition Model.cpp:192
storm::dd::Bdd< Type > const & getReachableStates() const
Retrieves the reachable states of the model.
Definition Model.cpp:102
virtual uint_fast64_t getNumberOfStates() const override
Returns the number of states of the model.
Definition Model.cpp:77
Base class for all nondeterministic symbolic models.
virtual std::set< storm::expressions::Variable > const & getNondeterminismVariables() const override
Retrieves the meta variables used to encode the nondeterminism in the model.
void setUpperBound(ValueType const &value)
Sets an upper bound for the solution that can potentially be used by the solver.
void setUpperBounds(std::vector< ValueType > const &values)
Sets upper bounds for the solution that can potentially be used by the solver.
virtual std::unique_ptr< MinMaxLinearEquationSolver< ValueType, SolutionType > > create(Environment const &env) const override
MinMaxLinearEquationSolverRequirements getRequirements(Environment const &env, bool hasUniqueSolution=false, bool hasNoEndComponents=false, boost::optional< storm::solver::OptimizationDirection > const &direction=boost::none, bool hasInitialScheduler=false, bool trackScheduler=false) const
Retrieves the requirements of the solver that would be created when calling create() right now.
A class representing the interface that all min-max linear equation solvers shall implement.
std::string getEnabledRequirementsAsString() const
Returns a string that enumerates the enabled requirements.
std::unique_ptr< Multiplier< ValueType > > create(Environment const &env, storm::storage::SparseMatrix< ValueType > const &matrix)
A bit vector that is internally represented as a vector of 64-bit values.
Definition BitVector.h:18
bool full() const
Retrieves whether all bits are set in this bit vector.
bool empty() const
Retrieves whether no bits are set to true in this bit vector.
void set(uint_fast64_t index, bool value=true)
Sets the given truth value at the given index.
uint_fast64_t getNumberOfSetBits() const
Returns the number of bits that are set to true in this bit vector.
bool get(uint_fast64_t index) const
Retrieves the truth value of the bit at the given index and performs a bound check.
bool empty() const
Checks if the decomposition is empty.
std::size_t size() const
Retrieves the number of blocks of this decomposition.
This class represents the decomposition of a nondeterministic model into its maximal end components.
This class defines which action is chosen in a particular state of a non-deterministic model.
Definition Scheduler.h:18
SchedulerChoice< ValueType > const & getChoice(uint_fast64_t modelState, uint_fast64_t memoryState=0) const
Gets the choice defined by the scheduler for the given model and memory state.
Definition Scheduler.cpp:86
A class that holds a possibly non-square matrix in the compressed row storage format.
index_type getRowGroupCount() const
Returns the number of row groups in the matrix.
std::vector< index_type > const & getRowGroupIndices() const
Returns the grouping of rows of this matrix.
std::vector< value_type > getConstrainedRowGroupSumVector(storm::storage::BitVector const &rowGroupConstraint, storm::storage::BitVector const &columnConstraint) const
Computes a vector whose entries represent the sums of selected columns for all rows in selected row g...
storm::storage::SparseMatrix< value_type > transpose(bool joinGroups=false, bool keepZeros=false) const
Transposes the matrix.
A class that provides convenience operations to display run times.
Definition Stopwatch.h:14
MilisecondType getTimeInMilliseconds() const
Gets the measured time in milliseconds.
Definition Stopwatch.cpp:21
void start()
Start stopwatch (again) and start measuring time.
Definition Stopwatch.cpp:48
void stop()
Stop stopwatch and add measured time to total time.
Definition Stopwatch.cpp:42
#define STORM_LOG_INFO(message)
Definition logging.h:29
#define STORM_LOG_DEBUG(message)
Definition logging.h:23
#define STORM_LOG_THROW(cond, exception, message)
Definition macros.h:30
storm::storage::BitVector computeTargetStatesForReachabilityRewardsFromExplicitRepresentation(storm::storage::SparseMatrix< ValueType > const &transitionMatrix)
void eliminateEndComponentsAndExtendedStatesUntilProbabilities(std::pair< storm::storage::SparseMatrix< ValueType >, std::vector< ValueType > > &explicitRepresentation, SolverRequirementsData< ValueType > &solverRequirementsData, storm::storage::BitVector const &targetStates)
void eliminateEndComponentsAndTargetStatesReachabilityRewards(std::pair< storm::storage::SparseMatrix< ValueType >, std::vector< ValueType > > &explicitRepresentation, SolverRequirementsData< ValueType > &solverRequirementsData, storm::storage::BitVector const &targetStates, bool computeOneStepTargetProbabilities)
void setUpperRewardBounds(storm::solver::MinMaxLinearEquationSolver< ValueType > &solver, storm::OptimizationDirection const &direction, storm::storage::SparseMatrix< ValueType > const &submatrix, std::vector< ValueType > const &choiceRewards, std::vector< ValueType > const &oneStepTargetProbabilities)
std::vector< ValueType > computeOneStepTargetProbabilitiesFromExtendedExplicitRepresentation(storm::storage::SparseMatrix< ValueType > const &extendedMatrix, storm::storage::BitVector const &properMaybeStates, storm::storage::BitVector const &targetStates)
void eliminateExtendedStatesFromExplicitRepresentation(std::pair< storm::storage::SparseMatrix< ValueType >, std::vector< ValueType > > &explicitRepresentation, boost::optional< std::vector< uint64_t > > &scheduler, storm::storage::BitVector const &properMaybeStates)
std::vector< uint64_t > computeValidInitialSchedulerForReachabilityRewards(storm::storage::SparseMatrix< ValueType > const &transitionMatrix, storm::storage::BitVector const &properMaybeStates, storm::storage::BitVector const &targetStates)
std::vector< uint64_t > computeValidInitialSchedulerForUntilProbabilities(storm::storage::SparseMatrix< ValueType > const &transitionMatrix, std::vector< ValueType > const &b)
std::pair< storm::storage::BitVector, storm::storage::BitVector > performProb01Max(storm::storage::SparseMatrix< T > const &transitionMatrix, std::vector< uint_fast64_t > const &nondeterministicChoiceIndices, storm::storage::SparseMatrix< T > const &backwardTransitions, storm::storage::BitVector const &phiStates, storm::storage::BitVector const &psiStates)
Definition graph.cpp:835
void computeSchedulerProbGreater0E(storm::storage::SparseMatrix< T > const &transitionMatrix, storm::storage::SparseMatrix< T > const &backwardTransitions, storm::storage::BitVector const &phiStates, storm::storage::BitVector const &psiStates, storm::storage::Scheduler< SchedulerValueType > &scheduler, boost::optional< storm::storage::BitVector > const &rowFilter)
Computes a scheduler for the ProbGreater0E-States such that in the induced system the given psiStates...
Definition graph.cpp:542
storm::storage::BitVector performProbGreater0A(storm::storage::SparseMatrix< T > const &transitionMatrix, std::vector< uint_fast64_t > const &nondeterministicChoiceIndices, storm::storage::SparseMatrix< T > const &backwardTransitions, storm::storage::BitVector const &phiStates, storm::storage::BitVector const &psiStates, bool useStepBound, uint_fast64_t maximalSteps, boost::optional< storm::storage::BitVector > const &choiceConstraint)
Computes the sets of states that have probability greater 0 of satisfying phi until psi under any pos...
Definition graph.cpp:857
void computeSchedulerProb1E(storm::storage::BitVector const &prob1EStates, storm::storage::SparseMatrix< T > const &transitionMatrix, storm::storage::SparseMatrix< T > const &backwardTransitions, storm::storage::BitVector const &phiStates, storm::storage::BitVector const &psiStates, storm::storage::Scheduler< SchedulerValueType > &scheduler, boost::optional< storm::storage::BitVector > const &rowFilter)
Computes a scheduler for the given prob1EStates such that in the induced system the given psiStates a...
Definition graph.cpp:624
storm::storage::BitVector performProb1A(storm::models::sparse::NondeterministicModel< T, RM > const &model, storm::storage::SparseMatrix< T > const &backwardTransitions, storm::storage::BitVector const &phiStates, storm::storage::BitVector const &psiStates)
Computes the sets of states that have probability 1 of satisfying phi until psi under all possible re...
Definition graph.cpp:997
storm::storage::BitVector performProbGreater0E(storm::storage::SparseMatrix< T > const &backwardTransitions, storm::storage::BitVector const &phiStates, storm::storage::BitVector const &psiStates, bool useStepBound, uint_fast64_t maximalSteps)
Computes the sets of states that have probability greater 0 of satisfying phi until psi under at leas...
Definition graph.cpp:689
storm::storage::BitVector performProb0E(storm::models::sparse::NondeterministicModel< T, RM > const &model, storm::storage::SparseMatrix< T > const &backwardTransitions, storm::storage::BitVector const &phiStates, storm::storage::BitVector const &psiStates)
Computes the sets of states that have probability 0 of satisfying phi until psi under at least one po...
Definition graph.cpp:976
std::pair< storm::storage::BitVector, storm::storage::BitVector > performProb01Min(storm::storage::SparseMatrix< T > const &transitionMatrix, std::vector< uint_fast64_t > const &nondeterministicChoiceIndices, storm::storage::SparseMatrix< T > const &backwardTransitions, storm::storage::BitVector const &phiStates, storm::storage::BitVector const &psiStates)
Definition graph.cpp:1079
storm::storage::BitVector performProb1E(storm::storage::SparseMatrix< T > const &transitionMatrix, std::vector< uint_fast64_t > const &nondeterministicChoiceIndices, storm::storage::SparseMatrix< T > const &backwardTransitions, storm::storage::BitVector const &phiStates, storm::storage::BitVector const &psiStates, boost::optional< storm::storage::BitVector > const &choiceConstraint)
Computes the sets of states that have probability 1 of satisfying phi until psi under at least one po...
Definition graph.cpp:757
bool isZero(ValueType const &a)
Definition constants.cpp:41
LabParser.cpp.
Definition cli.cpp:18
boost::optional< std::vector< uint64_t > > initialScheduler
boost::optional< SparseMdpEndComponentInformation< ValueType > > ecInformation
boost::optional< std::vector< ValueType > > oneStepTargetProbabilities