Storm
A Modern Probabilistic Model Checker
Loading...
Searching...
No Matches
SparseDtmcEliminationModelChecker.cpp
Go to the documentation of this file.
2
3#include <algorithm>
4#include <chrono>
5#include <random>
6
26#include "storm/utility/graph.h"
30
31namespace storm {
32namespace modelchecker {
33
35
36template<typename SparseDtmcModelType>
41
42template<typename SparseDtmcModelType>
52
53template<typename SparseDtmcModelType>
56 storm::logic::StateFormula const& stateFormula = checkTask.getFormula();
57 std::unique_ptr<CheckResult> subResultPointer = this->check(stateFormula);
58 storm::storage::BitVector const& psiStates = subResultPointer->asExplicitQualitativeCheckResult().getTruthValuesVector();
59
60 storm::storage::SparseMatrix<ValueType> const& transitionMatrix = this->getModel().getTransitionMatrix();
61 uint_fast64_t numberOfStates = transitionMatrix.getRowCount();
62 if (psiStates.empty()) {
63 return std::unique_ptr<CheckResult>(
64 new ExplicitQuantitativeCheckResult<ValueType>(std::vector<ValueType>(numberOfStates, storm::utility::zero<ValueType>())));
65 }
66 if (psiStates.full()) {
67 return std::unique_ptr<CheckResult>(
68 new ExplicitQuantitativeCheckResult<ValueType>(std::vector<ValueType>(numberOfStates, storm::utility::one<ValueType>())));
69 }
70
71 storm::storage::BitVector const& initialStates = this->getModel().getInitialStates();
72 STORM_LOG_THROW(initialStates.getNumberOfSetBits() == 1, storm::exceptions::IllegalArgumentException,
73 "Input model is required to have exactly one initial state.");
74 STORM_LOG_THROW(checkTask.isOnlyInitialStatesRelevantSet(), storm::exceptions::IllegalArgumentException,
75 "Cannot compute long-run probabilities for all states.");
76
77 storm::storage::SparseMatrix<ValueType> backwardTransitions = this->getModel().getBackwardTransitions();
78 storm::storage::BitVector maybeStates =
79 storm::utility::graph::performProbGreater0(backwardTransitions, storm::storage::BitVector(transitionMatrix.getRowCount(), true), psiStates);
80
81 std::vector<ValueType> result(transitionMatrix.getRowCount(), storm::utility::zero<ValueType>());
82
83 // Determine whether we need to perform some further computation.
84 bool furtherComputationNeeded = true;
85 if (checkTask.isOnlyInitialStatesRelevantSet() && initialStates.isDisjointFrom(maybeStates)) {
86 STORM_LOG_DEBUG("The long-run probability for all initial states was found in a preprocessing step.");
87 furtherComputationNeeded = false;
88 }
89 if (maybeStates.empty()) {
90 STORM_LOG_DEBUG("The long-run probability for all states was found in a preprocessing step.");
91 furtherComputationNeeded = false;
92 }
93
94 if (furtherComputationNeeded) {
95 if (checkTask.isOnlyInitialStatesRelevantSet()) {
96 // Determine the set of states that is reachable from the initial state without jumping over a target state.
98 transitionMatrix, initialStates, storm::storage::BitVector(numberOfStates, true), storm::storage::BitVector(numberOfStates, false));
99
100 // Subtract from the maybe states the set of states that is not reachable (on a path from the initial to a target state).
101 maybeStates &= reachableStates;
102 }
103
104 std::vector<ValueType> stateValues(maybeStates.size(), storm::utility::zero<ValueType>());
105 storm::utility::vector::setVectorValues(stateValues, psiStates, storm::utility::one<ValueType>());
106 result =
107 computeLongRunValues(transitionMatrix, backwardTransitions, initialStates, maybeStates, checkTask.isOnlyInitialStatesRelevantSet(), stateValues);
108 }
109
110 // Construct check result based on whether we have computed values for all states or just the initial states.
111 std::unique_ptr<CheckResult> checkResult(new ExplicitQuantitativeCheckResult<ValueType>(result));
112 if (checkTask.isOnlyInitialStatesRelevantSet()) {
113 // If we computed the results for the initial states only, we need to filter the result to only
114 // communicate these results.
115 checkResult->filter(ExplicitQualitativeCheckResult(initialStates));
116 }
117 return checkResult;
118}
119
120template<typename SparseDtmcModelType>
123 // Do some sanity checks to establish some required properties.
124 RewardModelType const& rewardModel = this->getModel().getRewardModel(checkTask.isRewardModelSet() ? checkTask.getRewardModel() : "");
125 STORM_LOG_THROW(!rewardModel.empty(), storm::exceptions::IllegalArgumentException, "Input model does not have a reward model.");
126
127 storm::storage::BitVector const& initialStates = this->getModel().getInitialStates();
128 STORM_LOG_THROW(initialStates.getNumberOfSetBits() == 1, storm::exceptions::IllegalArgumentException,
129 "Input model is required to have exactly one initial state.");
130 STORM_LOG_THROW(checkTask.isOnlyInitialStatesRelevantSet(), storm::exceptions::IllegalArgumentException,
131 "Cannot compute long-run probabilities for all states.");
132
133 storm::storage::SparseMatrix<ValueType> const& transitionMatrix = this->getModel().getTransitionMatrix();
134 uint_fast64_t numberOfStates = transitionMatrix.getRowCount();
135
136 // Get the state-reward values from the reward model.
137 std::vector<ValueType> stateRewardValues = rewardModel.getTotalRewardVector(this->getModel().getTransitionMatrix());
138
139 storm::storage::BitVector maybeStates(stateRewardValues.size());
140 uint_fast64_t index = 0;
141 for (auto const& value : stateRewardValues) {
142 if (value != storm::utility::zero<ValueType>()) {
143 maybeStates.set(index, true);
144 }
145 ++index;
146 }
147
148 storm::storage::SparseMatrix<ValueType> backwardTransitions = this->getModel().getBackwardTransitions();
149
150 storm::storage::BitVector allStates(numberOfStates, true);
151 maybeStates = storm::utility::graph::performProbGreater0(backwardTransitions, allStates, maybeStates);
152
153 std::vector<ValueType> result(numberOfStates, storm::utility::zero<ValueType>());
154
155 // Determine whether we need to perform some further computation.
156 bool furtherComputationNeeded = true;
157 if (checkTask.isOnlyInitialStatesRelevantSet() && initialStates.isDisjointFrom(maybeStates)) {
158 furtherComputationNeeded = false;
159 }
160
161 if (furtherComputationNeeded) {
162 if (checkTask.isOnlyInitialStatesRelevantSet()) {
163 // Determine the set of states that is reachable from the initial state without jumping over a target state.
165 transitionMatrix, initialStates, storm::storage::BitVector(numberOfStates, true), storm::storage::BitVector(numberOfStates, false));
166
167 // Subtract from the maybe states the set of states that is not reachable (on a path from the initial to a target state).
168 maybeStates &= reachableStates;
169 }
170
171 result = computeLongRunValues(transitionMatrix, backwardTransitions, initialStates, maybeStates, checkTask.isOnlyInitialStatesRelevantSet(),
172 stateRewardValues);
173 }
174
175 // Construct check result based on whether we have computed values for all states or just the initial states.
176 std::unique_ptr<CheckResult> checkResult(new ExplicitQuantitativeCheckResult<ValueType>(result));
177 if (checkTask.isOnlyInitialStatesRelevantSet()) {
178 // If we computed the results for the initial states only, we need to filter the result to only
179 // communicate these results.
180 checkResult->filter(ExplicitQualitativeCheckResult(initialStates));
181 }
182 return checkResult;
183}
184
185template<typename SparseDtmcModelType>
186std::vector<typename SparseDtmcEliminationModelChecker<SparseDtmcModelType>::ValueType>
188 storm::storage::SparseMatrix<ValueType> const& backwardTransitions,
189 storm::storage::BitVector const& initialStates,
190 storm::storage::BitVector const& maybeStates,
191 bool computeResultsForInitialStatesOnly, std::vector<ValueType>& stateValues) {
192 std::chrono::high_resolution_clock::time_point totalTimeStart = std::chrono::high_resolution_clock::now();
193
194 // Start by decomposing the DTMC into its BSCCs.
195 std::chrono::high_resolution_clock::time_point sccDecompositionStart = std::chrono::high_resolution_clock::now();
197 transitionMatrix, storm::storage::StronglyConnectedComponentDecompositionOptions().onlyBottomSccs());
198 auto sccDecompositionEnd = std::chrono::high_resolution_clock::now();
199
200 std::chrono::high_resolution_clock::time_point conversionStart = std::chrono::high_resolution_clock::now();
201
202 // Then, we convert the reduced matrix to a more flexible format to be able to perform state elimination more easily.
203 storm::storage::FlexibleSparseMatrix<ValueType> flexibleMatrix(transitionMatrix);
204 flexibleMatrix.filterEntries(maybeStates, maybeStates);
205 storm::storage::FlexibleSparseMatrix<ValueType> flexibleBackwardTransitions(backwardTransitions);
206 flexibleBackwardTransitions.filterEntries(maybeStates, maybeStates);
207 auto conversionEnd = std::chrono::high_resolution_clock::now();
208
209 std::chrono::high_resolution_clock::time_point modelCheckingStart = std::chrono::high_resolution_clock::now();
210
213 boost::optional<std::vector<uint_fast64_t>> distanceBasedPriorities;
215 distanceBasedPriorities = getDistanceBasedPriorities(transitionMatrix, backwardTransitions, initialStates, stateValues,
217 }
218
219 uint_fast64_t numberOfStates = transitionMatrix.getRowCount();
220 storm::storage::BitVector regularStatesInBsccs(numberOfStates);
221 storm::storage::BitVector relevantBsccs(bsccDecomposition.size());
222 storm::storage::BitVector bsccRepresentativesAsBitVector(numberOfStates);
223 std::vector<storm::storage::sparse::state_type> bsccRepresentatives;
224 uint_fast64_t currentIndex = 0;
225 for (auto const& bscc : bsccDecomposition) {
226 // Since all states in an SCC can reach all other states, we only need to check whether an arbitrary
227 // state is a maybe state.
228 if (maybeStates.get(*bscc.cbegin())) {
229 relevantBsccs.set(currentIndex);
230 bsccRepresentatives.push_back(*bscc.cbegin());
231 bsccRepresentativesAsBitVector.set(*bscc.cbegin(), true);
232 for (auto const& state : bscc) {
233 regularStatesInBsccs.set(state, true);
234 }
235 }
236 ++currentIndex;
237 }
238 regularStatesInBsccs &= ~bsccRepresentativesAsBitVector;
239
240 // Compute the average time to stay in each state for all states in BSCCs.
241 std::vector<ValueType> averageTimeInStates(stateValues.size(), storm::utility::one<ValueType>());
242
243 // First, we eliminate all states in BSCCs (except for the representative states).
244 std::shared_ptr<StatePriorityQueue> priorityQueue =
245 createStatePriorityQueue(distanceBasedPriorities, flexibleMatrix, flexibleBackwardTransitions, stateValues, regularStatesInBsccs);
246 storm::solver::stateelimination::MultiValueStateEliminator<ValueType> stateEliminator(flexibleMatrix, flexibleBackwardTransitions, priorityQueue,
247 stateValues, averageTimeInStates);
248
249 while (priorityQueue->hasNext()) {
250 storm::storage::sparse::state_type state = priorityQueue->pop();
251 stateEliminator.eliminateState(state, true);
252#ifdef STORM_DEV
253 STORM_LOG_ASSERT(checkConsistent(flexibleMatrix, flexibleBackwardTransitions), "The forward and backward transition matrices became inconsistent.");
254#endif
255 }
256
257 // Now, we set the values of all states in BSCCs to that of the representative value (and clear the
258 // transitions of the representative states while doing so).
259 auto representativeIt = bsccRepresentatives.begin();
260 for (auto sccIndex : relevantBsccs) {
261 // We only need to set the values for all states of the BSCC if we are not computing the values for the
262 // initial states only.
263 ValueType bsccValue = stateValues[*representativeIt] / averageTimeInStates[*representativeIt];
264 auto const& bscc = bsccDecomposition[sccIndex];
265 if (!computeResultsForInitialStatesOnly) {
266 for (auto const& state : bscc) {
267 stateValues[state] = bsccValue;
268 }
269 } else {
270 for (auto const& state : bscc) {
271 stateValues[state] = storm::utility::zero<ValueType>();
272 }
273 stateValues[*representativeIt] = bsccValue;
274 }
275
276 FlexibleRowType& representativeForwardRow = flexibleMatrix.getRow(*representativeIt);
277 representativeForwardRow.clear();
278 representativeForwardRow.shrink_to_fit();
279
280 FlexibleRowType& representativeBackwardRow = flexibleBackwardTransitions.getRow(*representativeIt);
281 auto it = representativeBackwardRow.begin(), ite = representativeBackwardRow.end();
282 for (; it != ite; ++it) {
283 if (it->getColumn() == *representativeIt) {
284 break;
285 }
286 }
287 representativeBackwardRow.erase(it);
288
289 ++representativeIt;
290 }
291
292 // If there are states remaining that are not in BSCCs, we need to eliminate them now.
293 storm::storage::BitVector remainingStates = maybeStates & ~regularStatesInBsccs;
294
295 // Set the value initial value of all states not in a BSCC to zero, because a) any previous value would
296 // incorrectly influence the result and b) the value have been erroneously changed for the predecessors of
297 // BSCCs by the previous state elimination.
298 for (auto state : remainingStates) {
299 if (!bsccRepresentativesAsBitVector.get(state)) {
300 stateValues[state] = storm::utility::zero<ValueType>();
301 }
302 }
303
304 // We only need to eliminate the remaining states if there was some BSCC that has a non-zero value, i.e.
305 // that consists of maybe states.
306 if (!relevantBsccs.empty()) {
307 performOrdinaryStateElimination(flexibleMatrix, flexibleBackwardTransitions, remainingStates, initialStates, computeResultsForInitialStatesOnly,
308 stateValues, distanceBasedPriorities);
309 }
310
311 std::chrono::high_resolution_clock::time_point modelCheckingEnd = std::chrono::high_resolution_clock::now();
312 std::chrono::high_resolution_clock::time_point totalTimeEnd = std::chrono::high_resolution_clock::now();
313
315 std::chrono::high_resolution_clock::duration sccDecompositionTime = sccDecompositionEnd - sccDecompositionStart;
316 std::chrono::milliseconds sccDecompositionTimeInMilliseconds = std::chrono::duration_cast<std::chrono::milliseconds>(sccDecompositionTime);
317 std::chrono::high_resolution_clock::duration conversionTime = conversionEnd - conversionStart;
318 std::chrono::milliseconds conversionTimeInMilliseconds = std::chrono::duration_cast<std::chrono::milliseconds>(conversionTime);
319 std::chrono::high_resolution_clock::duration modelCheckingTime = modelCheckingEnd - modelCheckingStart;
320 std::chrono::milliseconds modelCheckingTimeInMilliseconds = std::chrono::duration_cast<std::chrono::milliseconds>(modelCheckingTime);
321 std::chrono::high_resolution_clock::duration totalTime = totalTimeEnd - totalTimeStart;
322 std::chrono::milliseconds totalTimeInMilliseconds = std::chrono::duration_cast<std::chrono::milliseconds>(totalTime);
323
325 STORM_PRINT_AND_LOG("Time breakdown:\n");
326 STORM_PRINT_AND_LOG(" * time for SCC decomposition: " << sccDecompositionTimeInMilliseconds.count() << "ms\n");
327 STORM_PRINT_AND_LOG(" * time for conversion: " << conversionTimeInMilliseconds.count() << "ms\n");
328 STORM_PRINT_AND_LOG(" * time for checking: " << modelCheckingTimeInMilliseconds.count() << "ms\n");
329 STORM_PRINT_AND_LOG("------------------------------------------\n");
330 STORM_PRINT_AND_LOG(" * total time: " << totalTimeInMilliseconds.count() << "ms\n");
331 }
332
333 // Now, we return the value for the only initial state.
334 STORM_LOG_DEBUG("Simplifying and returning result.");
335 for (auto& value : stateValues) {
336 value = storm::utility::simplify(value);
337 }
338 return stateValues;
339}
340
341template<typename SparseDtmcModelType>
344 storm::logic::BoundedUntilFormula const& pathFormula = checkTask.getFormula();
345
346 STORM_LOG_THROW(!pathFormula.hasLowerBound() && pathFormula.hasUpperBound(), storm::exceptions::InvalidPropertyException,
347 "Formula needs to have single upper time bound.");
348 STORM_LOG_THROW(pathFormula.hasIntegerUpperBound(), storm::exceptions::InvalidPropertyException, "Formula needs to have discrete upper time bound.");
349
350 // Retrieve the appropriate bitvectors by model checking the subformulas.
351 std::unique_ptr<CheckResult> leftResultPointer = this->check(pathFormula.getLeftSubformula());
352 std::unique_ptr<CheckResult> rightResultPointer = this->check(pathFormula.getRightSubformula());
353 storm::storage::BitVector const& phiStates = leftResultPointer->asExplicitQualitativeCheckResult().getTruthValuesVector();
354 storm::storage::BitVector const& psiStates = rightResultPointer->asExplicitQualitativeCheckResult().getTruthValuesVector();
355
356 // Start by determining the states that have a non-zero probability of reaching the target states within the
357 // time bound.
359 this->getModel().getBackwardTransitions(), phiStates, psiStates, true, pathFormula.getUpperBound<uint64_t>());
360 statesWithProbabilityGreater0 &= ~psiStates;
361
362 // Determine whether we need to perform some further computation.
363 bool furtherComputationNeeded = true;
364 if (checkTask.isOnlyInitialStatesRelevantSet() && this->getModel().getInitialStates().isDisjointFrom(statesWithProbabilityGreater0)) {
365 STORM_LOG_DEBUG("The probability for all initial states was found in a preprocessing step.");
366 furtherComputationNeeded = false;
367 } else if (statesWithProbabilityGreater0.empty()) {
368 STORM_LOG_DEBUG("The probability for all states was found in a preprocessing step.");
369 furtherComputationNeeded = false;
370 }
371
372 storm::storage::SparseMatrix<ValueType> const& transitionMatrix = this->getModel().getTransitionMatrix();
373 storm::storage::BitVector const& initialStates = this->getModel().getInitialStates();
374
375 std::vector<ValueType> result(transitionMatrix.getRowCount(), storm::utility::zero<ValueType>());
376
377 if (furtherComputationNeeded) {
378 uint_fast64_t timeBound = pathFormula.getUpperBound<uint64_t>();
379
380 if (checkTask.isOnlyInitialStatesRelevantSet()) {
381 // Determine the set of states that is reachable from the initial state without jumping over a target state.
382 storm::storage::BitVector reachableStates =
383 storm::utility::graph::getReachableStates(transitionMatrix, initialStates, phiStates, psiStates, true, timeBound);
384
385 // Subtract from the maybe states the set of states that is not reachable (on a path from the initial to a target state).
386 statesWithProbabilityGreater0 &= reachableStates;
387 }
388
389 // We then build the submatrix that only has the transitions of the maybe states.
391 transitionMatrix.getSubmatrix(true, statesWithProbabilityGreater0, statesWithProbabilityGreater0, true);
392
393 std::vector<uint_fast64_t> distancesFromInitialStates;
394 storm::storage::BitVector relevantStates;
395 if (checkTask.isOnlyInitialStatesRelevantSet()) {
396 // Determine the set of initial states of the sub-model.
397 storm::storage::BitVector subInitialStates = this->getModel().getInitialStates() % statesWithProbabilityGreater0;
398
399 // Precompute the distances of the relevant states to the initial states.
400 distancesFromInitialStates = storm::utility::graph::getDistances(submatrix, subInitialStates, statesWithProbabilityGreater0);
401
402 // Set all states to be relevant for later use.
403 relevantStates = storm::storage::BitVector(statesWithProbabilityGreater0.getNumberOfSetBits(), true);
404 }
405
406 // Create the vector of one-step probabilities to go to target states.
407 std::vector<ValueType> b = transitionMatrix.getConstrainedRowSumVector(statesWithProbabilityGreater0, psiStates);
408
409 // Create the vector with which to multiply.
410 std::vector<ValueType> subresult(b);
411 std::vector<ValueType> tmp(subresult.size());
412
413 // Subtract one from the time bound because initializing the sub-result to b already accounts for one step.
414 --timeBound;
415
416 // Perform matrix-vector multiplications until the time-bound is met.
417 for (uint_fast64_t timeStep = 0; timeStep < timeBound; ++timeStep) {
418 submatrix.multiplyWithVector(subresult, tmp);
419 storm::utility::vector::addVectors(tmp, b, subresult);
420
421 // If we are computing the results for the initial states only, we can use the minimal distance from
422 // each state to the initial states to determine whether we still need to consider the values for
423 // these states. If not, we can null-out all their probabilities.
424 if (checkTask.isOnlyInitialStatesRelevantSet()) {
425 for (auto state : relevantStates) {
426 if (distancesFromInitialStates[state] > (timeBound - timeStep)) {
427 for (auto& element : submatrix.getRow(state)) {
428 element.setValue(storm::utility::zero<ValueType>());
429 }
430 b[state] = storm::utility::zero<ValueType>();
431 relevantStates.set(state, false);
432 }
433 }
434 }
435 }
436
437 // Set the values of the resulting vector accordingly.
438 storm::utility::vector::setVectorValues(result, statesWithProbabilityGreater0, subresult);
439 }
440 storm::utility::vector::setVectorValues<ValueType>(result, psiStates, storm::utility::one<ValueType>());
441
442 // Construct check result based on whether we have computed values for all states or just the initial states.
443 std::unique_ptr<CheckResult> checkResult(new ExplicitQuantitativeCheckResult<ValueType>(result));
444 if (checkTask.isOnlyInitialStatesRelevantSet()) {
445 // If we computed the results for the initial (and prob 0 and prob1) states only, we need to filter the
446 // result to only communicate these results.
447 checkResult->filter(ExplicitQualitativeCheckResult(this->getModel().getInitialStates() | psiStates));
448 }
449 return checkResult;
450}
451
452template<typename SparseDtmcModelType>
455 storm::logic::UntilFormula const& pathFormula = checkTask.getFormula();
456
457 // Retrieve the appropriate bitvectors by model checking the subformulas.
458 std::unique_ptr<CheckResult> leftResultPointer = this->check(pathFormula.getLeftSubformula());
459 std::unique_ptr<CheckResult> rightResultPointer = this->check(pathFormula.getRightSubformula());
460 storm::storage::BitVector const& phiStates = leftResultPointer->asExplicitQualitativeCheckResult().getTruthValuesVector();
461 storm::storage::BitVector const& psiStates = rightResultPointer->asExplicitQualitativeCheckResult().getTruthValuesVector();
462
463 return computeUntilProbabilities(this->getModel().getTransitionMatrix(), this->getModel().getBackwardTransitions(), this->getModel().getInitialStates(),
464 phiStates, psiStates, checkTask.isOnlyInitialStatesRelevantSet());
465}
466
467template<typename SparseDtmcModelType>
469 storm::storage::SparseMatrix<ValueType> const& probabilityMatrix, storm::storage::SparseMatrix<ValueType> const& backwardTransitions,
470 storm::storage::BitVector const& initialStates, storm::storage::BitVector const& phiStates, storm::storage::BitVector const& psiStates,
471 bool computeForInitialStatesOnly) {
472 // Then, compute the subset of states that has a probability of 0 or 1, respectively.
473 std::pair<storm::storage::BitVector, storm::storage::BitVector> statesWithProbability01 =
474 storm::utility::graph::performProb01(backwardTransitions, phiStates, psiStates);
475 storm::storage::BitVector statesWithProbability0 = statesWithProbability01.first;
476 storm::storage::BitVector statesWithProbability1 = statesWithProbability01.second;
477 storm::storage::BitVector maybeStates = ~(statesWithProbability0 | statesWithProbability1);
478
479 // Determine whether we need to perform some further computation.
480 bool furtherComputationNeeded = true;
481 if (computeForInitialStatesOnly && initialStates.isDisjointFrom(maybeStates)) {
482 STORM_LOG_DEBUG("The probability for all initial states was found in a preprocessing step.");
483 furtherComputationNeeded = false;
484 } else if (maybeStates.empty()) {
485 STORM_LOG_DEBUG("The probability for all states was found in a preprocessing step.");
486 furtherComputationNeeded = false;
487 }
488
489 std::vector<ValueType> result(maybeStates.size());
490 if (furtherComputationNeeded) {
491 // If we compute the results for the initial states only, we can cut off all maybe state that are not
492 // reachable from them.
493 if (computeForInitialStatesOnly) {
494 // Determine the set of states that is reachable from the initial state without jumping over a target state.
495 storm::storage::BitVector reachableStates =
496 storm::utility::graph::getReachableStates(probabilityMatrix, initialStates, maybeStates, statesWithProbability1);
497
498 // Subtract from the maybe states the set of states that is not reachable (on a path from the initial to a target state).
499 maybeStates &= reachableStates;
500 }
501
502 // Create a vector for the probabilities to go to a state with probability 1 in one step.
503 std::vector<ValueType> oneStepProbabilities = probabilityMatrix.getConstrainedRowSumVector(maybeStates, statesWithProbability1);
504
505 // Determine the set of initial states of the sub-model.
506 storm::storage::BitVector newInitialStates = initialStates % maybeStates;
507
508 // We then build the submatrix that only has the transitions of the maybe states.
509 storm::storage::SparseMatrix<ValueType> submatrix = probabilityMatrix.getSubmatrix(false, maybeStates, maybeStates);
510 storm::storage::SparseMatrix<ValueType> submatrixTransposed = submatrix.transpose();
511
512 std::vector<ValueType> subresult = computeReachabilityValues(submatrix, oneStepProbabilities, submatrixTransposed, newInitialStates,
513 computeForInitialStatesOnly, oneStepProbabilities);
514 storm::utility::vector::setVectorValues<ValueType>(result, maybeStates, subresult);
515 }
516
517 // Construct full result.
518 storm::utility::vector::setVectorValues<ValueType>(result, statesWithProbability0, storm::utility::zero<ValueType>());
519 storm::utility::vector::setVectorValues<ValueType>(result, statesWithProbability1, storm::utility::one<ValueType>());
520
521 if (computeForInitialStatesOnly) {
522 // If we computed the results for the initial (and prob 0 and prob1) states only, we need to filter the
523 // result to only communicate these results.
524 std::unique_ptr<ExplicitQuantitativeCheckResult<ValueType>> checkResult = std::make_unique<ExplicitQuantitativeCheckResult<ValueType>>();
525 for (auto state : ~maybeStates | initialStates) {
526 (*checkResult)[state] = result[state];
527 }
528 return std::move(checkResult); // move() required by, e.g., clang 3.8
529 }
530 return std::make_unique<ExplicitQuantitativeCheckResult<ValueType>>(result);
531}
532
533template<typename SparseDtmcModelType>
536 storm::logic::EventuallyFormula const& eventuallyFormula = checkTask.getFormula();
537
538 // Retrieve the appropriate bitvectors by model checking the subformulas.
539 std::unique_ptr<CheckResult> subResultPointer = this->check(eventuallyFormula.getSubformula());
540 storm::storage::BitVector trueStates(this->getModel().getNumberOfStates(), true);
541 storm::storage::BitVector const& targetStates = subResultPointer->asExplicitQualitativeCheckResult().getTruthValuesVector();
542
543 // Do some sanity checks to establish some required properties.
544 RewardModelType const& rewardModel = this->getModel().getRewardModel(checkTask.isRewardModelSet() ? checkTask.getRewardModel() : "");
545
546 STORM_LOG_THROW(!rewardModel.empty(), storm::exceptions::IllegalArgumentException, "Input model does not have a reward model.");
547 return computeReachabilityRewards(
548 this->getModel().getTransitionMatrix(), this->getModel().getBackwardTransitions(), this->getModel().getInitialStates(), targetStates,
549 [&](uint_fast64_t numberOfRows, storm::storage::SparseMatrix<ValueType> const& transitionMatrix, storm::storage::BitVector const& maybeStates) {
550 return rewardModel.getTotalRewardVector(numberOfRows, transitionMatrix, maybeStates);
551 },
553}
554
555template<typename SparseDtmcModelType>
557 storm::storage::SparseMatrix<ValueType> const& probabilityMatrix, storm::storage::SparseMatrix<ValueType> const& backwardTransitions,
558 storm::storage::BitVector const& initialStates, storm::storage::BitVector const& targetStates, std::vector<ValueType>& stateRewardValues,
559 bool computeForInitialStatesOnly) {
560 return computeReachabilityRewards(
561 probabilityMatrix, backwardTransitions, initialStates, targetStates,
562 [&](uint_fast64_t numberOfRows, storm::storage::SparseMatrix<ValueType> const&, storm::storage::BitVector const& maybeStates) {
563 std::vector<ValueType> result(numberOfRows);
564 storm::utility::vector::selectVectorValues(result, maybeStates, stateRewardValues);
565 return result;
566 },
567 computeForInitialStatesOnly);
568}
569
570template<typename SparseDtmcModelType>
572 storm::storage::SparseMatrix<ValueType> const& probabilityMatrix, storm::storage::SparseMatrix<ValueType> const& backwardTransitions,
573 storm::storage::BitVector const& initialStates, storm::storage::BitVector const& targetStates,
574 std::function<std::vector<ValueType>(uint_fast64_t, storm::storage::SparseMatrix<ValueType> const&, storm::storage::BitVector const&)> const&
575 totalStateRewardVectorGetter,
576 bool computeForInitialStatesOnly) {
577 uint_fast64_t numberOfStates = probabilityMatrix.getRowCount();
578
579 // Compute the subset of states that has a reachability reward less than infinity.
580 storm::storage::BitVector trueStates(numberOfStates, true);
581 storm::storage::BitVector infinityStates = storm::utility::graph::performProb1(backwardTransitions, trueStates, targetStates);
582 infinityStates.complement();
583 storm::storage::BitVector maybeStates = ~targetStates & ~infinityStates;
584
585 // Determine whether we need to perform some further computation.
586 bool furtherComputationNeeded = true;
587 if (computeForInitialStatesOnly) {
588 if (initialStates.isSubsetOf(infinityStates)) {
589 STORM_LOG_DEBUG("The reward of all initial states was found in a preprocessing step.");
590 furtherComputationNeeded = false;
591 }
592 if (initialStates.isSubsetOf(targetStates)) {
593 STORM_LOG_DEBUG("The reward of all initial states was found in a preprocessing step.");
594 furtherComputationNeeded = false;
595 }
596 }
597
598 std::vector<ValueType> result(maybeStates.size());
599 if (furtherComputationNeeded) {
600 // If we compute the results for the initial states only, we can cut off all maybe state that are not
601 // reachable from them.
602 if (computeForInitialStatesOnly) {
603 // Determine the set of states that is reachable from the initial state without jumping over a target state.
604 storm::storage::BitVector reachableStates = storm::utility::graph::getReachableStates(probabilityMatrix, initialStates, maybeStates, targetStates);
605
606 // Subtract from the maybe states the set of states that is not reachable (on a path from the initial to a target state).
607 maybeStates &= reachableStates;
608 }
609
610 // Determine the set of initial states of the sub-model.
611 storm::storage::BitVector newInitialStates = initialStates % maybeStates;
612
613 // We then build the submatrix that only has the transitions of the maybe states.
614 storm::storage::SparseMatrix<ValueType> submatrix = probabilityMatrix.getSubmatrix(false, maybeStates, maybeStates);
615 storm::storage::SparseMatrix<ValueType> submatrixTransposed = submatrix.transpose();
616
617 // Project the state reward vector to all maybe-states.
618 std::vector<ValueType> stateRewardValues = totalStateRewardVectorGetter(submatrix.getRowCount(), probabilityMatrix, maybeStates);
619
620 std::vector<ValueType> subresult =
621 computeReachabilityValues(submatrix, stateRewardValues, submatrixTransposed, newInitialStates, computeForInitialStatesOnly,
622 probabilityMatrix.getConstrainedRowSumVector(maybeStates, targetStates));
623 storm::utility::vector::setVectorValues<ValueType>(result, maybeStates, subresult);
624 }
625
626 // Construct full result.
627 storm::utility::vector::setVectorValues<ValueType>(result, infinityStates, storm::utility::infinity<ValueType>());
628 storm::utility::vector::setVectorValues<ValueType>(result, targetStates, storm::utility::zero<ValueType>());
629 if (computeForInitialStatesOnly) {
630 // If we computed the results for the initial (and inf) states only, we need to filter the result to
631 // only communicate these results.
632 std::unique_ptr<ExplicitQuantitativeCheckResult<ValueType>> checkResult = std::make_unique<ExplicitQuantitativeCheckResult<ValueType>>();
633 for (auto state : ~maybeStates | initialStates) {
634 (*checkResult)[state] = result[state];
635 }
636 return std::move(checkResult); // move() required by, e.g., clang 3.8
637 }
638 return std::make_unique<ExplicitQuantitativeCheckResult<ValueType>>(result);
639}
640
641template<typename SparseDtmcModelType>
644 storm::logic::ConditionalFormula const& conditionalFormula = checkTask.getFormula();
645
646 // Retrieve the appropriate bitvectors by model checking the subformulas.
647 STORM_LOG_THROW(conditionalFormula.getSubformula().isEventuallyFormula(), storm::exceptions::InvalidPropertyException, "Expected 'eventually' formula.");
648 STORM_LOG_THROW(conditionalFormula.getConditionFormula().isEventuallyFormula(), storm::exceptions::InvalidPropertyException,
649 "Expected 'eventually' formula.");
650
651 std::unique_ptr<CheckResult> leftResultPointer = this->check(conditionalFormula.getSubformula().asEventuallyFormula().getSubformula());
652 std::unique_ptr<CheckResult> rightResultPointer = this->check(conditionalFormula.getConditionFormula().asEventuallyFormula().getSubformula());
653 storm::storage::BitVector phiStates = leftResultPointer->asExplicitQualitativeCheckResult().getTruthValuesVector();
654 storm::storage::BitVector psiStates = rightResultPointer->asExplicitQualitativeCheckResult().getTruthValuesVector();
655 storm::storage::BitVector trueStates(this->getModel().getNumberOfStates(), true);
656
657 // Do some sanity checks to establish some required properties.
658 // STORM_LOG_WARN_COND(storm::settings::getModule<storm::settings::modules::EliminationSettings>().getEliminationMethod() ==
659 // storm::settings::modules::EliminationSettings::EliminationMethod::State, "The chosen elimination method is not available for computing conditional
660 // probabilities. Falling back to regular state elimination.");
661 STORM_LOG_THROW(this->getModel().getInitialStates().getNumberOfSetBits() == 1, storm::exceptions::IllegalArgumentException,
662 "Input model is required to have exactly one initial state.");
663 STORM_LOG_THROW(checkTask.isOnlyInitialStatesRelevantSet(), storm::exceptions::IllegalArgumentException,
664 "Cannot compute conditional probabilities for all states.");
665 storm::storage::sparse::state_type initialState = *this->getModel().getInitialStates().begin();
666
667 storm::storage::SparseMatrix<ValueType> backwardTransitions = this->getModel().getBackwardTransitions();
668
669 // Compute the 'true' psi states, i.e. those psi states that can be reached without passing through another psi state first.
670 psiStates = storm::utility::graph::getReachableStates(this->getModel().getTransitionMatrix(), this->getModel().getInitialStates(), trueStates, psiStates) &
671 psiStates;
672
673 std::pair<storm::storage::BitVector, storm::storage::BitVector> statesWithProbability01 =
674 storm::utility::graph::performProb01(backwardTransitions, trueStates, psiStates);
675 storm::storage::BitVector statesWithProbabilityGreater0 = ~statesWithProbability01.first;
676 storm::storage::BitVector statesWithProbability1 = std::move(statesWithProbability01.second);
677
678 STORM_LOG_THROW(this->getModel().getInitialStates().isSubsetOf(statesWithProbabilityGreater0), storm::exceptions::InvalidPropertyException,
679 "The condition of the conditional probability has zero probability.");
680
681 // If the initial state is known to have probability 1 of satisfying the condition, we can apply regular model checking.
682 if (this->getModel().getInitialStates().isSubsetOf(statesWithProbability1)) {
683 STORM_LOG_INFO("The condition holds with probability 1, so the regular reachability probability is computed.");
684 std::shared_ptr<storm::logic::BooleanLiteralFormula> trueFormula = std::make_shared<storm::logic::BooleanLiteralFormula>(true);
685 std::shared_ptr<storm::logic::UntilFormula> untilFormula =
686 std::make_shared<storm::logic::UntilFormula>(trueFormula, conditionalFormula.getSubformula().asSharedPointer());
687 return this->computeUntilProbabilities(env, *untilFormula);
688 }
689
690 // From now on, we know the condition does not have a trivial probability in the initial state.
691
692 // Compute the states that can be reached on a path that has a psi state in it.
693 storm::storage::BitVector statesWithPsiPredecessor =
694 storm::utility::graph::performProbGreater0(this->getModel().getTransitionMatrix(), trueStates, psiStates);
695 storm::storage::BitVector statesReachingPhi = storm::utility::graph::performProbGreater0(backwardTransitions, trueStates, phiStates);
696
697 // The set of states we need to consider are those that have a non-zero probability to satisfy the condition or are on some path that has a psi state in it.
698 storm::storage::BitVector maybeStates = statesWithProbabilityGreater0 | (statesWithPsiPredecessor & statesReachingPhi);
699
700 // Determine the set of initial states of the sub-DTMC.
701 storm::storage::BitVector newInitialStates = this->getModel().getInitialStates() % maybeStates;
702
703 // Create a dummy vector for the one-step probabilities.
704 std::vector<ValueType> oneStepProbabilities(maybeStates.getNumberOfSetBits(), storm::utility::zero<ValueType>());
705
706 // We then build the submatrix that only has the transitions of the maybe states.
707 storm::storage::SparseMatrix<ValueType> submatrix = this->getModel().getTransitionMatrix().getSubmatrix(false, maybeStates, maybeStates);
708 storm::storage::SparseMatrix<ValueType> submatrixTransposed = submatrix.transpose();
709
710 // The states we want to eliminate are those that are tagged with "maybe" but are not a phi or psi state.
711 phiStates = phiStates % maybeStates;
712
713 // If there are no phi states in the reduced model, the conditional probability is trivially zero.
714 if (phiStates.empty()) {
715 return std::unique_ptr<CheckResult>(new ExplicitQuantitativeCheckResult<ValueType>(initialState, storm::utility::zero<ValueType>()));
716 }
717
718 psiStates = psiStates % maybeStates;
719
720 // Keep only the states that we do not eliminate in the maybe states.
721 maybeStates = phiStates | psiStates;
722
723 storm::storage::BitVector statesToEliminate = ~maybeStates & ~newInitialStates;
724
725 // Before starting the model checking process, we assign priorities to states so we can use them to
726 // impose ordering constraints later.
727 boost::optional<std::vector<uint_fast64_t>> distanceBasedPriorities;
731 distanceBasedPriorities = getDistanceBasedPriorities(submatrix, submatrixTransposed, newInitialStates, oneStepProbabilities,
733 }
734
735 storm::storage::FlexibleSparseMatrix<ValueType> flexibleMatrix(submatrix);
736 storm::storage::FlexibleSparseMatrix<ValueType> flexibleBackwardTransitions(submatrixTransposed, true);
737
738 std::shared_ptr<StatePriorityQueue> statePriorities =
739 createStatePriorityQueue(distanceBasedPriorities, flexibleMatrix, flexibleBackwardTransitions, oneStepProbabilities, statesToEliminate);
740
741 STORM_LOG_INFO("Computing conditional probilities.\n");
742 uint_fast64_t numberOfStatesToEliminate = statePriorities->size();
743 STORM_LOG_INFO("Eliminating " << numberOfStatesToEliminate << " states using the state elimination technique.\n");
744 performPrioritizedStateElimination(statePriorities, flexibleMatrix, flexibleBackwardTransitions, oneStepProbabilities, this->getModel().getInitialStates(),
745 true);
746
748 storm::solver::stateelimination::ConditionalStateEliminator<ValueType>(flexibleMatrix, flexibleBackwardTransitions, oneStepProbabilities, phiStates,
749 psiStates);
750
751 // Eliminate the transitions going into the initial state (if there are any).
752 if (!flexibleBackwardTransitions.getRow(*newInitialStates.begin()).empty()) {
753 stateEliminator.eliminateState(*newInitialStates.begin(), false);
754 }
755
756 // Now we need to basically eliminate all chains of not-psi states after phi states and chains of not-phi
757 // states after psi states.
758 for (auto const& trans1 : flexibleMatrix.getRow(*newInitialStates.begin())) {
759 auto initialStateSuccessor = trans1.getColumn();
760
761 STORM_LOG_TRACE("Exploring successor " << initialStateSuccessor << " of the initial state.");
762
763 if (phiStates.get(initialStateSuccessor)) {
764 STORM_LOG_TRACE("Is a phi state.");
765
766 // If the state is both a phi and a psi state, we do not need to eliminate chains.
767 if (psiStates.get(initialStateSuccessor)) {
768 continue;
769 }
770
771 // At this point, we know that the state satisfies phi and not psi.
772 // This means, we must compute the probability to reach psi states, which in turn means that we need
773 // to eliminate all chains of non-psi states between the current state and psi states.
774 bool hasNonPsiSuccessor = true;
775 while (hasNonPsiSuccessor) {
776 stateEliminator.setFilterPhi();
777 hasNonPsiSuccessor = false;
778
779 // Only treat the state if it has an outgoing transition other than a self-loop.
780 auto const currentRow = flexibleMatrix.getRow(initialStateSuccessor);
781 if (currentRow.size() > 1 || (!currentRow.empty() && currentRow.front().getColumn() != initialStateSuccessor)) {
782 for (auto const& element : currentRow) {
783 // If any of the successors is a phi state, we eliminate it (wrt. all its phi predecessors).
784 if (!psiStates.get(element.getColumn())) {
785 FlexibleRowType const& successorRow = flexibleMatrix.getRow(element.getColumn());
786 // Eliminate the successor only if there possibly is a psi state reachable through it.
787 if (successorRow.size() > 1 || (!successorRow.empty() && successorRow.front().getColumn() != element.getColumn())) {
788 STORM_LOG_TRACE("Found non-psi successor " << element.getColumn() << " that needs to be eliminated.");
789 stateEliminator.eliminateState(element.getColumn(), false);
790 hasNonPsiSuccessor = true;
791 }
792 }
793 }
794 STORM_LOG_ASSERT(!flexibleMatrix.getRow(initialStateSuccessor).empty(), "(1) New transitions expected to be non-empty.");
795 }
796 }
797 stateEliminator.unsetFilter();
798 } else {
799 STORM_LOG_ASSERT(psiStates.get(initialStateSuccessor), "Expected psi state.");
800 STORM_LOG_TRACE("Is a psi state.");
801
802 // At this point, we know that the state satisfies psi and not phi.
803 // This means, we must compute the probability to reach phi states, which in turn means that we need
804 // to eliminate all chains of non-phi states between the current state and phi states.
805
806 bool hasNonPhiSuccessor = true;
807 while (hasNonPhiSuccessor) {
808 stateEliminator.setFilterPsi();
809 hasNonPhiSuccessor = false;
810
811 // Only treat the state if it has an outgoing transition other than a self-loop.
812 auto const currentRow = flexibleMatrix.getRow(initialStateSuccessor);
813 if (currentRow.size() > 1 || (!currentRow.empty() && currentRow.front().getColumn() != initialStateSuccessor)) {
814 for (auto const& element : currentRow) {
815 // If any of the successors is a psi state, we eliminate it (wrt. all its psi predecessors).
816 if (!phiStates.get(element.getColumn())) {
817 FlexibleRowType const& successorRow = flexibleMatrix.getRow(element.getColumn());
818 if (successorRow.size() > 1 || (!successorRow.empty() && successorRow.front().getColumn() != element.getColumn())) {
819 STORM_LOG_TRACE("Found non-phi successor " << element.getColumn() << " that needs to be eliminated.");
820 stateEliminator.eliminateState(element.getColumn(), false);
821 hasNonPhiSuccessor = true;
822 }
823 }
824 }
825 }
826 }
827 stateEliminator.unsetFilter();
828 }
829 }
830
831 ValueType numerator = storm::utility::zero<ValueType>();
832 ValueType denominator = storm::utility::zero<ValueType>();
833
834 for (auto const& trans1 : flexibleMatrix.getRow(*newInitialStates.begin())) {
835 auto initialStateSuccessor = trans1.getColumn();
836 if (phiStates.get(initialStateSuccessor)) {
837 if (psiStates.get(initialStateSuccessor)) {
838 numerator += trans1.getValue();
839 denominator += trans1.getValue();
840 } else {
841 ValueType additiveTerm = storm::utility::zero<ValueType>();
842 for (auto const& trans2 : flexibleMatrix.getRow(initialStateSuccessor)) {
843 if (psiStates.get(trans2.getColumn())) {
844 additiveTerm += trans2.getValue();
845 }
846 }
847 additiveTerm *= trans1.getValue();
848 numerator += additiveTerm;
849 denominator += additiveTerm;
850 }
851 } else {
852 STORM_LOG_ASSERT(psiStates.get(initialStateSuccessor), "Expected psi state.");
853 denominator += trans1.getValue();
854 ValueType additiveTerm = storm::utility::zero<ValueType>();
855 for (auto const& trans2 : flexibleMatrix.getRow(initialStateSuccessor)) {
856 if (phiStates.get(trans2.getColumn())) {
857 additiveTerm += trans2.getValue();
858 }
859 }
860 numerator += trans1.getValue() * additiveTerm;
861 }
862 }
863
864 return std::unique_ptr<CheckResult>(new ExplicitQuantitativeCheckResult<ValueType>(initialState, numerator / denominator));
865}
866
867template<typename SparseDtmcModelType>
869 std::shared_ptr<StatePriorityQueue>& priorityQueue, storm::storage::FlexibleSparseMatrix<ValueType>& transitionMatrix,
870 storm::storage::FlexibleSparseMatrix<ValueType>& backwardTransitions, std::vector<ValueType>& values, storm::storage::BitVector const& initialStates,
871 bool computeResultsForInitialStatesOnly) {
872 storm::solver::stateelimination::PrioritizedStateEliminator<ValueType> stateEliminator(transitionMatrix, backwardTransitions, priorityQueue, values);
873
874 while (priorityQueue->hasNext()) {
875 storm::storage::sparse::state_type state = priorityQueue->pop();
876 bool removeForwardTransitions = computeResultsForInitialStatesOnly && !initialStates.get(state);
877 stateEliminator.eliminateState(state, removeForwardTransitions);
878 if (removeForwardTransitions) {
879 values[state] = storm::utility::zero<ValueType>();
880 }
881#ifdef STORM_DEV
882 STORM_LOG_ASSERT(checkConsistent(transitionMatrix, backwardTransitions), "The forward and backward transition matrices became inconsistent.");
883#endif
884 }
885}
886
887template<typename SparseDtmcModelType>
888void SparseDtmcEliminationModelChecker<SparseDtmcModelType>::performOrdinaryStateElimination(
890 storm::storage::BitVector const& subsystem, storm::storage::BitVector const& initialStates, bool computeResultsForInitialStatesOnly,
891 std::vector<ValueType>& values, boost::optional<std::vector<uint_fast64_t>> const& distanceBasedPriorities) {
892 std::shared_ptr<StatePriorityQueue> statePriorities =
893 createStatePriorityQueue(distanceBasedPriorities, transitionMatrix, backwardTransitions, values, subsystem);
894
895 std::size_t numberOfStatesToEliminate = statePriorities->size();
896 STORM_LOG_DEBUG("Eliminating " << numberOfStatesToEliminate << " states using the state elimination technique.\n");
897 performPrioritizedStateElimination(statePriorities, transitionMatrix, backwardTransitions, values, initialStates, computeResultsForInitialStatesOnly);
898 STORM_LOG_DEBUG("Eliminated " << numberOfStatesToEliminate << " states.\n");
899}
900
901template<typename SparseDtmcModelType>
902uint_fast64_t SparseDtmcEliminationModelChecker<SparseDtmcModelType>::performHybridStateElimination(
905 storm::storage::BitVector const& initialStates, bool computeResultsForInitialStatesOnly, std::vector<ValueType>& values,
906 boost::optional<std::vector<uint_fast64_t>> const& distanceBasedPriorities) {
907 // When using the hybrid technique, we recursively treat the SCCs up to some size.
908 std::vector<storm::storage::sparse::state_type> entryStateQueue;
909 STORM_LOG_DEBUG("Eliminating " << subsystem.size() << " states using the hybrid elimination technique.\n");
910 uint_fast64_t maximalDepth = treatScc(transitionMatrix, values, initialStates, subsystem, initialStates, forwardTransitions, backwardTransitions, false, 0,
912 computeResultsForInitialStatesOnly, distanceBasedPriorities);
913
914 // If the entry states were to be eliminated last, we need to do so now.
916 STORM_LOG_DEBUG("Eliminating " << entryStateQueue.size() << " entry states as a last step.");
917 std::vector<storm::storage::sparse::state_type> sortedStates(entryStateQueue.begin(), entryStateQueue.end());
918 std::shared_ptr<StatePriorityQueue> queuePriorities = std::make_shared<StaticStatePriorityQueue>(sortedStates);
919 performPrioritizedStateElimination(queuePriorities, transitionMatrix, backwardTransitions, values, initialStates, computeResultsForInitialStatesOnly);
920 }
921 STORM_LOG_DEBUG("Eliminated " << subsystem.size() << " states.\n");
922 return maximalDepth;
923}
924
925template<typename SparseDtmcModelType>
926std::vector<typename SparseDtmcEliminationModelChecker<SparseDtmcModelType>::ValueType>
927SparseDtmcEliminationModelChecker<SparseDtmcModelType>::computeReachabilityValues(storm::storage::SparseMatrix<ValueType> const& transitionMatrix,
928 std::vector<ValueType>& values,
929 storm::storage::SparseMatrix<ValueType> const& backwardTransitions,
930 storm::storage::BitVector const& initialStates,
931 bool computeResultsForInitialStatesOnly,
932 std::vector<ValueType> const& oneStepProbabilitiesToTarget) {
933 // Then, we convert the reduced matrix to a more flexible format to be able to perform state elimination more easily.
934 storm::storage::FlexibleSparseMatrix<ValueType> flexibleMatrix(transitionMatrix);
935 storm::storage::FlexibleSparseMatrix<ValueType> flexibleBackwardTransitions(backwardTransitions);
936
939 boost::optional<std::vector<uint_fast64_t>> distanceBasedPriorities;
941 distanceBasedPriorities = getDistanceBasedPriorities(transitionMatrix, backwardTransitions, initialStates, oneStepProbabilitiesToTarget,
943 }
944
945 // Create a bit vector that represents the subsystem of states we still have to eliminate.
946 storm::storage::BitVector subsystem = storm::storage::BitVector(transitionMatrix.getRowCount(), true);
947
950 performOrdinaryStateElimination(flexibleMatrix, flexibleBackwardTransitions, subsystem, initialStates, computeResultsForInitialStatesOnly, values,
951 distanceBasedPriorities);
954 uint64_t maximalDepth = performHybridStateElimination(transitionMatrix, flexibleMatrix, flexibleBackwardTransitions, subsystem, initialStates,
955 computeResultsForInitialStatesOnly, values, distanceBasedPriorities);
956 STORM_LOG_TRACE("Maximal depth of decomposition was " << maximalDepth << ".");
957 }
958
959 STORM_LOG_ASSERT(flexibleMatrix.empty(), "Not all transitions were eliminated.");
960 STORM_LOG_ASSERT(flexibleBackwardTransitions.empty(), "Not all transitions were eliminated.");
961
962 // Now, we return the value for the only initial state.
963 STORM_LOG_DEBUG("Simplifying and returning result.");
964 for (auto& value : values) {
965 value = storm::utility::simplify(value);
966 }
967 return values;
968}
969
970template<typename SparseDtmcModelType>
971uint_fast64_t SparseDtmcEliminationModelChecker<SparseDtmcModelType>::treatScc(
972 storm::storage::FlexibleSparseMatrix<ValueType>& matrix, std::vector<ValueType>& values, storm::storage::BitVector const& entryStates,
973 storm::storage::BitVector const& scc, storm::storage::BitVector const& initialStates, storm::storage::SparseMatrix<ValueType> const& forwardTransitions,
974 storm::storage::FlexibleSparseMatrix<ValueType>& backwardTransitions, bool eliminateEntryStates, uint_fast64_t level, uint_fast64_t maximalSccSize,
975 std::vector<storm::storage::sparse::state_type>& entryStateQueue, bool computeResultsForInitialStatesOnly,
976 boost::optional<std::vector<uint_fast64_t>> const& distanceBasedPriorities) {
977 uint_fast64_t maximalDepth = level;
978
979 // If the SCCs are large enough, we try to split them further.
980 if (scc.getNumberOfSetBits() > maximalSccSize) {
981 STORM_LOG_TRACE("SCC is large enough (" << scc.getNumberOfSetBits() << " states) to be decomposed further.");
982
983 // Here, we further decompose the SCC into sub-SCCs.
984 storm::storage::BitVector nonEntrySccStates = scc & ~entryStates;
986 forwardTransitions, storm::storage::StronglyConnectedComponentDecompositionOptions().subsystem(nonEntrySccStates));
987 STORM_LOG_TRACE("Decomposed SCC into " << decomposition.size() << " sub-SCCs.");
988
989 // Store a bit vector of remaining SCCs so we can be flexible when it comes to the order in which
990 // we eliminate the SCCs.
991 storm::storage::BitVector remainingSccs(decomposition.size(), true);
992
993 // First, get rid of the trivial SCCs.
994 storm::storage::BitVector statesInTrivialSccs(matrix.getRowCount());
995 for (uint_fast64_t sccIndex = 0; sccIndex < decomposition.size(); ++sccIndex) {
996 storm::storage::StronglyConnectedComponent const& scc = decomposition.getBlock(sccIndex);
997 if (scc.isTrivial()) {
998 // Put the only state of the trivial SCC into the set of states to eliminate.
999 statesInTrivialSccs.set(*scc.begin(), true);
1000 remainingSccs.set(sccIndex, false);
1001 }
1002 }
1003
1004 std::shared_ptr<StatePriorityQueue> statePriorities =
1005 createStatePriorityQueue(distanceBasedPriorities, matrix, backwardTransitions, values, statesInTrivialSccs);
1006 STORM_LOG_TRACE("Eliminating " << statePriorities->size() << " trivial SCCs.");
1007 performPrioritizedStateElimination(statePriorities, matrix, backwardTransitions, values, initialStates, computeResultsForInitialStatesOnly);
1008 STORM_LOG_TRACE("Eliminated all trivial SCCs.");
1009
1010 // And then recursively treat the remaining sub-SCCs.
1011 STORM_LOG_TRACE("Eliminating " << remainingSccs.getNumberOfSetBits() << " remaining SCCs on level " << level << ".");
1012 for (auto sccIndex : remainingSccs) {
1013 storm::storage::StronglyConnectedComponent const& newScc = decomposition.getBlock(sccIndex);
1014
1015 // Rewrite SCC into bit vector and subtract it from the remaining states.
1016 storm::storage::BitVector newSccAsBitVector(forwardTransitions.getRowCount(), newScc.begin(), newScc.end());
1017
1018 // Determine the set of entry states of the SCC.
1019 storm::storage::BitVector entryStates(forwardTransitions.getRowCount());
1020 for (auto const& state : newScc) {
1021 for (auto const& predecessor : backwardTransitions.getRow(state)) {
1022 if (predecessor.getValue() != storm::utility::zero<ValueType>() && !newSccAsBitVector.get(predecessor.getColumn())) {
1023 entryStates.set(state);
1024 }
1025 }
1026 }
1027
1028 // Recursively descend in SCC-hierarchy.
1029 uint_fast64_t depth =
1030 treatScc(matrix, values, entryStates, newSccAsBitVector, initialStates, forwardTransitions, backwardTransitions,
1031 eliminateEntryStates || !storm::settings::getModule<storm::settings::modules::EliminationSettings>().isEliminateEntryStatesLastSet(),
1032 level + 1, maximalSccSize, entryStateQueue, computeResultsForInitialStatesOnly, distanceBasedPriorities);
1033 maximalDepth = std::max(maximalDepth, depth);
1034 }
1035 } else {
1036 // In this case, we perform simple state elimination in the current SCC.
1037 STORM_LOG_TRACE("SCC of size " << scc.getNumberOfSetBits() << " is small enough to be eliminated directly.");
1038 std::shared_ptr<StatePriorityQueue> statePriorities =
1039 createStatePriorityQueue(distanceBasedPriorities, matrix, backwardTransitions, values, scc & ~entryStates);
1040 performPrioritizedStateElimination(statePriorities, matrix, backwardTransitions, values, initialStates, computeResultsForInitialStatesOnly);
1041 STORM_LOG_TRACE("Eliminated all states of SCC.");
1042 }
1043
1044 // Finally, eliminate the entry states (if we are required to do so).
1045 if (eliminateEntryStates) {
1046 STORM_LOG_TRACE("Finally, eliminating entry states.");
1047 std::shared_ptr<StatePriorityQueue> naivePriorities = createStatePriorityQueue(entryStates);
1048 performPrioritizedStateElimination(naivePriorities, matrix, backwardTransitions, values, initialStates, computeResultsForInitialStatesOnly);
1049 STORM_LOG_TRACE("Eliminated/added entry states.");
1050 } else {
1051 STORM_LOG_TRACE("Finally, adding entry states to queue.");
1052 for (auto state : entryStates) {
1053 entryStateQueue.push_back(state);
1054 }
1055 }
1056
1057 return maximalDepth;
1058}
1059
1060template<typename SparseDtmcModelType>
1061bool SparseDtmcEliminationModelChecker<SparseDtmcModelType>::checkConsistent(storm::storage::FlexibleSparseMatrix<ValueType>& transitionMatrix,
1063 for (uint_fast64_t forwardIndex = 0; forwardIndex < transitionMatrix.getRowCount(); ++forwardIndex) {
1064 for (auto const& forwardEntry : transitionMatrix.getRow(forwardIndex)) {
1065 if (forwardEntry.getColumn() == forwardIndex) {
1066 continue;
1067 }
1068
1069 bool foundCorrespondingElement = false;
1070 for (auto const& backwardEntry : backwardTransitions.getRow(forwardEntry.getColumn())) {
1071 if (backwardEntry.getColumn() == forwardIndex) {
1072 foundCorrespondingElement = true;
1073 }
1074 }
1075
1076 if (!foundCorrespondingElement) {
1077 return false;
1078 }
1079 }
1080 }
1081 return true;
1082}
1083
1084template class SparseDtmcEliminationModelChecker<storm::models::sparse::Dtmc<double>>;
1085
1086#ifdef STORM_HAVE_CARL
1087template class SparseDtmcEliminationModelChecker<storm::models::sparse::Dtmc<storm::RationalNumber>>;
1088template class SparseDtmcEliminationModelChecker<storm::models::sparse::Dtmc<storm::RationalFunction>>;
1089#endif
1090} // namespace modelchecker
1091} // namespace storm
Formula const & getRightSubformula() const
Formula const & getLeftSubformula() const
storm::expressions::Expression const & getUpperBound(unsigned i=0) const
bool hasIntegerUpperBound(unsigned i=0) const
Formula const & getConditionFormula() const
Formula const & getSubformula() const
EventuallyFormula & asEventuallyFormula()
Definition Formula.cpp:333
bool isInFragment(FragmentSpecification const &fragment) const
Definition Formula.cpp:196
virtual bool isEventuallyFormula() const
Definition Formula.cpp:88
std::shared_ptr< Formula const > asSharedPointer()
Definition Formula.cpp:548
FragmentSpecification & setOnlyEventuallyFormuluasInConditionalFormulasAllowed(bool newValue)
FragmentSpecification & setCumulativeRewardFormulasAllowed(bool newValue)
FragmentSpecification & setConditionalProbabilityFormulasAllowed(bool newValue)
FragmentSpecification & setLongRunAverageProbabilitiesAllowed(bool newValue)
FragmentSpecification & setInstantaneousFormulasAllowed(bool newValue)
FragmentSpecification & setNestedOperatorsAllowed(bool newValue)
Formula const & getSubformula() const
bool isRewardModelSet() const
Retrieves whether a reward model was set.
Definition CheckTask.h:190
std::string const & getRewardModel() const
Retrieves the reward model over which to perform the checking (if set).
Definition CheckTask.h:197
FormulaType const & getFormula() const
Retrieves the formula from this task.
Definition CheckTask.h:140
bool isOnlyInitialStatesRelevantSet() const
Retrieves whether only the initial states are relevant in the computation.
Definition CheckTask.h:204
virtual std::unique_ptr< CheckResult > computeBoundedUntilProbabilities(Environment const &env, CheckTask< storm::logic::BoundedUntilFormula, ValueType > const &checkTask) override
virtual std::unique_ptr< CheckResult > computeUntilProbabilities(Environment const &env, CheckTask< storm::logic::UntilFormula, ValueType > const &checkTask) override
virtual std::unique_ptr< CheckResult > computeLongRunAverageRewards(Environment const &env, CheckTask< storm::logic::LongRunAverageRewardFormula, ValueType > const &checkTask) override
SparseDtmcEliminationModelChecker(storm::models::sparse::Dtmc< ValueType > const &model)
Creates an elimination-based model checker for the given model.
virtual std::unique_ptr< CheckResult > computeLongRunAverageProbabilities(Environment const &env, CheckTask< storm::logic::StateFormula, ValueType > const &checkTask) override
virtual std::unique_ptr< CheckResult > computeReachabilityRewards(Environment const &env, CheckTask< storm::logic::EventuallyFormula, ValueType > const &checkTask) override
storm::storage::FlexibleSparseMatrix< ValueType >::row_type FlexibleRowType
virtual std::unique_ptr< CheckResult > computeConditionalProbabilities(Environment const &env, CheckTask< storm::logic::ConditionalFormula, ValueType > const &checkTask) override
virtual bool canHandle(CheckTask< storm::logic::Formula, ValueType > const &checkTask) const override
This class represents a discrete-time Markov chain.
Definition Dtmc.h:14
EliminationOrder
An enum that contains all available state elimination orders.
void eliminateState(storm::storage::sparse::state_type state, bool removeForwardTransitions)
A bit vector that is internally represented as a vector of 64-bit values.
Definition BitVector.h:18
void complement()
Negates all bits in the bit vector.
bool isDisjointFrom(BitVector const &other) const
Checks whether none of the bits that are set in the current bit vector are also set in the given bit ...
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.
bool isSubsetOf(BitVector const &other) const
Checks whether all bits that are set in the current bit vector are also set in the given bit vector.
const_iterator begin() const
Returns an iterator to the indices of the set bits in the bit vector.
void set(uint_fast64_t index, bool value=true)
Sets the given truth value at the given index.
size_t size() const
Retrieves the number of bits this bit vector can store.
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.
The flexible sparse matrix is used during state elimination.
row_type & getRow(index_type)
Returns an object representing the given row.
index_type getRowCount() const
Returns the number of rows of the matrix.
A class that holds a possibly non-square matrix in the compressed row storage format.
const_rows getRow(index_type row) const
Returns an object representing the given row.
void multiplyWithVector(std::vector< value_type > const &vector, std::vector< value_type > &result, std::vector< value_type > const *summand=nullptr) const
Multiplies the matrix with the given vector and writes the result to the given result vector.
SparseMatrix getSubmatrix(bool useGroups, storm::storage::BitVector const &rowConstraint, storm::storage::BitVector const &columnConstraint, bool insertDiagonalEntries=false, storm::storage::BitVector const &makeZeroColumns=storm::storage::BitVector()) const
Creates a submatrix of the current matrix by dropping all rows and columns whose bits are not set to ...
std::vector< value_type > getConstrainedRowSumVector(storm::storage::BitVector const &rowConstraint, storm::storage::BitVector const &columnConstraint) const
Computes a vector whose i-th entry is the sum of the entries in the i-th selected row where only thos...
storm::storage::SparseMatrix< value_type > transpose(bool joinGroups=false, bool keepZeros=false) const
Transposes the matrix.
index_type getRowCount() const
Returns the number of rows of the matrix.
iterator begin()
Returns an iterator to the states in this SCC.
Definition StateBlock.cpp:5
iterator end()
Returns an iterator that points one past the end of the states in this SCC.
This class represents the decomposition of a graph-like structure into its strongly connected compone...
This class represents a strongly connected component, i.e., a set of states such that every state can...
bool isTrivial() const
Retrieves whether this SCC is trivial.
#define STORM_LOG_INFO(message)
Definition logging.h:29
#define STORM_LOG_DEBUG(message)
Definition logging.h:23
#define STORM_LOG_TRACE(message)
Definition logging.h:17
#define STORM_LOG_ASSERT(cond, message)
Definition macros.h:11
#define STORM_LOG_THROW(cond, exception, message)
Definition macros.h:30
#define STORM_PRINT_AND_LOG(message)
Definition macros.h:68
SFTBDDChecker::ValueType ValueType
Expression ite(Expression const &condition, Expression const &thenExpression, Expression const &elseExpression)
FragmentSpecification prctl()
SettingsType const & getModule()
Get module.
std::pair< storm::storage::BitVector, storm::storage::BitVector > performProb01(storm::models::sparse::DeterministicModel< T > const &model, storm::storage::BitVector const &phiStates, storm::storage::BitVector const &psiStates)
Computes the sets of states that have probability 0 or 1, respectively, of satisfying phi until psi i...
Definition graph.cpp:400
storm::storage::BitVector getReachableStates(storm::storage::SparseMatrix< T > const &transitionMatrix, storm::storage::BitVector const &initialStates, storm::storage::BitVector const &constraintStates, storm::storage::BitVector const &targetStates, bool useStepBound, uint_fast64_t maximalSteps, boost::optional< storm::storage::BitVector > const &choiceFilter)
Performs a forward depth-first search through the underlying graph structure to identify the states t...
Definition graph.cpp:48
storm::storage::BitVector performProbGreater0(storm::storage::SparseMatrix< T > const &backwardTransitions, storm::storage::BitVector const &phiStates, storm::storage::BitVector const &psiStates, bool useStepBound, uint_fast64_t maximalSteps)
Performs a backward depth-first search trough the underlying graph structure of the given model to de...
Definition graph.cpp:322
storm::storage::BitVector performProb1(storm::storage::SparseMatrix< T > const &backwardTransitions, storm::storage::BitVector const &, storm::storage::BitVector const &psiStates, storm::storage::BitVector const &statesWithProbabilityGreater0)
Computes the set of states of the given model for which all paths lead to the given set of target sta...
Definition graph.cpp:383
std::vector< uint_fast64_t > getDistances(storm::storage::SparseMatrix< T > const &transitionMatrix, storm::storage::BitVector const &initialStates, boost::optional< storm::storage::BitVector > const &subsystem)
Performs a breadth-first search through the underlying graph structure to compute the distance from a...
Definition graph.cpp:288
std::shared_ptr< StatePriorityQueue > createStatePriorityQueue(boost::optional< std::vector< uint_fast64_t > > const &distanceBasedStatePriorities, storm::storage::FlexibleSparseMatrix< ValueType > const &transitionMatrix, storm::storage::FlexibleSparseMatrix< ValueType > const &backwardTransitions, std::vector< ValueType > const &oneStepProbabilities, storm::storage::BitVector const &states)
bool eliminationOrderNeedsForwardDistances(storm::settings::modules::EliminationSettings::EliminationOrder const &order)
bool eliminationOrderNeedsDistances(storm::settings::modules::EliminationSettings::EliminationOrder const &order)
std::vector< uint_fast64_t > getDistanceBasedPriorities(storm::storage::SparseMatrix< ValueType > const &transitionMatrix, storm::storage::SparseMatrix< ValueType > const &transitionMatrixTransposed, storm::storage::BitVector const &initialStates, std::vector< ValueType > const &oneStepProbabilities, bool forward, bool reverse)
bool eliminationOrderNeedsReversedDistances(storm::settings::modules::EliminationSettings::EliminationOrder const &order)
void addVectors(std::vector< InValueType1 > const &firstOperand, std::vector< InValueType2 > const &secondOperand, std::vector< OutValueType > &target)
Adds the two given vectors and writes the result to the target vector.
Definition vector.h:443
void setVectorValues(std::vector< T > &vector, storm::storage::BitVector const &positions, std::vector< T > const &values)
Sets the provided values at the provided positions in the given vector.
Definition vector.h:83
void selectVectorValues(std::vector< T > &vector, storm::storage::BitVector const &positions, std::vector< T > const &values)
Selects the elements from a vector at the specified positions and writes them consecutively into anot...
Definition vector.h:189
ValueType simplify(ValueType value)
LabParser.cpp.
Definition cli.cpp:18