Storm
A Modern Probabilistic Model Checker
Loading...
Searching...
No Matches
StandardPcaaWeightVectorChecker.h
Go to the documentation of this file.
1#pragma once
2
16
17namespace storm {
18namespace modelchecker {
19namespace multiobjective {
20
27template<class SparseModelType>
29 public:
30 typedef typename SparseModelType::ValueType ValueType;
32 typename std::conditional<std::is_same<SparseModelType, storm::models::sparse::MarkovAutomaton<ValueType>>::value,
35
36 /*
37 * Creates a weight vextor checker.
38 *
39 * @param model The (preprocessed) model
40 * @param objectives The (preprocessed) objectives
41 * @param possibleECActions Overapproximation of the actions that are part of an EC
42 * @param possibleBottomStates The states for which it is posible to not collect further reward with prob. 1
43 *
44 */
46
52 virtual void check(Environment const& env, std::vector<ValueType> const& weightVector) override;
53
60 virtual std::vector<ValueType> getUnderApproximationOfInitialStateResults() const override;
61 virtual std::vector<ValueType> getOverApproximationOfInitialStateResults() const override;
62
69
70 protected:
72 virtual void initializeModelTypeSpecificData(SparseModelType const& model) = 0;
74 storm::storage::SparseMatrix<ValueType> const& transitions) const = 0;
76
77 void infiniteHorizonWeightedPhase(Environment const& env, std::vector<ValueType> const& weightedActionRewardVector,
78 boost::optional<std::vector<ValueType>> const& weightedStateRewardVector);
79
85 void unboundedWeightedPhase(Environment const& env, std::vector<ValueType> const& weightedRewardVector, std::vector<ValueType> const& weightVector);
86
91 void unboundedIndividualPhase(Environment const& env, std::vector<ValueType> const& weightVector);
92
102 virtual void boundedPhase(Environment const& env, std::vector<ValueType> const& weightVector, std::vector<ValueType>& weightedRewardVector) = 0;
103
104 void updateEcQuotient(std::vector<ValueType> const& weightedRewardVector);
105
106 void setBoundsToSolver(storm::solver::AbstractEquationSolver<ValueType>& solver, bool requiresLower, bool requiresUpper, uint64_t objIndex,
107 storm::storage::SparseMatrix<ValueType> const& transitions, storm::storage::BitVector const& rowsWithSumLessOne,
108 std::vector<ValueType> const& rewards) const;
109 void setBoundsToSolver(storm::solver::AbstractEquationSolver<ValueType>& solver, bool requiresLower, bool requiresUpper,
110 std::vector<ValueType> const& weightVector, storm::storage::BitVector const& objectiveFilter,
111 storm::storage::SparseMatrix<ValueType> const& transitions, storm::storage::BitVector const& rowsWithSumLessOne,
112 std::vector<ValueType> const& rewards) const;
113 void computeAndSetBoundsToSolver(storm::solver::AbstractEquationSolver<ValueType>& solver, bool requiresLower, bool requiresUpper,
114 storm::storage::SparseMatrix<ValueType> const& transitions, storm::storage::BitVector const& rowsWithSumLessOne,
115 std::vector<ValueType> const& rewards) const;
116
120 void transformEcqSolutionToOriginalModel(std::vector<ValueType> const& ecqSolution, std::vector<uint_fast64_t> const& ecqOptimalChoices,
121 std::map<uint64_t, uint64_t> const& ecqStateToOptimalMecMap, std::vector<ValueType>& originalSolution,
122 std::vector<uint_fast64_t>& originalOptimalChoices) const;
123
124 // Data regarding the given model
125 // The transition matrix of the considered model
127 // The initial state of the considered model
128 uint64_t initialState;
129 // Overapproximation of the set of choices that are part of an end component.
131 // The actions that have reward assigned for at least one objective without upper timeBound
133 // The states for which there is a scheduler yielding reward 0 for each total reward objective
135 // stores the state action rewards for each objective.
136 std::vector<std::vector<ValueType>> actionRewards;
137 // stores the state rewards for each objective.
138 // These are only relevant for LRA objectives for MAs (otherwise, they appear within the action rewards). For other objectives/models, the corresponding
139 // vector will be empty.
140 std::vector<std::vector<ValueType>> stateRewards;
141
142 // stores the indices of the objectives for which we need to compute the long run average values
144 // stores the indices of the objectives for which there is no upper time bound
146
147 // Memory for the solution of the most recent call of check(..)
148 // becomes true after the first call of check(..)
150 // The result for the weighted reward vector (for all states of the model)
151 std::vector<ValueType> weightedResult;
152 // The results for the individual objectives (w.r.t. all states of the model)
153 std::vector<std::vector<ValueType>> objectiveResults;
154 // Stores for each objective the distance between the computed result (w.r.t. the initial state) and an over/under approximation for the actual result.
155 // The distances are stored as a (possibly negative) offset that has to be added (+) to to the objectiveResults.
156 std::vector<ValueType> offsetsToUnderApproximation;
157 std::vector<ValueType> offsetsToOverApproximation;
158 // The scheduler choices that optimize the weighted rewards of undounded objectives.
159 std::vector<uint64_t> optimalChoices;
160
161 struct EcQuotient {
163 std::vector<uint_fast64_t> ecqToOriginalChoiceMapping;
164 std::vector<uint_fast64_t> originalToEcqStateMapping;
165 std::vector<storm::storage::FlatSetStateContainer> ecqToOriginalStateMapping;
167 storm::storage::BitVector origReward0Choices; // includes total and LRA rewards
168 storm::storage::BitVector origTotalReward0Choices; // considers just total rewards
170
171 std::vector<ValueType> auxStateValues;
172 std::vector<ValueType> auxChoiceValues;
173 };
174 boost::optional<EcQuotient> ecQuotient;
175
180 boost::optional<LraMecDecomposition> lraMecDecomposition;
181};
182
183} // namespace multiobjective
184} // namespace modelchecker
185} // namespace storm
Helper class for model checking queries that depend on the long run behavior of the (nondeterministic...
Helper class for model checking queries that depend on the long run behavior of the (nondeterministic...
Helper Class that takes preprocessed Pcaa data and a weight vector and ...
void unboundedWeightedPhase(Environment const &env, std::vector< ValueType > const &weightedRewardVector, std::vector< ValueType > const &weightVector)
Determines the scheduler that optimizes the weighted reward vector of the unbounded objectives.
virtual void boundedPhase(Environment const &env, std::vector< ValueType > const &weightVector, std::vector< ValueType > &weightedRewardVector)=0
For each time epoch (starting with the maximal stepBound occurring in the objectives),...
virtual storm::modelchecker::helper::SparseNondeterministicInfiniteHorizonHelper< ValueType > createNondetInfiniteHorizonHelper(storm::storage::SparseMatrix< ValueType > const &transitions) const =0
virtual std::vector< ValueType > getOverApproximationOfInitialStateResults() const override
void computeAndSetBoundsToSolver(storm::solver::AbstractEquationSolver< ValueType > &solver, bool requiresLower, bool requiresUpper, storm::storage::SparseMatrix< ValueType > const &transitions, storm::storage::BitVector const &rowsWithSumLessOne, std::vector< ValueType > const &rewards) const
virtual DeterministicInfiniteHorizonHelperType createDetInfiniteHorizonHelper(storm::storage::SparseMatrix< ValueType > const &transitions) const =0
virtual std::vector< ValueType > getUnderApproximationOfInitialStateResults() const override
Retrieves the results of the individual objectives at the initial state of the given model.
void transformEcqSolutionToOriginalModel(std::vector< ValueType > const &ecqSolution, std::vector< uint_fast64_t > const &ecqOptimalChoices, std::map< uint64_t, uint64_t > const &ecqStateToOptimalMecMap, std::vector< ValueType > &originalSolution, std::vector< uint_fast64_t > &originalOptimalChoices) const
Transforms the results of a min-max-solver that considers a reduced model (without end components) to...
void infiniteHorizonWeightedPhase(Environment const &env, std::vector< ValueType > const &weightedActionRewardVector, boost::optional< std::vector< ValueType > > const &weightedStateRewardVector)
typename std::conditional< std::is_same< SparseModelType, storm::models::sparse::MarkovAutomaton< ValueType > >::value, storm::modelchecker::helper::SparseNondeterministicInfiniteHorizonHelper< ValueType >, storm::modelchecker::helper::SparseDeterministicInfiniteHorizonHelper< ValueType > >::type DeterministicInfiniteHorizonHelperType
virtual void initializeModelTypeSpecificData(SparseModelType const &model)=0
virtual storm::storage::Scheduler< ValueType > computeScheduler() const override
Retrieves a scheduler that induces the current values Note that check(..) has to be called before ret...
virtual void check(Environment const &env, std::vector< ValueType > const &weightVector) override
void updateEcQuotient(std::vector< ValueType > const &weightedRewardVector)
void initialize(preprocessing::SparseMultiObjectivePreprocessorResult< SparseModelType > const &preprocessorResult)
void unboundedIndividualPhase(Environment const &env, std::vector< ValueType > const &weightVector)
Computes the values of the objectives that do not have a stepBound w.r.t.
void setBoundsToSolver(storm::solver::AbstractEquationSolver< ValueType > &solver, bool requiresLower, bool requiresUpper, uint64_t objIndex, storm::storage::SparseMatrix< ValueType > const &transitions, storm::storage::BitVector const &rowsWithSumLessOne, std::vector< ValueType > const &rewards) const
A bit vector that is internally represented as a vector of 64-bit values.
Definition BitVector.h:18
This class represents the decomposition of a nondeterministic model into its maximal end components.
This class defines which action is chosen in a particular state of a non-deterministic model.
Definition Scheduler.h:18
A class that holds a possibly non-square matrix in the compressed row storage format.
LabParser.cpp.
Definition cli.cpp:18