Storm 1.11.1.1
A Modern Probabilistic Model Checker
Loading...
Searching...
No Matches
SparseDeterministicInfiniteHorizonHelper.cpp
Go to the documentation of this file.
2
3#include <numeric>
4
10
14
16
20
24
27
28namespace storm {
29namespace modelchecker {
30namespace helper {
31
32template<typename ValueType>
37
38template<typename ValueType>
40 std::vector<ValueType> const& exitRates)
41 : SparseInfiniteHorizonHelper<ValueType, false>(transitionMatrix, exitRates) {}
42
43template<typename ValueType>
45 if (this->_longRunComponentDecomposition == nullptr) {
46 // The decomposition has not been provided or computed, yet.
47 this->_computedLongRunComponentDecomposition = std::make_unique<storm::storage::StronglyConnectedComponentDecomposition<ValueType>>(
49 this->_longRunComponentDecomposition = this->_computedLongRunComponentDecomposition.get();
50 }
51}
52
53template<typename ValueType>
55 ValueGetter const& actionValueGetter,
57 // For deterministic models, we compute the LRA for a BSCC
58
59 STORM_LOG_ASSERT(!this->isProduceSchedulerSet(), "Scheduler production enabled for deterministic model.");
60
61 auto trivialResult = computeLraForTrivialBscc(env, stateValueGetter, actionValueGetter, component);
62 if (trivialResult.first) {
63 return trivialResult.second;
64 }
65
66 // Solve nontrivial BSCC with the method specified in the settings
67 storm::solver::LraMethod method = env.solver().lra().getDetLraMethod();
69 method == storm::solver::LraMethod::ValueIteration) {
70 method = storm::solver::LraMethod::GainBiasEquations;
71 STORM_LOG_INFO("Selecting " << storm::solver::toString(method)
72 << " as the solution technique for long-run properties to guarantee exact results. If you want to override this, please "
73 "explicitly specify a different LRA method.");
74 } else if (env.solver().isForceSoundness() && env.solver().lra().isDetLraMethodSetFromDefault() && method != storm::solver::LraMethod::ValueIteration) {
75 method = storm::solver::LraMethod::ValueIteration;
76 STORM_LOG_INFO("Selecting " << storm::solver::toString(method)
77 << " as the solution technique for long-run properties to guarantee sound results. If you want to override this, please "
78 "explicitly specify a different LRA method.");
79 }
80 STORM_LOG_TRACE("Computing LRA for BSCC of size " << component.size() << " using '" << storm::solver::toString(method) << "'.");
81 if (method == storm::solver::LraMethod::ValueIteration) {
82 return computeLraForBsccVi(env, stateValueGetter, actionValueGetter, component);
83 } else if (method == storm::solver::LraMethod::LraDistributionEquations) {
84 // We only need the first element of the pair as the lra distribution is not relevant at this point.
85 return computeLraForBsccSteadyStateDistr(env, stateValueGetter, actionValueGetter, component).first;
86 }
87 STORM_LOG_WARN_COND(method == storm::solver::LraMethod::GainBiasEquations,
88 "Unsupported lra method selected. Defaulting to " << storm::solver::toString(storm::solver::LraMethod::GainBiasEquations) << ".");
89 // We don't need the bias values
90 return computeLraForBsccGainBias(env, stateValueGetter, actionValueGetter, component).first;
91}
92
93template<typename ValueType>
95 Environment const& env, ValueGetter const& stateValueGetter, ValueGetter const& actionValueGetter,
97 // For deterministic models, we can catch the case where all values are the same. This includes the special case where the BSCC consist only of just one
98 // state.
99 bool first = true;
100 ValueType val = storm::utility::zero<ValueType>();
101 for (auto const& element : component) {
102 auto state = internal::getComponentElementState(element);
104 "Unexpected choice index at state " << state << " of deterministic model.");
105 ValueType curr =
106 stateValueGetter(state) + (this->isContinuousTime() ? (*this->_exitRates)[state] * actionValueGetter(state) : actionValueGetter(state));
107 if (first) {
108 val = curr;
109 first = false;
110 } else if (val != curr) {
111 return {false, storm::utility::zero<ValueType>()};
112 }
113 }
114 // All values are the same
115 return {true, val};
116}
117
118template<>
120 Environment const& /*env*/, ValueGetter const& /*stateValueGetter*/, ValueGetter const& /*actionValueGetter*/,
122 STORM_LOG_THROW(false, storm::exceptions::NotSupportedException, "The requested Method for LRA computation is not supported for parametric models.");
123}
124template<typename ValueType>
126 ValueGetter const& actionValueGetter,
128 // Collect parameters of the computation
129 ValueType aperiodicFactor = storm::utility::convertNumber<ValueType>(env.solver().lra().getAperiodicFactor());
130
131 // Now create a helper and perform the algorithm
132 if (this->isContinuousTime()) {
133 // We assume a CTMC (with deterministic timed states and no instant states)
136 viHelper(bscc, this->_transitionMatrix, aperiodicFactor, this->_markovianStates, this->_exitRates);
137 return viHelper.performValueIteration(env, stateValueGetter, actionValueGetter, this->_exitRates);
138 } else {
139 // We assume a DTMC (with deterministic timed states and no instant states)
142 viHelper(bscc, this->_transitionMatrix, aperiodicFactor);
143 return viHelper.performValueIteration(env, stateValueGetter, actionValueGetter);
144 }
145}
146
147template<typename ValueType>
149 Environment const& env, ValueGetter const& stateValuesGetter, ValueGetter const& actionValuesGetter,
151 // we want that the returned vector is sorted as the bscc. So let's assert that the bscc is sorted ascendingly.
152 STORM_LOG_ASSERT(std::is_sorted(bscc.begin(), bscc.end()), "Expected that bsccs are sorted.");
153
154 // We build the equation system as in Line 3 of Algorithm 3 from
155 // Kretinsky, Meggendorfer: Efficient Strategy Iteration for Mean Payoff in Markov Decision Processes (ATVA 2017)
156 // https://doi.org/10.1007/978-3-319-68167-2_25
157 // The first variable corresponds to the gain of the bscc whereas the subsequent variables yield the bias for each state s_1, s_2, ....
158 // No bias variable for s_0 is needed since it is always set to zero, yielding an nxn equation system matrix
159 // To make this work for CTMC, we could uniformize the model. This preserves LRA and ensures that we can compute the
160 // LRA as for a DTMC (the soujourn time in each state is the same). If we then multiply the equations with the uniformization rate,
161 // the uniformization rate cancels out. Hence, we obtain the equation system below.
162
163 // Get a mapping from global state indices to local ones.
164 std::unordered_map<uint64_t, uint64_t> toLocalIndexMap;
165 uint64_t localIndex = 0;
166 for (auto const& globalIndex : bscc) {
167 toLocalIndexMap[globalIndex] = localIndex;
168 ++localIndex;
169 }
170
171 // Prepare an environment for the underlying equation solver
172 auto subEnv = env;
173 if (subEnv.solver().getLinearEquationSolverType() == storm::solver::EquationSolverType::Topological) {
174 // Topological solver does not make any sense since the BSCC is connected.
175 subEnv.solver().setLinearEquationSolverType(subEnv.solver().topological().getUnderlyingEquationSolverType(),
176 subEnv.solver().topological().isUnderlyingEquationSolverTypeSetFromDefault());
177 }
178 subEnv.solver().setLinearEquationSolverPrecision(env.solver().lra().getPrecision(), env.solver().lra().getRelativeTerminationCriterion());
179
180 // Build the equation system matrix and vector.
182 bool isEquationSystemFormat =
185 std::vector<ValueType> eqSysVector;
186 eqSysVector.reserve(bscc.size());
187 // The first row asserts that the weighted bias variables and the reward at s_0 sum up to the gain
188 uint64_t row = 0;
189 ValueType entryValue;
190 for (auto const& globalState : bscc) {
191 ValueType rateAtState = this->_exitRates ? (*this->_exitRates)[globalState] : storm::utility::one<ValueType>();
192 // Coefficient for the gain variable
193 if (isEquationSystemFormat) {
194 // '1-0' in row 0 and -(-1) in other rows
195 builder.addNextValue(row, 0, storm::utility::one<ValueType>());
196 } else if (row > 0) {
197 // No coeficient in row 0, othwerise substract the gain
198 builder.addNextValue(row, 0, -storm::utility::one<ValueType>());
199 }
200 // Compute weighted sum over successor state. As this is a BSCC, each successor state will again be in the BSCC.
201 if (row > 0) {
202 if (isEquationSystemFormat) {
203 builder.addDiagonalEntry(row, rateAtState);
204 } else if (!storm::utility::isOne(rateAtState)) {
205 builder.addDiagonalEntry(row, storm::utility::one<ValueType>() - rateAtState);
206 }
207 }
208 for (auto const& entry : this->_transitionMatrix.getRow(globalState)) {
209 uint64_t col = toLocalIndexMap[entry.getColumn()];
210 if (col == 0) {
211 // Skip transition to state_0. This corresponds to setting the bias of state_0 to zero
212 continue;
213 }
214 entryValue = entry.getValue() * rateAtState;
215 if (isEquationSystemFormat) {
216 entryValue = -entryValue;
217 }
218 builder.addNextValue(row, col, entryValue);
219 }
220 eqSysVector.push_back(stateValuesGetter(globalState) + rateAtState * actionValuesGetter(globalState));
221 ++row;
222 }
223
224 // Create a linear equation solver
225 auto solver = linearEquationSolverFactory.create(subEnv, builder.build());
226 // Check solver requirements.
227 auto requirements = solver->getRequirements(subEnv);
228 STORM_LOG_THROW(!requirements.hasEnabledCriticalRequirement(), storm::exceptions::UnmetRequirementException,
229 "Solver requirements " + requirements.getEnabledRequirementsAsString() + " not checked.");
230 // Todo: Find bounds on the bias variables. Just inserting the maximal value from the vector probably does not work.
231
232 std::vector<ValueType> eqSysSol(bscc.size(), storm::utility::zero<ValueType>());
233 // Take the mean of the rewards as an initial guess for the gain
234 // eqSysSol.front() = std::accumulate(eqSysVector.begin(), eqSysVector.end(), storm::utility::zero<ValueType>()) / storm::utility::convertNumber<ValueType,
235 // uint64_t>(bscc.size());
236 solver->solveEquations(subEnv, eqSysSol, eqSysVector);
237
238 ValueType gain = eqSysSol.front();
239 // insert bias value for state 0
240 eqSysSol.front() = storm::utility::zero<ValueType>();
241 // Return the gain and the bias values
242 return std::pair<ValueType, std::vector<ValueType>>(std::move(gain), std::move(eqSysSol));
243}
244
245template<typename ValueType>
248 // We catch the (easy) case where the BSCC is a singleton.
249 if (bscc.size() == 1) {
250 return {storm::utility::one<ValueType>()};
251 }
254 if (env.solver().isForceSoundness()) {
256 } else {
258 }
259 }
262 }
263
265 return computeSteadyStateDistrForBsccEqSys(env, bscc);
266 } else {
268 "Unexpected algorithm for steady state distribution computation.");
269 return computeSteadyStateDistrForBsccEVTs(env, bscc);
270 }
271}
272
273template<typename ValueType>
276 // Computes steady state distributions by computing EVTs on a slightly modified system. See https://arxiv.org/abs/2401.10638 for more information.
277 storm::storage::BitVector bsccAsBitVector(this->_transitionMatrix.getColumnCount(), false);
278 bsccAsBitVector.set(bscc.begin(), bscc.end(), true);
279
280 // Remove one state from the BSCC (so it becomes an SCC) and compute visiting times on that SCC.
281 uint64_t const proxyState = *bsccAsBitVector.rbegin();
282 bsccAsBitVector.set(proxyState, false);
283 auto visittimesHelper = this->isContinuousTime() ? SparseDeterministicVisitingTimesHelper<ValueType>(this->_transitionMatrix, *this->_exitRates)
284 : SparseDeterministicVisitingTimesHelper<ValueType>(this->_transitionMatrix);
285 this->createBackwardTransitions();
286 visittimesHelper.provideBackwardTransitions(*this->_backwardTransitions);
287 std::vector<ValueType> initialValues;
288 initialValues.reserve(bscc.size() - 1);
289 auto row = this->_transitionMatrix.getRow(proxyState);
290 auto entryIt = row.begin();
291 auto const entryItEnd = row.end();
292 for (auto state : bsccAsBitVector) {
293 if (entryIt != entryItEnd && state == entryIt->getColumn()) {
294 initialValues.push_back(entryIt->getValue());
295 ++entryIt;
296 } else {
297 initialValues.push_back(storm::utility::zero<ValueType>());
298 }
299 }
300 STORM_LOG_ASSERT(entryIt == entryItEnd || entryIt->getColumn() == proxyState, "Unexpected matrix row.");
301 auto evtEnv = env;
302 if (evtEnv.solver().isForceSoundness()) {
303 auto prec = evtEnv.solver().getPrecisionOfLinearEquationSolver(evtEnv.solver().getLinearEquationSolverType());
304 if (prec.first.is_initialized()) {
305 auto requiredPrecision = *prec.first;
306 // We need to adapt the precision. We focus on propagation of relative errors. Absolut precision is then also fine as we're computing with values in
307 // [0,1].
308 //
309 // We are going to normalize a vector of EVTs. Assuming that the EVT vector has been computed with relative precision eps, the sum of all
310 // EVTs still is eps-precise w.r.t. the exact (unknown) sum of the EVTs. Let x and y be the computed EVT of a given state and the sum of all
311 // computed EVTs, respectively. We have relative errors of delta_x = (x-x')/x' and delta_y = (y-y')/y' respectively, where x' and y' are the exact
312 // values. Note that x = x' * (1+delta_x) and |delta_x| <= eps.
313 //
314 // The relative error in the normalized value x/y = (x'/y')*((1+delta_x)/(1+delta_y)) = (x'/y')*(1+((delta_x-delta_y)/(1+\delta_y))) can be upper
315 // bounded by 2*eps/(1-eps). We set eps so that this term is equal to requiredPrecision
316 storm::RationalNumber eps = requiredPrecision / (storm::utility::convertNumber<storm::RationalNumber, uint64_t>(2) + requiredPrecision);
317 evtEnv.solver().setLinearEquationSolverPrecision(eps, prec.second);
318 }
319 }
320 auto visitingTimes = visittimesHelper.computeExpectedVisitingTimes(evtEnv, bsccAsBitVector, initialValues);
321 visitingTimes.push_back(storm::utility::one<ValueType>()); // Add the value for the proxy state
322 bsccAsBitVector.set(proxyState, true);
323
324 ValueType sumOfVisitingTimes = storm::utility::zero<ValueType>();
325 if (this->isContinuousTime()) {
326 auto resultIt = visitingTimes.begin();
327 for (auto state : bsccAsBitVector) {
328 *resultIt /= (*this->_exitRates)[state];
329 sumOfVisitingTimes += *resultIt;
330 ++resultIt;
331 }
332 } else {
333 sumOfVisitingTimes = std::accumulate(visitingTimes.begin(), visitingTimes.end(), storm::utility::zero<ValueType>());
334 }
335 storm::utility::vector::scaleVectorInPlace(visitingTimes, storm::utility::one<ValueType>() / sumOfVisitingTimes);
336 return visitingTimes;
337}
338
339template<typename ValueType>
341 Environment const& inputEnv, storm::storage::StronglyConnectedComponent const& bscc) {
342 // We want that the returned values are sorted properly. Let's assert that strongly connected components use a sorted container type
343 STORM_LOG_ASSERT(std::is_sorted(bscc.begin(), bscc.end()), "Expected that bsccs are sorted.");
344
345 // Prepare an environment for the underlying linear equation solver
346 auto env = inputEnv;
347 if (env.solver().getLinearEquationSolverType() == storm::solver::EquationSolverType::Topological) {
348 // Topological solver does not make any sense since the BSCC is connected.
349 env.solver().setLinearEquationSolverType(env.solver().topological().getUnderlyingEquationSolverType(),
350 env.solver().topological().isUnderlyingEquationSolverTypeSetFromDefault());
351 }
352
353 STORM_LOG_WARN_COND(!env.solver().isForceSoundness(),
354 "Sound computations are not properly implemented for this computation. You might get incorrect results.");
355
356 // Let A be ab auxiliary Matrix with A[s,s] = R(s,s) - r(s) & A[s,s'] = R(s,s') for s,s' in BSCC and s!=s'.
357 // We build and solve the equation system for
358 // x*A=0 & x_0+...+x_n=1 <=> A^t*x=0=x-x & x_0+...+x_n=1 <=> (1+A^t)*x = x & 1-x_0-...-x_n-1=x_n
359 // Then, x[i] will be the fraction of the time we are in state i.
360
361 // Get a mapping from global state indices to local ones as well as a bitvector containing states within the BSCC.
362 std::unordered_map<uint64_t, uint64_t> toLocalIndexMap;
363 storm::storage::BitVector bsccStates(this->_transitionMatrix.getRowCount(), false);
364 uint64_t localIndex = 0;
365 for (auto const& globalIndex : bscc) {
366 bsccStates.set(globalIndex, true);
367 toLocalIndexMap[globalIndex] = localIndex;
368 ++localIndex;
369 }
370
371 // Build the auxiliary Matrix A.
372 auto auxMatrix = this->_transitionMatrix.getSubmatrix(false, bsccStates, bsccStates, true); // add diagonal entries!
373 uint64_t row = 0;
374 for (auto const& globalIndex : bscc) {
375 ValueType rateAtState = this->_exitRates ? (*this->_exitRates)[globalIndex] : storm::utility::one<ValueType>();
376 for (auto& entry : auxMatrix.getRow(row)) {
377 if (entry.getColumn() == row) {
378 // This value is non-zero since we have a BSCC with more than one state
379 entry.setValue(rateAtState * (entry.getValue() - storm::utility::one<ValueType>()));
380 } else if (this->isContinuousTime()) {
381 entry.setValue(entry.getValue() * rateAtState);
382 }
383 }
384 ++row;
385 }
386 assert(row == auxMatrix.getRowCount());
387
388 // We need to consider A^t. This will not delete diagonal entries since they are non-zero.
389 auxMatrix = auxMatrix.transpose();
390
391 // Check whether we need the fixpoint characterization
393 bool isFixpointFormat = linearEquationSolverFactory.getEquationProblemFormat(env) == storm::solver::LinearEquationSolverProblemFormat::FixedPointSystem;
394 if (isFixpointFormat) {
395 // Add a 1 on the diagonal
396 for (row = 0; row < auxMatrix.getRowCount(); ++row) {
397 for (auto& entry : auxMatrix.getRow(row)) {
398 if (entry.getColumn() == row) {
399 entry.setValue(storm::utility::one<ValueType>() + entry.getValue());
400 }
401 }
402 }
403 }
404
405 // We now build the equation system matrix.
406 // We can drop the last row of A and add ones in this row instead to assert that the variables sum up to one
407 // Phase 1: replace the existing entries of the last row with ones
408 uint64_t col = 0;
409 uint64_t lastRow = auxMatrix.getRowCount() - 1;
410 for (auto& entry : auxMatrix.getRow(lastRow)) {
411 entry.setColumn(col);
412 if (isFixpointFormat) {
413 if (col == lastRow) {
414 entry.setValue(storm::utility::zero<ValueType>());
415 } else {
416 entry.setValue(-storm::utility::one<ValueType>());
417 }
418 } else {
419 entry.setValue(storm::utility::one<ValueType>());
420 }
421 ++col;
422 }
423 storm::storage::SparseMatrixBuilder<ValueType> builder(std::move(auxMatrix));
424 for (; col <= lastRow; ++col) {
425 if (isFixpointFormat) {
426 if (col != lastRow) {
427 builder.addNextValue(lastRow, col, -storm::utility::one<ValueType>());
428 }
429 } else {
430 builder.addNextValue(lastRow, col, storm::utility::one<ValueType>());
431 }
432 }
433
434 std::vector<ValueType> bsccEquationSystemRightSide(bscc.size(), storm::utility::zero<ValueType>());
435 bsccEquationSystemRightSide.back() = storm::utility::one<ValueType>();
436
437 // Create a linear equation solver
438 auto solver = linearEquationSolverFactory.create(env, builder.build());
439 solver->setBounds(storm::utility::zero<ValueType>(), storm::utility::one<ValueType>());
440 // Check solver requirements.
441 auto requirements = solver->getRequirements(env);
442 requirements.clearLowerBounds();
443 requirements.clearUpperBounds();
444 STORM_LOG_THROW(!requirements.hasEnabledCriticalRequirement(), storm::exceptions::UnmetRequirementException,
445 "Solver requirements " + requirements.getEnabledRequirementsAsString() + " not checked.");
446
447 std::vector<ValueType> steadyStateDistr(bscc.size(), storm::utility::one<ValueType>() / storm::utility::convertNumber<ValueType, uint64_t>(bscc.size()));
448 solver->solveEquations(env, steadyStateDistr, bsccEquationSystemRightSide);
449
450 // As a last step, we normalize these values to counter numerical inaccuracies a bit.
451 // This is only reasonable in non-exact mode.
452 if (!env.solver().isForceExact()) {
453 ValueType sum = std::accumulate(steadyStateDistr.begin(), steadyStateDistr.end(), storm::utility::zero<ValueType>());
454 storm::utility::vector::scaleVectorInPlace<ValueType, ValueType>(steadyStateDistr, storm::utility::one<ValueType>() / sum);
455 }
456
457 return steadyStateDistr;
458}
459
460template<typename ValueType>
462 Environment const& env, ValueGetter const& stateValuesGetter, ValueGetter const& actionValuesGetter,
464 // Get the LRA distribution for the provided bscc.
465 auto steadyStateDistr = computeSteadyStateDistrForBscc(env, bscc);
466
467 // Calculate final LRA Value
468 ValueType result = storm::utility::zero<ValueType>();
469 auto solIt = steadyStateDistr.begin();
470 for (auto const& globalState : bscc) {
471 if (this->isContinuousTime()) {
472 result += (*solIt) * (stateValuesGetter(globalState) + (*this->_exitRates)[globalState] * actionValuesGetter(globalState));
473 } else {
474 result += (*solIt) * (stateValuesGetter(globalState) + actionValuesGetter(globalState));
475 }
476 ++solIt;
477 }
478 assert(solIt == steadyStateDistr.end());
479
480 return std::pair<ValueType, std::vector<ValueType>>(std::move(result), std::move(steadyStateDistr));
481}
482
483template<typename ValueType>
484std::pair<storm::storage::SparseMatrix<ValueType>, std::vector<ValueType>> SparseDeterministicInfiniteHorizonHelper<ValueType>::buildSspMatrixVector(
485 std::vector<ValueType> const& bsccLraValues, std::vector<uint64_t> const& inputStateToBsccIndexMap, storm::storage::BitVector const& statesNotInComponent,
486 bool asEquationSystem) {
487 // Create SSP Matrix.
488 // In contrast to the version for nondeterministic models, we eliminate the auxiliary states representing each BSCC on the fly
489
490 // Probability mass that would lead to a BSCC will be considered in the rhs of the equation system
491 auto sspMatrix = this->_transitionMatrix.getSubmatrix(false, statesNotInComponent, statesNotInComponent, asEquationSystem);
492 if (asEquationSystem) {
493 sspMatrix.convertToEquationSystem();
494 }
495
496 // Create the SSP right-hand-side
497 std::vector<ValueType> rhs;
498 rhs.reserve(sspMatrix.getRowCount());
499 for (auto state : statesNotInComponent) {
500 ValueType stateValue = storm::utility::zero<ValueType>();
501 for (auto const& transition : this->_transitionMatrix.getRow(state)) {
502 if (!statesNotInComponent.get(transition.getColumn())) {
503 // This transition leads to a BSCC!
504 stateValue += transition.getValue() * bsccLraValues[inputStateToBsccIndexMap[transition.getColumn()]];
505 }
506 }
507 rhs.push_back(std::move(stateValue));
508 }
509
510 return std::make_pair(std::move(sspMatrix), std::move(rhs));
511}
512
513template<typename ValueType>
515 std::vector<ValueType> const& componentLraValues) {
516 STORM_LOG_ASSERT(this->_longRunComponentDecomposition != nullptr, "Decomposition not computed, yet.");
517
518 // For fast transition rewriting, we build a mapping from the input state indices to the state indices of a new transition matrix
519 // which redirects all transitions leading to a former BSCC state to a new (imaginary) auxiliary state.
520 // Each auxiliary state gets assigned the value of that BSCC and we compute expected rewards (aka stochastic shortest path, SSP) on that new system.
521 // For efficiency reasons, we actually build the system where the auxiliary states are already eliminated.
522
523 // First gather the states that are part of a component
524 // and create a mapping from states that lie in a component to the corresponding component index.
525 storm::storage::BitVector statesInComponents(this->_transitionMatrix.getRowGroupCount());
526 std::vector<uint64_t> stateIndexMap(this->_transitionMatrix.getRowGroupCount(), std::numeric_limits<uint64_t>::max());
527 for (uint64_t currentComponentIndex = 0; currentComponentIndex < this->_longRunComponentDecomposition->size(); ++currentComponentIndex) {
528 for (auto const& element : (*this->_longRunComponentDecomposition)[currentComponentIndex]) {
529 uint64_t state = internal::getComponentElementState(element);
530 statesInComponents.set(state);
531 stateIndexMap[state] = currentComponentIndex;
532 }
533 }
534 // Map the non-component states to their index in the SSP. Note that the order of these states will be preserved.
535 uint64_t numberOfNonComponentStates = 0;
536 storm::storage::BitVector statesNotInComponent = ~statesInComponents;
537 for (auto nonComponentState : statesNotInComponent) {
538 stateIndexMap[nonComponentState] = numberOfNonComponentStates;
539 ++numberOfNonComponentStates;
540 }
541
542 // The next step is to create the equation system solving the SSP (unless the whole system consists of BSCCs)
543 std::vector<ValueType> sspValues;
544 if (numberOfNonComponentStates > 0) {
546 bool isEqSysFormat = linearEquationSolverFactory.getEquationProblemFormat(env) == storm::solver::LinearEquationSolverProblemFormat::EquationSystem;
547 auto sspMatrixVector = buildSspMatrixVector(componentLraValues, stateIndexMap, statesNotInComponent, isEqSysFormat);
548 std::unique_ptr<storm::solver::LinearEquationSolver<ValueType>> solver = linearEquationSolverFactory.create(env, sspMatrixVector.first);
549 auto lowerUpperBounds = std::minmax_element(componentLraValues.begin(), componentLraValues.end());
550 solver->setBounds(*lowerUpperBounds.first, *lowerUpperBounds.second);
551 // Check solver requirements
552 auto requirements = solver->getRequirements(env);
553 requirements.clearUpperBounds();
554 requirements.clearLowerBounds();
555 STORM_LOG_THROW(!requirements.hasEnabledCriticalRequirement(), storm::exceptions::UnmetRequirementException,
556 "Solver requirements " + requirements.getEnabledRequirementsAsString() + " not checked.");
557 sspValues.assign(sspMatrixVector.first.getRowCount(),
558 (*lowerUpperBounds.first + *lowerUpperBounds.second) / storm::utility::convertNumber<ValueType, uint64_t>(2));
559 solver->solveEquations(env, sspValues, sspMatrixVector.second);
560 }
561
562 // Prepare result vector.
563 std::vector<ValueType> result(this->_transitionMatrix.getRowGroupCount());
564 for (uint64_t state = 0; state < stateIndexMap.size(); ++state) {
565 if (statesNotInComponent.get(state)) {
566 result[state] = sspValues[stateIndexMap[state]];
567 } else {
568 result[state] = componentLraValues[stateIndexMap[state]];
569 }
570 }
571 return result;
572}
573
574template<typename ValueType>
576 createDecomposition();
577 STORM_LOG_THROW(this->_longRunComponentDecomposition->size() <= 1, storm::exceptions::InvalidOperationException, "");
578 return computeLongRunAverageStateDistribution(env, [](uint64_t) { return storm::utility::zero<ValueType>(); });
579}
580
581template<typename ValueType>
583 uint64_t const& initialState) {
584 STORM_LOG_ASSERT(initialState < this->_transitionMatrix.getRowGroupCount(),
585 "Invlid initial state index: " << initialState << ". Have only " << this->_transitionMatrix.getRowGroupCount() << " states.");
586 return computeLongRunAverageStateDistribution(env, [&initialState](uint64_t stateIndex) {
587 return initialState == stateIndex ? storm::utility::one<ValueType>() : storm::utility::zero<ValueType>();
588 });
589}
590
591template<typename ValueType>
593 Environment const& env, ValueGetter const& initialDistributionGetter) {
594 createDecomposition();
595
596 Environment subEnv = env;
597 if (subEnv.solver().isForceSoundness()) {
598 // We need to adapt the precision. We focus on propagation of relative errors. Absolut precision is then also fine as we're dealing with values in
599 // [0,1]
601 if (prec.first.is_initialized()) {
602 auto requiredPrecision = *prec.first;
603 // We are going to multiply two numbers x and y that each have relative errors of delta_x = (x-x')/x' and delta_y = (y-y')/y' respectively.
604 // Here, x' and y' are the exact values. Note that x = x' * (1+delta_x) and |delta_x| <= eps
605 // The result x*y= x' * y' * (1+delta_x) * (1+delta_y) will have a relative error of delta_x + delta_y + delta_x*\delta_y <= (2eps * eps^2)
606 // We set eps such that (2eps * eps^2) <= requiredPrecision
607 storm::RationalNumber eps = storm::utility::sqrt<RationalNumber>(storm::utility::one<storm::RationalNumber>() + requiredPrecision) -
608 storm::utility::one<storm::RationalNumber>();
609 subEnv.solver().setLinearEquationSolverPrecision(eps, prec.second);
610 STORM_LOG_INFO("Precision for BSCC reachability and BSCC steady state distribution analysis set to " << storm::utility::convertNumber<double>(eps)
611 << ".");
612 }
613 }
614 // Compute for each BSCC get the probability with which we reach that BSCC
615 auto bsccReachProbs = computeBsccReachabilityProbabilities(subEnv, initialDistributionGetter);
616 // We are now ready to compute the resulting lra distribution
617 std::vector<ValueType> steadyStateDistr(this->_transitionMatrix.getRowGroupCount(), storm::utility::zero<ValueType>());
618 for (uint64_t currentComponentIndex = 0; currentComponentIndex < this->_longRunComponentDecomposition->size(); ++currentComponentIndex) {
619 auto const& component = (*this->_longRunComponentDecomposition)[currentComponentIndex];
620 // Compute distribution for current bscc
621 auto bsccDistr = this->computeSteadyStateDistrForBscc(subEnv, component);
622 // Scale with probability to reach that bscc
623 auto const& scalingFactor = bsccReachProbs[currentComponentIndex];
624 if (!storm::utility::isOne(scalingFactor)) {
625 storm::utility::vector::scaleVectorInPlace(bsccDistr, scalingFactor);
626 }
627 // Set the values in the result vector
628 auto bsccDistrIt = bsccDistr.begin();
629 for (auto const& element : component) {
630 uint64_t state = internal::getComponentElementState(element);
631 steadyStateDistr[state] = *bsccDistrIt;
632 ++bsccDistrIt;
633 }
634 STORM_LOG_ASSERT(bsccDistrIt == bsccDistr.end(), "Unexpected number of entries in bscc distribution");
635 }
636 return steadyStateDistr;
637}
638
639template<typename ValueType>
641 std::vector<ValueType> const& toBsccProbabilities) {
643}
644
645template<>
647 std::vector<storm::RationalFunction> const& /*toBsccProbabilities*/) {
648 STORM_LOG_THROW(false, storm::exceptions::NotSupportedException,
649 "Computing upper bounds for expected visiting times over rational functions is not supported.");
650}
651
652template<typename ValueType>
654 ValueGetter const& initialDistributionGetter) {
655 STORM_LOG_ASSERT(this->_longRunComponentDecomposition != nullptr, "Decomposition not computed, yet.");
656
657 // Compute for each BSCC get the probability with which we reach that BSCC
658 std::vector<ValueType> bsccReachProbs;
659 if (auto numBSCCs = this->_longRunComponentDecomposition->size(); numBSCCs <= 1) {
660 STORM_LOG_ASSERT(numBSCCs == 1, "Found 0 BSCCs in a Markov chain. This should not be possible.");
661 bsccReachProbs = std::vector<ValueType>({storm::utility::one<ValueType>()});
662 } else {
664 bsccReachProbs = computeBsccReachabilityProbabilitiesClassic(env, initialDistributionGetter);
665 } else {
666 bsccReachProbs = computeBsccReachabilityProbabilitiesEVTs(env, initialDistributionGetter);
667 }
668 }
669
670 // As a last step, we normalize these values to counter inaccuracies a bit.
671 // This is only reasonable in non-exact mode and can invalidate accuracy guarantees in sound mode.
672 if (!env.solver().isForceExact() && !env.solver().isForceSoundness()) {
673 ValueType sum = std::accumulate(bsccReachProbs.begin(), bsccReachProbs.end(), storm::utility::zero<ValueType>());
674 storm::utility::vector::scaleVectorInPlace<ValueType, ValueType>(bsccReachProbs, storm::utility::one<ValueType>() / sum);
675 }
676 return bsccReachProbs;
677}
678
679template<typename ValueType>
681 Environment const& env, ValueGetter const& initialDistributionGetter) {
682 // Solve a linear equation system for each BSCC
683
684 // Get the states that do not lie on any BSCC
685 storm::storage::BitVector nonBsccStates(this->_transitionMatrix.getRowCount(), true);
686 for (uint64_t currentComponentIndex = 0; currentComponentIndex < this->_longRunComponentDecomposition->size(); ++currentComponentIndex) {
687 for (auto const& element : (*this->_longRunComponentDecomposition)[currentComponentIndex]) {
688 nonBsccStates.set(internal::getComponentElementState(element), false);
689 }
690 }
691 bool const hasNonBsccStates = !nonBsccStates.empty();
692
693 // Set-up a linear equation solver (unless all states are on a BSCC)
694 std::unique_ptr<storm::solver::LinearEquationSolver<ValueType>> solver;
695 if (hasNonBsccStates) {
697 bool isEquationSystemFormat =
699 auto subMatrix = this->_transitionMatrix.getSubmatrix(false, nonBsccStates, nonBsccStates, isEquationSystemFormat);
700 if (isEquationSystemFormat) {
701 subMatrix.convertToEquationSystem();
702 }
703 solver = linearEquationSolverFactory.create(env, std::move(subMatrix));
704 solver->setBounds(storm::utility::zero<ValueType>(), storm::utility::one<ValueType>());
705 // Check solver requirements.
706 auto requirements = solver->getRequirements(env);
707 requirements.clearLowerBounds();
708 requirements.clearUpperBounds();
709 STORM_LOG_THROW(!requirements.hasEnabledCriticalRequirement(), storm::exceptions::UnmetRequirementException,
710 "Solver requirements " + requirements.getEnabledRequirementsAsString() + " not checked.");
711 solver->setCachingEnabled(true);
712 }
713
714 // Run over all BSCCs
715 std::vector<ValueType> bsccReachProbs(this->_longRunComponentDecomposition->size(), storm::utility::zero<ValueType>());
716 for (uint64_t currentComponentIndex = 0; currentComponentIndex < this->_longRunComponentDecomposition->size(); ++currentComponentIndex) {
717 auto const& bscc = (*this->_longRunComponentDecomposition)[currentComponentIndex];
718 auto& bsccVal = bsccReachProbs[currentComponentIndex];
719 // Deal with initial states within the BSCC
720 for (auto const& element : (*this->_longRunComponentDecomposition)[currentComponentIndex]) {
721 uint64_t state = internal::getComponentElementState(element);
722 bsccVal += initialDistributionGetter(state);
723 }
724 // Add reachability probabilities for initial states outside of this BSCC
725 if (hasNonBsccStates) {
726 storm::storage::BitVector bsccAsBitVector(this->_transitionMatrix.getColumnCount(), false);
727 bsccAsBitVector.set(bscc.begin(), bscc.end(), true);
728 // set up and solve the equation system for this BSCC
729 std::vector<ValueType> eqSysRhs;
730 eqSysRhs.reserve(nonBsccStates.getNumberOfSetBits());
731 for (auto state : nonBsccStates) {
732 eqSysRhs.push_back(this->_transitionMatrix.getConstrainedRowSum(state, bsccAsBitVector));
733 }
734 std::vector<ValueType> eqSysSolution(eqSysRhs.size());
735 solver->solveEquations(env, eqSysSolution, eqSysRhs);
736 // Sum up reachability probabilities over initial states
737 uint64_t subsysState = 0;
738 for (auto globalState : nonBsccStates) {
739 bsccVal += initialDistributionGetter(globalState) * eqSysSolution[subsysState];
740 ++subsysState;
741 }
742 }
743 }
744 return bsccReachProbs;
745}
746
747template<typename ValueType>
749 Environment const& env, ValueGetter const& initialDistributionGetter) {
750 // Get the expected number of times we visit each non-BSCC state
751 // See https://arxiv.org/abs/2401.10638 for more information.
752 // We deliberately exclude the exit rates here as we want to make this computation on the induced DTMC to get the expected number of times
753 // that a successor state is chosen probabilistically.
754 auto visittimesHelper = SparseDeterministicVisitingTimesHelper<ValueType>(this->_transitionMatrix);
755 this->createBackwardTransitions();
756 visittimesHelper.provideBackwardTransitions(*this->_backwardTransitions);
757 auto expVisitTimes = visittimesHelper.computeExpectedVisitingTimes(env, initialDistributionGetter);
758 // Then use the expected visiting times to compute BSCC reachability probabilities
759 storm::storage::BitVector nonBsccStates(this->_transitionMatrix.getRowCount(), true);
760 for (uint64_t currentComponentIndex = 0; currentComponentIndex < this->_longRunComponentDecomposition->size(); ++currentComponentIndex) {
761 for (auto const& element : (*this->_longRunComponentDecomposition)[currentComponentIndex]) {
762 nonBsccStates.set(internal::getComponentElementState(element), false);
763 }
764 }
765 std::vector<ValueType> bsccReachProbs(this->_longRunComponentDecomposition->size(), storm::utility::zero<ValueType>());
766 for (uint64_t currentComponentIndex = 0; currentComponentIndex < this->_longRunComponentDecomposition->size(); ++currentComponentIndex) {
767 auto& bsccVal = bsccReachProbs[currentComponentIndex];
768 for (auto const& element : (*this->_longRunComponentDecomposition)[currentComponentIndex]) {
769 uint64_t state = internal::getComponentElementState(element);
770 bsccVal += initialDistributionGetter(state);
771 for (auto const& pred : this->_backwardTransitions->getRow(state)) {
772 if (nonBsccStates.get(pred.getColumn())) {
773 bsccVal += pred.getValue() * expVisitTimes[pred.getColumn()];
774 }
775 }
776 }
777 }
778 return bsccReachProbs;
779}
780
784
785} // namespace helper
786} // namespace modelchecker
787} // namespace storm
SolverEnvironment & solver()
ModelCheckerEnvironment & modelchecker()
storm::RationalNumber const & getAperiodicFactor() const
storm::RationalNumber const & getPrecision() const
storm::solver::LraMethod const & getDetLraMethod() const
SteadyStateDistributionAlgorithm getSteadyStateDistributionAlgorithm() const
void setLinearEquationSolverType(storm::solver::EquationSolverType const &value, bool isSetFromDefault=false)
storm::solver::EquationSolverType const & getLinearEquationSolverType() const
void setLinearEquationSolverPrecision(boost::optional< storm::RationalNumber > const &newPrecision, boost::optional< bool > const &relativePrecision=boost::none)
LongRunAverageSolverEnvironment & lra()
std::pair< boost::optional< storm::RationalNumber >, boost::optional< bool > > getPrecisionOfLinearEquationSolver(storm::solver::EquationSolverType const &solverType) const
static std::vector< ValueType > computeUpperBoundOnExpectedVisitingTimes(storm::storage::SparseMatrix< ValueType > const &transitionMatrix, storm::storage::SparseMatrix< ValueType > const &backwardTransitions, std::vector< ValueType > const &oneStepTargetProbabilities)
Computes for each state an upper bound for the maximal expected times each state is visited.
Helper class for model checking queries that depend on the long run behavior of the (nondeterministic...
SparseDeterministicInfiniteHorizonHelper(storm::storage::SparseMatrix< ValueType > const &transitionMatrix)
Initializes the helper for a discrete time model (i.e.
ValueType computeLraForBsccVi(Environment const &env, ValueGetter const &stateValuesGetter, ValueGetter const &actionValuesGetter, storm::storage::StronglyConnectedComponent const &bscc)
As computeLraForComponent but uses value iteration as a solution method (independent of what is set i...
std::pair< ValueType, std::vector< ValueType > > computeLraForBsccGainBias(Environment const &env, ValueGetter const &stateValuesGetter, ValueGetter const &actionValuesGetter, storm::storage::StronglyConnectedComponent const &bscc)
As computeLraForComponent but solves a linear equation system encoding gain and bias (independent of ...
virtual ValueType computeLraForComponent(Environment const &env, ValueGetter const &stateValuesGetter, ValueGetter const &actionValuesGetter, storm::storage::StronglyConnectedComponent const &component) override
std::pair< ValueType, std::vector< ValueType > > computeLraForBsccSteadyStateDistr(Environment const &env, ValueGetter const &stateValuesGetter, ValueGetter const &actionValuesGetter, storm::storage::StronglyConnectedComponent const &bscc)
As computeLraForComponent but does the computation by computing the long run average (steady state) d...
std::pair< storm::storage::SparseMatrix< ValueType >, std::vector< ValueType > > buildSspMatrixVector(std::vector< ValueType > const &bsccLraValues, std::vector< uint64_t > const &inputStateToBsccIndexMap, storm::storage::BitVector const &statesNotInComponent, bool asEquationSystem)
std::pair< bool, ValueType > computeLraForTrivialBscc(Environment const &env, ValueGetter const &stateValuesGetter, ValueGetter const &actionValuesGetter, storm::storage::StronglyConnectedComponent const &bscc)
SparseInfiniteHorizonHelper< ValueType, true >::ValueGetter ValueGetter
Function mapping from indices to values.
virtual std::vector< ValueType > buildAndSolveSsp(Environment const &env, std::vector< ValueType > const &mecLraValues) override
std::vector< ValueType > computeLongRunAverageStateDistribution(Environment const &env)
Computes the long run average state distribution, i.e., a vector that assigns for each state s the av...
std::vector< ValueType > computeSteadyStateDistrForBsccEVTs(Environment const &env, storm::storage::StronglyConnectedComponent const &bscc)
std::vector< ValueType > computeSteadyStateDistrForBscc(Environment const &env, storm::storage::StronglyConnectedComponent const &bscc)
Computes the long run average (steady state) distribution for the given BSCC.
std::vector< ValueType > computeBsccReachabilityProbabilitiesClassic(Environment const &env, ValueGetter const &initialDistributionGetter)
std::vector< ValueType > computeBsccReachabilityProbabilitiesEVTs(Environment const &env, ValueGetter const &initialDistributionGetter)
std::vector< ValueType > computeBsccReachabilityProbabilities(Environment const &env, ValueGetter const &initialDistributionGetter)
Computes for each BSCC the probability to reach that SCC assuming the given distribution over initial...
std::vector< ValueType > computeSteadyStateDistrForBsccEqSys(Environment const &env, storm::storage::StronglyConnectedComponent const &bscc)
Helper class for computing for each state the expected number of times to visit that state assuming a...
Helper class for model checking queries that depend on the long run behavior of the (nondeterministic...
Helper class that performs iterations of the value iteration method.
Definition LraViHelper.h:42
virtual std::unique_ptr< LinearEquationSolver< ValueType > > create(Environment const &env) const override
Creates an equation solver with the current settings, but without a matrix.
virtual LinearEquationSolverProblemFormat getEquationProblemFormat(Environment const &env) const
Retrieves the problem format that the solver expects if it was created with the current settings.
A bit vector that is internally represented as a vector of 64-bit values.
Definition BitVector.h:16
const_reverse_iterator rbegin() const
Returns a reverse iterator to the indices of the set bits in the bit vector.
bool empty() const
Retrieves whether no bits are set to true in this bit vector.
uint64_t getNumberOfSetBits() const
Returns the number of bits that are set to true in this bit vector.
void set(uint64_t index, bool value=true)
Sets the given truth value at the given index.
bool get(uint64_t index) const
Retrieves the truth value of the bit at the given index and performs a bound check.
A class that can be used to build a sparse matrix by adding value by value.
void addNextValue(index_type row, index_type column, value_type const &value)
Sets the matrix entry at the given row and column to the given value.
void addDiagonalEntry(index_type row, ValueType const &value)
Makes sure that a diagonal entry will be inserted at the given row.
SparseMatrix< value_type > build(index_type overriddenRowCount=0, index_type overriddenColumnCount=0, index_type overriddenRowGroupCount=0)
A class that holds a possibly non-square matrix in the compressed row storage format.
std::size_t size() const
Retrieves the number of states in this SCC.
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 a strongly connected component, i.e., a set of states such that every state can...
#define STORM_LOG_INFO(message)
Definition logging.h:24
#define STORM_LOG_TRACE(message)
Definition logging.h:12
#define STORM_LOG_ASSERT(cond, message)
Definition macros.h:11
#define STORM_LOG_WARN_COND(cond, message)
Definition macros.h:38
#define STORM_LOG_THROW(cond, exception, message)
Definition macros.h:30
uint64_t getComponentElementState(typename storm::storage::StronglyConnectedComponent::value_type const &element)
Auxiliary functions that deal with the different kinds of components (MECs on potentially nondetermin...
uint64_t const * getComponentElementChoicesBegin(typename storm::storage::StronglyConnectedComponent::value_type const &element)
std::vector< ValueType > computeUpperBoundsForExpectedVisitingTimes(storm::storage::SparseMatrix< ValueType > const &nonBsccMatrix, std::vector< ValueType > const &toBsccProbabilities)
std::string toString(GurobiSolverMethod const &method)
Yields a string representation of the GurobiSolverMethod.
void scaleVectorInPlace(std::vector< ValueType1 > &target, ValueType2 const &factor)
Multiplies each element of the given vector with the given factor and writes the result into the vect...
Definition vector.h:451
bool isOne(ValueType const &a)
Definition constants.cpp:33
StronglyConnectedComponentDecompositionOptions & onlyBottomSccs(bool value=true)
Sets if only bottom SCCs, i.e. SCCs in which all states have no way of leaving the SCC),...