Storm 1.11.1.1
A Modern Probabilistic Model Checker
Loading...
Searching...
No Matches
SymbolicDtmcPrctlHelper.cpp
Go to the documentation of this file.
2
11#include "storm/utility/graph.h"
12
13namespace storm {
14namespace modelchecker {
15namespace helper {
16
17template<storm::dd::DdType DdType, typename ValueType>
20 storm::dd::Bdd<DdType> const& phiStates, storm::dd::Bdd<DdType> const& psiStates, bool qualitative,
21 boost::optional<storm::dd::Add<DdType, ValueType>> const& startValues) {
22 // We need to identify the states which have to be taken out of the matrix, i.e. all states that have
23 // probability 0 and 1 of satisfying the until-formula.
24 STORM_LOG_TRACE("Found " << phiStates.getNonZeroCount() << " phi states and " << psiStates.getNonZeroCount() << " psi states.");
25 std::pair<storm::dd::Bdd<DdType>, storm::dd::Bdd<DdType>> statesWithProbability01 =
26 storm::utility::graph::performProb01(model, transitionMatrix.notZero(), phiStates, psiStates);
27 storm::dd::Bdd<DdType> maybeStates = !statesWithProbability01.first && !statesWithProbability01.second && model.getReachableStates();
28
29 STORM_LOG_INFO("Preprocessing: " << statesWithProbability01.first.getNonZeroCount() << " states with probability 0, "
30 << statesWithProbability01.second.getNonZeroCount() << " with probability 1 (" << maybeStates.getNonZeroCount()
31 << " states remaining).");
32
33 // Check whether we need to compute exact probabilities for some states.
34 if (qualitative) {
35 // Set the values for all maybe-states to 0.5 to indicate that their probability values are neither 0 nor 1.
36 return statesWithProbability01.second.template toAdd<ValueType>() +
37 maybeStates.template toAdd<ValueType>() * model.getManager().getConstant(storm::utility::convertNumber<ValueType>(0.5));
38 } else {
39 // If there are maybe states, we need to solve an equation system.
40 if (!maybeStates.isZero()) {
41 return computeUntilProbabilities(env, model, transitionMatrix, maybeStates, statesWithProbability01.second,
42 startValues ? maybeStates.ite(startValues.get(), model.getManager().template getAddZero<ValueType>())
43 : model.getManager().template getAddZero<ValueType>());
44 } else {
45 return statesWithProbability01.second.template toAdd<ValueType>();
46 }
47 }
48}
49
50template<storm::dd::DdType DdType, typename ValueType>
53 storm::dd::Bdd<DdType> const& maybeStates, storm::dd::Bdd<DdType> const& statesWithProbability1,
54 boost::optional<storm::dd::Add<DdType, ValueType>> const& startValues) {
55 // Create the matrix and the vector for the equation system.
56 storm::dd::Add<DdType, ValueType> maybeStatesAdd = maybeStates.template toAdd<ValueType>();
57
58 // Start by cutting away all rows that do not belong to maybe states. Note that this leaves columns targeting
59 // non-maybe states in the matrix.
60 storm::dd::Add<DdType, ValueType> submatrix = transitionMatrix * maybeStatesAdd;
61
62 // Then compute the vector that contains the one-step probabilities to a state with probability 1 for all
63 // maybe states.
64 storm::dd::Add<DdType, ValueType> prob1StatesAsColumn =
65 statesWithProbability1.template toAdd<ValueType>().swapVariables(model.getRowColumnMetaVariablePairs());
66 storm::dd::Add<DdType, ValueType> subvector = submatrix * prob1StatesAsColumn;
67 subvector = subvector.sumAbstract(model.getColumnVariables());
68
69 // Finally cut away all columns targeting non-maybe states and convert the matrix into the matrix needed
70 // for solving the equation system (i.e. compute (I-A)).
71 submatrix *= maybeStatesAdd.swapVariables(model.getRowColumnMetaVariablePairs());
74 submatrix = (model.getRowColumnIdentity() * maybeStatesAdd) - submatrix;
75 }
76
77 // Solve the equation system.
78 std::unique_ptr<storm::solver::SymbolicLinearEquationSolver<DdType, ValueType>> solver = linearEquationSolverFactory.create(
79 env, submatrix, maybeStates, model.getRowVariables(), model.getColumnVariables(), model.getRowColumnMetaVariablePairs());
80 solver->setBounds(storm::utility::zero<ValueType>(), storm::utility::one<ValueType>());
81 storm::dd::Add<DdType, ValueType> result = solver->solveEquations(env, model.getManager().template getAddZero<ValueType>(), subvector);
82
83 return statesWithProbability1.template toAdd<ValueType>() + result;
84}
85
86template<storm::dd::DdType DdType, typename ValueType>
89 storm::dd::Bdd<DdType> const& psiStates, bool qualitative) {
91 computeUntilProbabilities(env, model, transitionMatrix, model.getReachableStates(), !psiStates && model.getReachableStates(), qualitative);
92 result = result.getDdManager().template getAddOne<ValueType>() - result;
93 return result;
94}
95
96template<storm::dd::DdType DdType, typename ValueType>
103
104template<storm::dd::DdType DdType, typename ValueType>
107 storm::dd::Bdd<DdType> const& phiStates, storm::dd::Bdd<DdType> const& psiStates, uint_fast64_t stepBound) {
108 // We need to identify the states which have to be taken out of the matrix, i.e. all states that have
109 // probability 0 or 1 of satisfying the until-formula.
110 storm::dd::Bdd<DdType> statesWithProbabilityGreater0 =
111 storm::utility::graph::performProbGreater0(model, transitionMatrix.notZero(), phiStates, psiStates, stepBound);
112 storm::dd::Bdd<DdType> maybeStates = statesWithProbabilityGreater0 && !psiStates && model.getReachableStates();
113
114 // If there are maybe states, we need to perform matrix-vector multiplications.
115 if (!maybeStates.isZero()) {
116 // Create the matrix and the vector for the equation system.
117 storm::dd::Add<DdType, ValueType> maybeStatesAdd = maybeStates.template toAdd<ValueType>();
118
119 // Start by cutting away all rows that do not belong to maybe states. Note that this leaves columns targeting
120 // non-maybe states in the matrix.
121 storm::dd::Add<DdType, ValueType> submatrix = transitionMatrix * maybeStatesAdd;
122
123 // Then compute the vector that contains the one-step probabilities to a state with probability 1 for all
124 // maybe states.
125 storm::dd::Add<DdType, ValueType> prob1StatesAsColumn = psiStates.template toAdd<ValueType>().swapVariables(model.getRowColumnMetaVariablePairs());
126 storm::dd::Add<DdType, ValueType> subvector = (submatrix * prob1StatesAsColumn).sumAbstract(model.getColumnVariables());
127
128 // Finally cut away all columns targeting non-maybe states.
129 submatrix *= maybeStatesAdd.swapVariables(model.getRowColumnMetaVariablePairs());
130
131 // Perform the matrix-vector multiplication.
133 std::unique_ptr<storm::solver::SymbolicLinearEquationSolver<DdType, ValueType>> solver = linearEquationSolverFactory.create(
134 env, submatrix, maybeStates, model.getRowVariables(), model.getColumnVariables(), model.getRowColumnMetaVariablePairs());
135 storm::dd::Add<DdType, ValueType> result = solver->multiply(model.getManager().template getAddZero<ValueType>(), &subvector, stepBound);
136
137 return psiStates.template toAdd<ValueType>() + result;
138 } else {
139 return psiStates.template toAdd<ValueType>();
140 }
141}
142
143template<storm::dd::DdType DdType, typename ValueType>
146 RewardModelType const& rewardModel, uint_fast64_t stepBound) {
147 // Only compute the result if the model has at least one reward this->getModel().
148 STORM_LOG_THROW(!rewardModel.empty(), storm::exceptions::InvalidPropertyException, "Missing reward model for formula. Skipping formula.");
149
150 // Compute the reward vector to add in each step based on the available reward models.
151 storm::dd::Add<DdType, ValueType> totalRewardVector = rewardModel.getTotalRewardVector(transitionMatrix, model.getColumnVariables());
152
153 // Perform the matrix-vector multiplication.
155 std::unique_ptr<storm::solver::SymbolicLinearEquationSolver<DdType, ValueType>> solver = linearEquationSolverFactory.create(
156 env, transitionMatrix, model.getReachableStates(), model.getRowVariables(), model.getColumnVariables(), model.getRowColumnMetaVariablePairs());
157 return solver->multiply(model.getManager().template getAddZero<ValueType>(), &totalRewardVector, stepBound);
158}
159
160template<storm::dd::DdType DdType, typename ValueType>
163 RewardModelType const& rewardModel, uint_fast64_t stepBound) {
164 // Only compute the result if the model has at least one reward this->getModel().
165 STORM_LOG_THROW(rewardModel.hasStateRewards(), storm::exceptions::InvalidPropertyException,
166 "Computing instantaneous rewards for a reward model that does not define any state-rewards. The result is trivially 0.");
167
168 // Perform the matrix-vector multiplication.
170 std::unique_ptr<storm::solver::SymbolicLinearEquationSolver<DdType, ValueType>> solver = linearEquationSolverFactory.create(
171 env, transitionMatrix, model.getReachableStates(), model.getRowVariables(), model.getColumnVariables(), model.getRowColumnMetaVariablePairs());
172 return solver->multiply(rewardModel.getStateRewardVector(), nullptr, stepBound);
173}
174
175template<storm::dd::DdType DdType, typename ValueType>
178 RewardModelType const& rewardModel, storm::dd::Bdd<DdType> const& targetStates, bool qualitative,
179 boost::optional<storm::dd::Add<DdType, ValueType>> const& startValues) {
180 // Only compute the result if there is at least one reward model.
181 STORM_LOG_THROW(!rewardModel.empty(), storm::exceptions::InvalidPropertyException, "Missing reward model for formula. Skipping formula.");
182
183 // Determine which states have a reward of infinity by definition.
184 storm::dd::Bdd<DdType> infinityStates = storm::utility::graph::performProb1(model, transitionMatrix.notZero(), model.getReachableStates(), targetStates);
185 infinityStates = !infinityStates && model.getReachableStates();
186 storm::dd::Bdd<DdType> maybeStates = (!targetStates && !infinityStates) && model.getReachableStates();
187
188 STORM_LOG_INFO("Preprocessing: " << infinityStates.getNonZeroCount() << " states with reward infinity, " << targetStates.getNonZeroCount()
189 << " target states (" << maybeStates.getNonZeroCount() << " states remaining).");
190
191 // Check whether we need to compute exact rewards for some states.
192 if (qualitative) {
193 // Set the values for all maybe-states to 1 to indicate that their reward values
194 // are neither 0 nor infinity.
195 return infinityStates.ite(model.getManager().getConstant(storm::utility::infinity<ValueType>()), model.getManager().template getAddZero<ValueType>()) +
196 maybeStates.template toAdd<ValueType>() * model.getManager().getConstant(storm::utility::one<ValueType>());
197 } else {
198 // If there are maybe states, we need to solve an equation system.
199 if (!maybeStates.isZero()) {
200 return computeReachabilityRewards(env, model, transitionMatrix, rewardModel, maybeStates, targetStates, infinityStates,
201 startValues ? maybeStates.ite(startValues.get(), model.getManager().template getAddZero<ValueType>())
202 : model.getManager().template getAddZero<ValueType>());
203 } else {
204 return infinityStates.ite(model.getManager().getConstant(storm::utility::infinity<ValueType>()),
205 model.getManager().getConstant(storm::utility::zero<ValueType>()));
206 }
207 }
208}
209
210template<storm::dd::DdType DdType, typename ValueType>
213 RewardModelType const& rewardModel, storm::dd::Bdd<DdType> const& maybeStates, storm::dd::Bdd<DdType> const& targetStates,
214 storm::dd::Bdd<DdType> const& infinityStates, boost::optional<storm::dd::Add<DdType, ValueType>> const& startValues) {
215 // Create the matrix and the vector for the equation system.
216 storm::dd::Add<DdType, ValueType> maybeStatesAdd = maybeStates.template toAdd<ValueType>();
217
218 // Start by cutting away all rows that do not belong to maybe states. Note that this leaves columns targeting
219 // non-maybe states in the matrix.
220 storm::dd::Add<DdType, ValueType> submatrix = transitionMatrix * maybeStatesAdd;
221
222 // Then compute the state reward vector to use in the computation.
223 storm::dd::Add<DdType, ValueType> subvector = rewardModel.getTotalRewardVector(maybeStatesAdd, submatrix, model.getColumnVariables());
224
225 // Finally cut away all columns targeting non-maybe states and convert the matrix into the matrix needed
226 // for solving the equation system (i.e. compute (I-A)).
227 submatrix *= maybeStatesAdd.swapVariables(model.getRowColumnMetaVariablePairs());
230 submatrix = (model.getRowColumnIdentity() * maybeStatesAdd) - submatrix;
231 }
232
233 // Solve the equation system.
234 std::unique_ptr<storm::solver::SymbolicLinearEquationSolver<DdType, ValueType>> solver = linearEquationSolverFactory.create(
235 env, submatrix, maybeStates, model.getRowVariables(), model.getColumnVariables(), model.getRowColumnMetaVariablePairs());
236 solver->setLowerBound(storm::utility::zero<ValueType>());
238 solver->solveEquations(env, startValues ? startValues.get() : maybeStatesAdd.getDdManager().template getAddZero<ValueType>(), subvector);
239
240 return infinityStates.ite(model.getManager().getConstant(storm::utility::infinity<ValueType>()), result);
241}
242
243template<storm::dd::DdType DdType, typename ValueType>
246 storm::dd::Bdd<DdType> const& targetStates, bool qualitative, boost::optional<storm::dd::Add<DdType, ValueType>> const& startValues) {
247 RewardModelType rewardModel(model.getManager().getConstant(storm::utility::one<ValueType>()), boost::none, boost::none);
248 return computeReachabilityRewards(env, model, transitionMatrix, rewardModel, targetStates, qualitative, startValues);
249}
250
253
256} // namespace helper
257} // namespace modelchecker
258} // 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:285
Add< LibraryType, ValueType > sumAbstract(std::set< storm::expressions::Variable > const &metaVariables) const
Sum-abstracts from the given meta variables.
Definition Add.cpp:171
Bdd< LibraryType > notZero() const
Computes a BDD that represents the function in which all assignments with a function value unequal to...
Definition Add.cpp:424
bool isZero() const
Retrieves whether this DD represents the constant zero function.
Definition Bdd.cpp:541
virtual uint_fast64_t getNonZeroCount() const override
Retrieves the number of encodings that are mapped to a non-zero value.
Definition Bdd.cpp:507
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:296
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:107
DdManager< LibraryType > & getDdManager() const
Retrieves the manager that is responsible for this DD.
Definition Dd.cpp:38
static storm::dd::Add< DdType, ValueType > computeUntilProbabilities(Environment const &env, storm::models::symbolic::Model< 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 storm::dd::Add< DdType, ValueType > computeCumulativeRewards(Environment const &env, storm::models::symbolic::Model< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, RewardModelType const &rewardModel, uint_fast64_t stepBound)
static storm::dd::Add< DdType, ValueType > computeNextProbabilities(Environment const &env, storm::models::symbolic::Model< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, storm::dd::Bdd< DdType > const &nextStates)
static storm::dd::Add< DdType, ValueType > computeBoundedUntilProbabilities(Environment const &env, storm::models::symbolic::Model< 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 storm::dd::Add< DdType, ValueType > computeReachabilityTimes(Environment const &env, storm::models::symbolic::Model< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, storm::dd::Bdd< DdType > const &targetStates, bool qualitative, boost::optional< storm::dd::Add< DdType, ValueType > > const &startValues=boost::none)
static storm::dd::Add< DdType, ValueType > computeReachabilityRewards(Environment const &env, storm::models::symbolic::Model< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, RewardModelType const &rewardModel, storm::dd::Bdd< DdType > const &targetStates, bool qualitative, boost::optional< storm::dd::Add< DdType, ValueType > > const &startValues=boost::none)
static storm::dd::Add< DdType, ValueType > computeGloballyProbabilities(Environment const &env, storm::models::symbolic::Model< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, storm::dd::Bdd< DdType > const &psiStates, bool qualitative)
static storm::dd::Add< DdType, ValueType > computeInstantaneousRewards(Environment const &env, storm::models::symbolic::Model< DdType, ValueType > const &model, storm::dd::Add< DdType, ValueType > const &transitionMatrix, RewardModelType const &rewardModel, uint_fast64_t stepBound)
Base class for all symbolic models.
Definition Model.h:42
storm::dd::DdManager< Type > & getManager() const
Retrieves the manager responsible for the DDs that represent this model.
Definition Model.cpp:87
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:192
storm::dd::Add< Type, ValueType > getRowColumnIdentity() const
Retrieves an ADD that represents the diagonal of the transition matrix.
Definition Model.cpp:233
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:218
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:187
storm::dd::Bdd< Type > const & getReachableStates() const
Retrieves the reachable states of the model.
Definition Model.cpp:97
bool hasStateRewards() const
Retrieves whether the reward model has state rewards.
storm::dd::Add< Type, ValueType > getTotalRewardVector(storm::dd::Add< Type, ValueType > const &transitionMatrix, std::set< storm::expressions::Variable > const &columnVariables) const
Creates a vector representing the complete reward vector based on the state-, state-action- and trans...
bool empty() const
Retrieves whether the reward model is empty.
storm::dd::Add< Type, ValueType > const & getStateRewardVector() const
Retrieves the state rewards of the reward model.
virtual std::unique_ptr< storm::solver::SymbolicLinearEquationSolver< DdType, ValueType > > create(Environment const &env) const override
LinearEquationSolverProblemFormat getEquationProblemFormat(Environment const &env) const
#define STORM_LOG_INFO(message)
Definition logging.h:24
#define STORM_LOG_TRACE(message)
Definition logging.h:12
#define STORM_LOG_THROW(cond, exception, message)
Definition macros.h:30
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:393
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:315
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:376