Storm
A Modern Probabilistic Model Checker
Loading...
Searching...
No Matches
SymbolicMdpPrctlHelper.cpp
Go to the documentation of this file.
2
4
8
10
12#include "storm/utility/graph.h"
13
15
18
22
23namespace storm {
24namespace modelchecker {
25namespace helper {
26
28
29template<storm::dd::DdType DdType, typename ValueType>
31 storm::dd::Add<DdType, ValueType> const& transitionMatrix, storm::dd::Bdd<DdType> const& maybeStates,
32 storm::dd::Bdd<DdType> const& targetStates) {
34
36 result = storm::utility::graph::computeSchedulerProbGreater0E(model, transitionMatrix.notZero(), maybeStates, targetStates);
37 } else if (type == EquationSystemType::ExpectedRewards) {
38 result = storm::utility::graph::computeSchedulerProb1E(model, transitionMatrix.notZero(), maybeStates, targetStates, maybeStates || targetStates);
39 }
40
41 return result;
42}
43
44template<storm::dd::DdType DdType, typename ValueType>
47 storm::dd::Add<DdType, ValueType> const& transitionMatrix, storm::dd::Bdd<DdType> const& maybeStates, storm::dd::Bdd<DdType> const& statesWithProbability1,
48 boost::optional<storm::dd::Add<DdType, ValueType>> const& startValues) {
49 // Create the matrix and the vector for the equation system.
50 storm::dd::Add<DdType, ValueType> maybeStatesAdd = maybeStates.template toAdd<ValueType>();
51
52 // Start by cutting away all rows that do not belong to maybe states. Note that this leaves columns targeting
53 // non-maybe states in the matrix.
54 storm::dd::Add<DdType, ValueType> submatrix = transitionMatrix * maybeStatesAdd;
55
56 // Then compute the vector that contains the one-step probabilities to a state with probability 1 for all
57 // maybe states.
58 storm::dd::Add<DdType, ValueType> prob1StatesAsColumn = statesWithProbability1.template toAdd<ValueType>();
59 prob1StatesAsColumn = prob1StatesAsColumn.swapVariables(model.getRowColumnMetaVariablePairs());
60 storm::dd::Add<DdType, ValueType> subvector = submatrix * prob1StatesAsColumn;
61 subvector = subvector.sumAbstract(model.getColumnVariables());
62
63 // Finally cut away all columns targeting non-maybe states.
64 submatrix *= maybeStatesAdd.swapVariables(model.getRowColumnMetaVariablePairs());
65
66 // Now solve the resulting equation system.
68 std::unique_ptr<storm::solver::SymbolicMinMaxLinearEquationSolver<DdType, ValueType>> solver =
69 linearEquationSolverFactory.create(submatrix, maybeStates, model.getIllegalMask() && maybeStates, model.getRowVariables(), model.getColumnVariables(),
71
72 if (storm::solver::minimize(dir)) {
73 solver->setHasUniqueSolution(true);
74 }
75
76 // Check requirements of solver.
77 storm::solver::MinMaxLinearEquationSolverRequirements requirements = solver->getRequirements(env, dir);
78 boost::optional<storm::dd::Bdd<DdType>> initialScheduler;
79 if (requirements.hasEnabledRequirement()) {
80 if (requirements.validInitialScheduler()) {
81 STORM_LOG_DEBUG("Computing valid scheduler, because the solver requires it.");
82 initialScheduler = computeValidSchedulerHint(EquationSystemType::UntilProbabilities, model, transitionMatrix, maybeStates, statesWithProbability1);
83 requirements.clearValidInitialScheduler();
84 }
85 requirements.clearBounds();
86 if (requirements.uniqueSolution()) {
87 // Check whether there are end components
88 if (storm::utility::graph::performProb0E(model, transitionMatrix.notZero(), maybeStates, !maybeStates && model.getReachableStates()).isZero()) {
89 requirements.clearUniqueSolution();
90 }
91 }
92 STORM_LOG_THROW(!requirements.hasEnabledCriticalRequirement(), storm::exceptions::UncheckedRequirementException,
93 "Solver requirements " + requirements.getEnabledRequirementsAsString() + " not checked.");
94 }
95 if (initialScheduler) {
96 solver->setInitialScheduler(initialScheduler.get());
97 }
98 solver->setBounds(storm::utility::zero<ValueType>(), storm::utility::one<ValueType>());
99 solver->setRequirementsChecked();
100
102 solver->solveEquations(env, dir, startValues ? startValues.get() : maybeStatesAdd.getDdManager().template getAddZero<ValueType>(), subvector);
103
105 model.getReachableStates(), statesWithProbability1.template toAdd<ValueType>() + result));
106}
107
108template<storm::dd::DdType DdType, typename ValueType>
111 storm::dd::Add<DdType, ValueType> const& transitionMatrix, storm::dd::Bdd<DdType> const& phiStates, storm::dd::Bdd<DdType> const& psiStates,
112 bool qualitative, boost::optional<storm::dd::Add<DdType, ValueType>> const& startValues) {
113 // We need to identify the states which have to be taken out of the matrix, i.e. all states that have
114 // probability 0 and 1 of satisfying the until-formula.
115 std::pair<storm::dd::Bdd<DdType>, storm::dd::Bdd<DdType>> statesWithProbability01;
116 if (dir == OptimizationDirection::Minimize) {
117 statesWithProbability01 = storm::utility::graph::performProb01Min(model, phiStates, psiStates);
118 } else {
119 statesWithProbability01 = storm::utility::graph::performProb01Max(model, phiStates, psiStates);
120 }
121
122 storm::dd::Bdd<DdType> maybeStates = !statesWithProbability01.first && !statesWithProbability01.second && model.getReachableStates();
123 STORM_LOG_INFO("Preprocessing: " << statesWithProbability01.first.getNonZeroCount() << " states with probability 0, "
124 << statesWithProbability01.second.getNonZeroCount() << " with probability 1 (" << maybeStates.getNonZeroCount()
125 << " states remaining).");
126
127 // Check whether we need to compute exact probabilities for some states.
128 if (qualitative) {
129 // Set the values for all maybe-states to 0.5 to indicate that their probability values are neither 0 nor 1.
131 model.getReachableStates(),
132 statesWithProbability01.second.template toAdd<ValueType>() +
133 maybeStates.template toAdd<ValueType>() * model.getManager().getConstant(storm::utility::convertNumber<ValueType>(0.5))));
134 } else {
135 // If there are maybe states, we need to solve an equation system.
136 if (!maybeStates.isZero()) {
137 return computeUntilProbabilities(env, dir, model, transitionMatrix, maybeStates, statesWithProbability01.second,
138 startValues ? maybeStates.ite(startValues.get(), model.getManager().template getAddZero<ValueType>())
139 : model.getManager().template getAddZero<ValueType>());
140 } else {
142 model.getReachableStates(), statesWithProbability01.second.template toAdd<ValueType>()));
143 }
144 }
145}
146
147template<storm::dd::DdType DdType, typename ValueType>
150 storm::dd::Add<DdType, ValueType> const& transitionMatrix, storm::dd::Bdd<DdType> const& psiStates, bool qualitative) {
151 std::unique_ptr<CheckResult> result =
152 computeUntilProbabilities(env, dir == OptimizationDirection::Minimize ? OptimizationDirection::Maximize : OptimizationDirection::Minimize, model,
153 transitionMatrix, model.getReachableStates(), !psiStates && model.getReachableStates(), qualitative);
154 result->asQuantitativeCheckResult<ValueType>().oneMinus();
155 return result;
156}
157
158template<storm::dd::DdType DdType, typename ValueType>
161 storm::dd::Add<DdType, ValueType> const& transitionMatrix, storm::dd::Bdd<DdType> const& nextStates) {
162 storm::dd::Add<DdType, ValueType> result = (transitionMatrix * nextStates.swapVariables(model.getRowColumnMetaVariablePairs()).template toAdd<ValueType>())
163 .sumAbstract(model.getColumnVariables());
164 if (dir == OptimizationDirection::Minimize) {
165 result = (result + model.getIllegalMask().template toAdd<ValueType>()).minAbstract(model.getNondeterminismVariables());
166 } else {
167 result = result.maxAbstract(model.getNondeterminismVariables());
168 }
169 return std::unique_ptr<CheckResult>(new SymbolicQuantitativeCheckResult<DdType, ValueType>(model.getReachableStates(), result));
170}
171
172template<storm::dd::DdType DdType, typename ValueType>
175 storm::dd::Add<DdType, ValueType> const& transitionMatrix, storm::dd::Bdd<DdType> const& phiStates, storm::dd::Bdd<DdType> const& psiStates,
176 uint_fast64_t stepBound) {
177 // We need to identify the states which have to be taken out of the matrix, i.e. all states that have
178 // probability 0 or 1 of satisfying the until-formula.
179 storm::dd::Bdd<DdType> statesWithProbabilityGreater0;
180 if (dir == OptimizationDirection::Minimize) {
181 statesWithProbabilityGreater0 = storm::utility::graph::performProbGreater0A(model, transitionMatrix.notZero(), phiStates, psiStates);
182 } else {
183 statesWithProbabilityGreater0 = storm::utility::graph::performProbGreater0E(model, transitionMatrix.notZero(), phiStates, psiStates);
184 }
185 storm::dd::Bdd<DdType> maybeStates = statesWithProbabilityGreater0 && !psiStates && model.getReachableStates();
186
187 // If there are maybe states, we need to perform matrix-vector multiplications.
188 if (!maybeStates.isZero()) {
189 // Create the matrix and the vector for the equation system.
190 storm::dd::Add<DdType, ValueType> maybeStatesAdd = maybeStates.template toAdd<ValueType>();
191
192 // Start by cutting away all rows that do not belong to maybe states. Note that this leaves columns targeting
193 // non-maybe states in the matrix.
194 storm::dd::Add<DdType, ValueType> submatrix = transitionMatrix * maybeStatesAdd;
195
196 // Then compute the vector that contains the one-step probabilities to a state with probability 1 for all
197 // maybe states.
198 storm::dd::Add<DdType, ValueType> prob1StatesAsColumn = psiStates.template toAdd<ValueType>().swapVariables(model.getRowColumnMetaVariablePairs());
199 storm::dd::Add<DdType, ValueType> subvector = (submatrix * prob1StatesAsColumn).sumAbstract(model.getColumnVariables());
200
201 // Finally cut away all columns targeting non-maybe states.
202 submatrix *= maybeStatesAdd.swapVariables(model.getRowColumnMetaVariablePairs());
203
205 std::unique_ptr<storm::solver::SymbolicMinMaxLinearEquationSolver<DdType, ValueType>> solver =
206 linearEquationSolverFactory.create(submatrix, maybeStates, model.getIllegalMask() && maybeStates, model.getRowVariables(),
208 storm::dd::Add<DdType, ValueType> result = solver->multiply(dir, model.getManager().template getAddZero<ValueType>(), &subvector, stepBound);
209
211 model.getReachableStates(), psiStates.template toAdd<ValueType>() + result));
212 } else {
213 return std::unique_ptr<CheckResult>(
214 new storm::modelchecker::SymbolicQuantitativeCheckResult<DdType, ValueType>(model.getReachableStates(), psiStates.template toAdd<ValueType>()));
215 }
216}
217
218template<storm::dd::DdType DdType, typename ValueType>
221 storm::dd::Add<DdType, ValueType> const& transitionMatrix, RewardModelType const& rewardModel, uint_fast64_t stepBound) {
222 // Only compute the result if the model has at least one reward this->getModel().
223 STORM_LOG_THROW(rewardModel.hasStateRewards(), storm::exceptions::InvalidPropertyException,
224 "Computing instantaneous rewards for a reward model that does not define any state-rewards. The result is trivially 0.");
225
226 // Perform the matrix-vector multiplication.
228 std::unique_ptr<storm::solver::SymbolicMinMaxLinearEquationSolver<DdType, ValueType>> solver =
229 linearEquationSolverFactory.create(transitionMatrix, model.getReachableStates(), model.getIllegalMask(), model.getRowVariables(),
231 storm::dd::Add<DdType, ValueType> result = solver->multiply(dir, rewardModel.getStateRewardVector(), nullptr, stepBound);
232
233 return std::unique_ptr<CheckResult>(new SymbolicQuantitativeCheckResult<DdType, ValueType>(model.getReachableStates(), result));
234}
235
236template<storm::dd::DdType DdType, typename ValueType>
239 storm::dd::Add<DdType, ValueType> const& transitionMatrix, RewardModelType const& rewardModel, uint_fast64_t stepBound) {
240 // Only compute the result if the model has at least one reward this->getModel().
241 STORM_LOG_THROW(!rewardModel.empty(), storm::exceptions::InvalidPropertyException, "Reward model for formula is empty. Skipping formula.");
242
243 // Compute the reward vector to add in each step based on the available reward models.
244 storm::dd::Add<DdType, ValueType> totalRewardVector = rewardModel.getTotalRewardVector(transitionMatrix, model.getColumnVariables());
245
246 // Perform the matrix-vector multiplication.
248 std::unique_ptr<storm::solver::SymbolicMinMaxLinearEquationSolver<DdType, ValueType>> solver =
249 linearEquationSolverFactory.create(model.getTransitionMatrix(), model.getReachableStates(), model.getIllegalMask(), model.getRowVariables(),
251 storm::dd::Add<DdType, ValueType> result = solver->multiply(dir, model.getManager().template getAddZero<ValueType>(), &totalRewardVector, stepBound);
252
253 return std::unique_ptr<CheckResult>(new SymbolicQuantitativeCheckResult<DdType, ValueType>(model.getReachableStates(), result));
254}
255
256template<storm::dd::DdType DdType, typename ValueType>
259 storm::dd::Add<DdType, ValueType> const& transitionMatrix, storm::dd::Bdd<DdType> const& transitionMatrixBdd, RewardModelType const& rewardModel,
260 storm::dd::Bdd<DdType> const& maybeStates, storm::dd::Bdd<DdType> const& targetStates, storm::dd::Bdd<DdType> const& infinityStates,
261 boost::optional<storm::dd::Add<DdType, ValueType>> const& startValues) {
262 // Create the matrix and the vector for the equation system.
263 storm::dd::Add<DdType, ValueType> maybeStatesAdd = maybeStates.template toAdd<ValueType>();
264
265 // Start by cutting away all rows that do not belong to maybe states. Note that this leaves columns targeting
266 // non-maybe states in the matrix.
267 storm::dd::Add<DdType, ValueType> submatrix = transitionMatrix * maybeStatesAdd;
268
269 // Then compute the state reward vector to use in the computation.
270 storm::dd::Add<DdType, ValueType> subvector = rewardModel.getTotalRewardVector(maybeStatesAdd, submatrix, model.getColumnVariables());
271
272 // Since we are cutting away target and infinity states, we need to account for this by giving
273 // choices the value infinity that have some successor contained in the infinity states.
274 storm::dd::Bdd<DdType> choicesWithInfinitySuccessor =
275 (maybeStates && transitionMatrixBdd && infinityStates.swapVariables(model.getRowColumnMetaVariablePairs())).existsAbstract(model.getColumnVariables());
276 subvector = choicesWithInfinitySuccessor.ite(model.getManager().template getInfinity<ValueType>(), subvector);
277
278 // Finally cut away all columns targeting non-maybe states.
279 submatrix *= maybeStatesAdd.swapVariables(model.getRowColumnMetaVariablePairs());
280
281 // Now solve the resulting equation system.
283 std::unique_ptr<storm::solver::SymbolicMinMaxLinearEquationSolver<DdType, ValueType>> solver =
284 linearEquationSolverFactory.create(submatrix, maybeStates, model.getIllegalMask() && maybeStates, model.getRowVariables(), model.getColumnVariables(),
286
287 if (storm::solver::maximize(dir)) {
288 solver->setHasUniqueSolution(true);
289 }
290
291 // Check requirements of solver.
292 storm::solver::MinMaxLinearEquationSolverRequirements requirements = solver->getRequirements(env, dir);
293 boost::optional<storm::dd::Bdd<DdType>> initialScheduler;
294 if (requirements.hasEnabledRequirement()) {
295 if (requirements.validInitialScheduler()) {
296 STORM_LOG_DEBUG("Computing valid scheduler, because the solver requires it.");
297 initialScheduler = computeValidSchedulerHint(EquationSystemType::ExpectedRewards, model, transitionMatrix, maybeStates, targetStates);
298 requirements.clearValidInitialScheduler();
299 }
300 requirements.clearLowerBounds();
301 if (requirements.uniqueSolution()) {
302 // Check whether there are end components
303 if (storm::utility::graph::performProb0E(model, transitionMatrixBdd, maybeStates, !maybeStates && model.getReachableStates()).isZero()) {
304 requirements.clearUniqueSolution();
305 }
306 }
307 STORM_LOG_THROW(!requirements.hasEnabledCriticalRequirement(), storm::exceptions::UncheckedRequirementException,
308 "Solver requirements " + requirements.getEnabledRequirementsAsString() + " not checked.");
309 }
310 if (initialScheduler) {
311 solver->setInitialScheduler(initialScheduler.get());
312 }
313 solver->setLowerBound(storm::utility::zero<ValueType>());
314 solver->setRequirementsChecked();
315
317 solver->solveEquations(env, dir, startValues ? startValues.get() : maybeStatesAdd.getDdManager().template getAddZero<ValueType>(), subvector);
318
320 model.getReachableStates(), infinityStates.ite(model.getManager().getConstant(storm::utility::infinity<ValueType>()), result)));
321}
322
323template<storm::dd::DdType DdType, typename ValueType>
326 storm::dd::Add<DdType, ValueType> const& transitionMatrix, RewardModelType const& rewardModel, storm::dd::Bdd<DdType> const& targetStates,
327 boost::optional<storm::dd::Add<DdType, ValueType>> const& startValues) {
328 // Only compute the result if there is at least one reward model.
329 STORM_LOG_THROW(!rewardModel.empty(), storm::exceptions::InvalidPropertyException, "Missing reward model for formula. Skipping formula.");
330
331 // Determine which states have a reward of infinity by definition.
332 storm::dd::Bdd<DdType> infinityStates;
333 storm::dd::Bdd<DdType> transitionMatrixBdd = transitionMatrix.notZero();
334 if (dir == OptimizationDirection::Minimize) {
336 model, transitionMatrixBdd, model.getReachableStates(), targetStates,
337 storm::utility::graph::performProbGreater0E(model, transitionMatrixBdd, model.getReachableStates(), targetStates));
338 } else {
340 model, transitionMatrixBdd, targetStates,
341 storm::utility::graph::performProbGreater0A(model, transitionMatrixBdd, model.getReachableStates(), targetStates));
342 }
343 infinityStates = !infinityStates && model.getReachableStates();
344
345 storm::dd::Bdd<DdType> maybeStates = (!targetStates && !infinityStates) && model.getReachableStates();
346
347 STORM_LOG_INFO("Preprocessing: " << infinityStates.getNonZeroCount() << " states with reward infinity, " << targetStates.getNonZeroCount()
348 << " target states (" << maybeStates.getNonZeroCount() << " states remaining).");
349
350 // If there are maybe states, we need to solve an equation system.
351 if (!maybeStates.isZero()) {
352 return computeReachabilityRewards(env, dir, model, transitionMatrix, transitionMatrixBdd, rewardModel, maybeStates, targetStates, infinityStates,
353 startValues ? maybeStates.ite(startValues.get(), model.getManager().template getAddZero<ValueType>())
354 : model.getManager().template getAddZero<ValueType>());
355 } else {
357 model.getReachableStates(),
358 infinityStates.ite(model.getManager().getConstant(storm::utility::infinity<ValueType>()), model.getManager().template getAddZero<ValueType>())));
359 }
360}
361
362template<storm::dd::DdType DdType, typename ValueType>
365 storm::dd::Add<DdType, ValueType> const& transitionMatrix, storm::dd::Bdd<DdType> const& targetStates,
366 boost::optional<storm::dd::Add<DdType, ValueType>> const& startValues) {
367 RewardModelType rewardModel(model.getManager().getConstant(storm::utility::one<ValueType>()), boost::none, boost::none);
368 return computeReachabilityRewards(env, dir, model, transitionMatrix, rewardModel, targetStates, startValues);
369}
370
373
375} // namespace helper
376} // namespace modelchecker
377} // namespace storm
Add< LibraryType, ValueType > swapVariables(std::vector< std::pair< storm::expressions::Variable, storm::expressions::Variable > > const &metaVariablePairs) const
Swaps the given pairs of meta variables in the ADD.
Definition Add.cpp:292
Add< LibraryType, ValueType > maxAbstract(std::set< storm::expressions::Variable > const &metaVariables) const
Max-abstracts from the given meta variables.
Definition Add.cpp:198
Add< LibraryType, ValueType > sumAbstract(std::set< storm::expressions::Variable > const &metaVariables) const
Sum-abstracts from the given meta variables.
Definition Add.cpp:178
Bdd< LibraryType > notZero() const
Computes a BDD that represents the function in which all assignments with a function value unequal to...
Definition Add.cpp:431
bool isZero() const
Retrieves whether this DD represents the constant zero function.
Definition Bdd.cpp:547
virtual uint_fast64_t getNonZeroCount() const override
Retrieves the number of encodings that are mapped to a non-zero value.
Definition Bdd.cpp:513
Bdd< LibraryType > swapVariables(std::vector< std::pair< storm::expressions::Variable, storm::expressions::Variable > > const &metaVariablePairs) const
Swaps the given pairs of meta variables in the BDD.
Definition Bdd.cpp:302
Bdd< LibraryType > ite(Bdd< LibraryType > const &thenBdd, Bdd< LibraryType > const &elseBdd) const
Performs an if-then-else with the given operands, i.e.
Definition Bdd.cpp:113
DdManager< LibraryType > & getDdManager() const
Retrieves the manager that is responsible for this DD.
Definition Dd.cpp:39
static std::unique_ptr< CheckResult > computeBoundedUntilProbabilities(Environment const &env, OptimizationDirection dir, storm::models::symbolic::NondeterministicModel< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, storm::dd::Bdd< DdType > const &phiStates, storm::dd::Bdd< DdType > const &psiStates, uint_fast64_t stepBound)
static std::unique_ptr< CheckResult > computeUntilProbabilities(Environment const &env, OptimizationDirection dir, storm::models::symbolic::NondeterministicModel< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, storm::dd::Bdd< DdType > const &phiStates, storm::dd::Bdd< DdType > const &psiStates, bool qualitative, boost::optional< storm::dd::Add< DdType, ValueType > > const &startValues=boost::none)
static std::unique_ptr< CheckResult > computeGloballyProbabilities(Environment const &env, OptimizationDirection dir, storm::models::symbolic::NondeterministicModel< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, storm::dd::Bdd< DdType > const &psiStates, bool qualitative)
static std::unique_ptr< CheckResult > computeReachabilityRewards(Environment const &env, OptimizationDirection dir, storm::models::symbolic::NondeterministicModel< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, RewardModelType const &rewardModel, storm::dd::Bdd< DdType > const &targetStates, boost::optional< storm::dd::Add< DdType, ValueType > > const &startValues=boost::none)
static std::unique_ptr< CheckResult > computeInstantaneousRewards(Environment const &env, OptimizationDirection dir, storm::models::symbolic::NondeterministicModel< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, RewardModelType const &rewardModel, uint_fast64_t stepBound)
storm::models::symbolic::NondeterministicModel< DdType, ValueType >::RewardModelType RewardModelType
static std::unique_ptr< CheckResult > computeReachabilityTimes(Environment const &env, OptimizationDirection dir, storm::models::symbolic::NondeterministicModel< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, storm::dd::Bdd< DdType > const &targetStates, boost::optional< storm::dd::Add< DdType, ValueType > > const &startValues=boost::none)
static std::unique_ptr< CheckResult > computeNextProbabilities(Environment const &env, OptimizationDirection dir, storm::models::symbolic::NondeterministicModel< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, storm::dd::Bdd< DdType > const &nextStates)
static std::unique_ptr< CheckResult > computeCumulativeRewards(Environment const &env, OptimizationDirection dir, storm::models::symbolic::NondeterministicModel< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, RewardModelType const &rewardModel, uint_fast64_t stepBound)
storm::dd::DdManager< Type > & getManager() const
Retrieves the manager responsible for the DDs that represent this model.
Definition Model.cpp:92
storm::dd::Add< Type, ValueType > const & getTransitionMatrix() const
Retrieves the matrix representing the transitions of the model.
Definition Model.cpp:177
std::set< storm::expressions::Variable > const & getColumnVariables() const
Retrieves the meta variables used to encode the columns of the transition matrix and the vector indic...
Definition Model.cpp:197
std::vector< std::pair< storm::expressions::Variable, storm::expressions::Variable > > const & getRowColumnMetaVariablePairs() const
Retrieves the pairs of row and column meta variables.
Definition Model.cpp:223
std::set< storm::expressions::Variable > const & getRowVariables() const
Retrieves the meta variables used to encode the rows of the transition matrix and the vector indices.
Definition Model.cpp:192
storm::dd::Bdd< Type > const & getReachableStates() const
Retrieves the reachable states of the model.
Definition Model.cpp:102
Base class for all nondeterministic symbolic models.
virtual std::set< storm::expressions::Variable > const & getNondeterminismVariables() const override
Retrieves the meta variables used to encode the nondeterminism in the model.
storm::dd::Bdd< Type > const & getIllegalMask() const
Retrieves a BDD characterizing all illegal nondeterminism encodings in the model.
virtual std::unique_ptr< storm::solver::SymbolicMinMaxLinearEquationSolver< DdType, ValueType > > create(storm::dd::Add< DdType, ValueType > const &A, storm::dd::Bdd< DdType > const &allRows, storm::dd::Bdd< DdType > const &illegalMask, std::set< storm::expressions::Variable > const &rowMetaVariables, std::set< storm::expressions::Variable > const &columnMetaVariables, std::set< storm::expressions::Variable > const &choiceVariables, std::vector< std::pair< storm::expressions::Variable, storm::expressions::Variable > > const &rowColumnMetaVariablePairs) const override
std::string getEnabledRequirementsAsString() const
Returns a string that enumerates the enabled requirements.
#define STORM_LOG_INFO(message)
Definition logging.h:29
#define STORM_LOG_DEBUG(message)
Definition logging.h:23
#define STORM_LOG_THROW(cond, exception, message)
Definition macros.h:30
std::vector< uint_fast64_t > computeValidSchedulerHint(Environment const &env, SemanticSolutionType const &type, storm::storage::SparseMatrix< ValueType > const &transitionMatrix, storm::storage::SparseMatrix< ValueType > const &backwardTransitions, storm::storage::BitVector const &maybeStates, storm::storage::BitVector const &filterStates, storm::storage::BitVector const &targetStates, boost::optional< storm::storage::BitVector > const &selectedChoices)
bool constexpr maximize(OptimizationDirection d)
bool constexpr minimize(OptimizationDirection d)
std::pair< storm::storage::BitVector, storm::storage::BitVector > performProb01Max(storm::storage::SparseMatrix< T > const &transitionMatrix, std::vector< uint_fast64_t > const &nondeterministicChoiceIndices, storm::storage::SparseMatrix< T > const &backwardTransitions, storm::storage::BitVector const &phiStates, storm::storage::BitVector const &psiStates)
Definition graph.cpp:835
void computeSchedulerProbGreater0E(storm::storage::SparseMatrix< T > const &transitionMatrix, storm::storage::SparseMatrix< T > const &backwardTransitions, storm::storage::BitVector const &phiStates, storm::storage::BitVector const &psiStates, storm::storage::Scheduler< SchedulerValueType > &scheduler, boost::optional< storm::storage::BitVector > const &rowFilter)
Computes a scheduler for the ProbGreater0E-States such that in the induced system the given psiStates...
Definition graph.cpp:542
storm::storage::BitVector performProbGreater0A(storm::storage::SparseMatrix< T > const &transitionMatrix, std::vector< uint_fast64_t > const &nondeterministicChoiceIndices, storm::storage::SparseMatrix< T > const &backwardTransitions, storm::storage::BitVector const &phiStates, storm::storage::BitVector const &psiStates, bool useStepBound, uint_fast64_t maximalSteps, boost::optional< storm::storage::BitVector > const &choiceConstraint)
Computes the sets of states that have probability greater 0 of satisfying phi until psi under any pos...
Definition graph.cpp:857
void computeSchedulerProb1E(storm::storage::BitVector const &prob1EStates, storm::storage::SparseMatrix< T > const &transitionMatrix, storm::storage::SparseMatrix< T > const &backwardTransitions, storm::storage::BitVector const &phiStates, storm::storage::BitVector const &psiStates, storm::storage::Scheduler< SchedulerValueType > &scheduler, boost::optional< storm::storage::BitVector > const &rowFilter)
Computes a scheduler for the given prob1EStates such that in the induced system the given psiStates a...
Definition graph.cpp:624
storm::storage::BitVector performProb1A(storm::models::sparse::NondeterministicModel< T, RM > const &model, storm::storage::SparseMatrix< T > const &backwardTransitions, storm::storage::BitVector const &phiStates, storm::storage::BitVector const &psiStates)
Computes the sets of states that have probability 1 of satisfying phi until psi under all possible re...
Definition graph.cpp:997
storm::storage::BitVector performProbGreater0E(storm::storage::SparseMatrix< T > const &backwardTransitions, storm::storage::BitVector const &phiStates, storm::storage::BitVector const &psiStates, bool useStepBound, uint_fast64_t maximalSteps)
Computes the sets of states that have probability greater 0 of satisfying phi until psi under at leas...
Definition graph.cpp:689
storm::storage::BitVector performProb0E(storm::models::sparse::NondeterministicModel< T, RM > const &model, storm::storage::SparseMatrix< T > const &backwardTransitions, storm::storage::BitVector const &phiStates, storm::storage::BitVector const &psiStates)
Computes the sets of states that have probability 0 of satisfying phi until psi under at least one po...
Definition graph.cpp:976
std::pair< storm::storage::BitVector, storm::storage::BitVector > performProb01Min(storm::storage::SparseMatrix< T > const &transitionMatrix, std::vector< uint_fast64_t > const &nondeterministicChoiceIndices, storm::storage::SparseMatrix< T > const &backwardTransitions, storm::storage::BitVector const &phiStates, storm::storage::BitVector const &psiStates)
Definition graph.cpp:1079
storm::storage::BitVector performProb1E(storm::storage::SparseMatrix< T > const &transitionMatrix, std::vector< uint_fast64_t > const &nondeterministicChoiceIndices, storm::storage::SparseMatrix< T > const &backwardTransitions, storm::storage::BitVector const &phiStates, storm::storage::BitVector const &psiStates, boost::optional< storm::storage::BitVector > const &choiceConstraint)
Computes the sets of states that have probability 1 of satisfying phi until psi under at least one po...
Definition graph.cpp:757
LabParser.cpp.
Definition cli.cpp:18