Storm
A Modern Probabilistic Model Checker
Loading...
Searching...
No Matches
BinaryPomdpTransformer.cpp
Go to the documentation of this file.
2#include <queue>
3
8
9namespace storm {
10namespace transformer {
11
12template<typename ValueType>
14 // Intentionally left empty
15}
16
17template<typename ValueType>
19 bool keepStateValuations) const {
20 auto data = transformTransitions(pomdp, transformSimple);
22 components.stateLabeling = transformStateLabeling(pomdp, data);
23 for (auto const& rewModel : pomdp.getRewardModels()) {
24 components.rewardModels.emplace(rewModel.first, transformRewardModel(pomdp, rewModel.second, data));
25 }
26 components.transitionMatrix = std::move(data.simpleMatrix);
27 components.observabilityClasses = std::move(data.simpleObservations);
28 if (keepStateValuations && pomdp.hasStateValuations()) {
29 components.stateValuations = pomdp.getStateValuations().blowup(data.simpleStateToOriginalState);
30 }
32 result.transformedPomdp = std::make_shared<storm::models::sparse::Pomdp<ValueType>>(std::move(components), true);
33 result.transformedStateToOriginalStateMap = data.simpleStateToOriginalState;
34 return result;
35}
36
42
43 uint64_t origState;
44 uint64_t firstRow;
45 uint64_t endRow;
48
49 uint64_t size() const {
50 return endRow - firstRow;
51 }
52
53 std::vector<BinaryPomdpTransformerRowGroup> split() const {
54 assert(size() > 1);
55 uint64_t midRow = firstRow + size() / 2;
56 std::vector<BinaryPomdpTransformerRowGroup> res;
57 res.emplace_back(origState, firstRow, midRow, origStateObservation);
59 newAuxStateId.resize(auxStateId.size() + 1, false);
60 res.back().auxStateId = newAuxStateId;
61 res.emplace_back(origState, midRow, endRow, origStateObservation);
62 newAuxStateId.set(auxStateId.size(), true);
63 res.back().auxStateId = newAuxStateId;
64 return res;
65 }
66};
67
71 return lhs.auxStateId < rhs.auxStateId;
72 } else {
74 }
75 }
76};
77
78template<typename ValueType>
79typename BinaryPomdpTransformer<ValueType>::TransformationData BinaryPomdpTransformer<ValueType>::transformTransitions(
80 storm::models::sparse::Pomdp<ValueType> const& pomdp, bool transformSimple) const {
81 auto const& matrix = pomdp.getTransitionMatrix();
82
83 // Initialize a FIFO Queue that stores the start and the end of each row group
84 std::queue<BinaryPomdpTransformerRowGroup> queue;
85 for (uint64_t state = 0; state < matrix.getRowGroupCount(); ++state) {
86 queue.emplace(state, matrix.getRowGroupIndices()[state], matrix.getRowGroupIndices()[state + 1], pomdp.getObservation(state));
87 }
88
89 std::vector<uint32_t> newObservations;
90 std::map<BinaryPomdpTransformerRowGroup, uint32_t, BinaryPomdpTransformerRowGroupCompare> observationMap;
91 storm::storage::SparseMatrixBuilder<ValueType> builder(0, 0, 0, true, true);
92 uint64_t currRow = 0;
93 std::vector<uint64_t> origRowToSimpleRowMap(pomdp.getNumberOfChoices(), std::numeric_limits<uint64_t>::max());
94 uint64_t currAuxState = queue.size();
95 std::vector<uint64_t> origStates;
96
97 while (!queue.empty()) {
98 auto group = std::move(queue.front());
99 queue.pop();
100
101 // Get the observation
102 uint64_t newObservation = observationMap.insert(std::make_pair(group, observationMap.size())).first->second;
103 newObservations.push_back(newObservation);
104
105 // Add matrix entries
106 builder.newRowGroup(currRow);
107 if (group.size() == 1) {
108 // Insert the row directly
109 for (auto const& entry : matrix.getRow(group.firstRow)) {
110 builder.addNextValue(currRow, entry.getColumn(), entry.getValue());
111 }
112 origRowToSimpleRowMap[group.firstRow] = currRow;
113 ++currRow;
114 } else if (group.size() > 1) {
115 // Split the row group into two equally large parts
116 for (auto& splittedGroup : group.split()) {
117 // Check whether we can insert the row now or whether an auxiliary state is needed
118 if (splittedGroup.size() == 1 && (!transformSimple || matrix.getRow(splittedGroup.firstRow).getNumberOfEntries() < 2)) {
119 for (auto const& entry : matrix.getRow(splittedGroup.firstRow)) {
120 builder.addNextValue(currRow, entry.getColumn(), entry.getValue());
121 }
122 origRowToSimpleRowMap[splittedGroup.firstRow] = currRow;
123 ++currRow;
124 } else {
125 queue.push(std::move(splittedGroup));
126 builder.addNextValue(currRow, currAuxState, storm::utility::one<ValueType>());
127 ++currAuxState;
128 ++currRow;
129 }
130 }
131 }
132 // Nothing to be done if group has size zero
133 origStates.push_back(group.origState);
134 }
135
136 TransformationData result;
137 result.simpleMatrix = builder.build(currRow, currAuxState, currAuxState);
138 result.simpleObservations = std::move(newObservations);
139 result.originalToSimpleChoiceMap = std::move(origRowToSimpleRowMap);
140 result.simpleStateToOriginalState = std::move(origStates);
141 return result;
142}
143
144template<typename ValueType>
145storm::models::sparse::StateLabeling BinaryPomdpTransformer<ValueType>::transformStateLabeling(storm::models::sparse::Pomdp<ValueType> const& pomdp,
146 TransformationData const& data) const {
147 storm::models::sparse::StateLabeling labeling(data.simpleMatrix.getRowGroupCount());
148 for (auto const& labelName : pomdp.getStateLabeling().getLabels()) {
149 storm::storage::BitVector newStates = pomdp.getStateLabeling().getStates(labelName);
150 newStates.resize(data.simpleMatrix.getRowGroupCount(), false);
151 if (labelName != "init") {
152 for (uint64_t newState = pomdp.getNumberOfStates(); newState < data.simpleMatrix.getRowGroupCount(); ++newState) {
153 newStates.set(newState, newStates[data.simpleStateToOriginalState[newState]]);
154 }
155 }
156 labeling.addLabel(labelName, std::move(newStates));
157 }
158 return labeling;
159}
160
161template<typename ValueType>
162storm::models::sparse::StandardRewardModel<ValueType> BinaryPomdpTransformer<ValueType>::transformRewardModel(
164 TransformationData const& data) const {
165 std::optional<std::vector<ValueType>> stateRewards, actionRewards;
166 if (rewardModel.hasStateRewards()) {
167 stateRewards = rewardModel.getStateRewardVector();
168 stateRewards.value().resize(data.simpleMatrix.getRowGroupCount(), storm::utility::zero<ValueType>());
169 }
170 if (rewardModel.hasStateActionRewards()) {
171 actionRewards = std::vector<ValueType>(data.simpleMatrix.getRowCount(), storm::utility::zero<ValueType>());
172 for (uint64_t pomdpChoice = 0; pomdpChoice < pomdp.getNumberOfChoices(); ++pomdpChoice) {
173 STORM_LOG_ASSERT(data.originalToSimpleChoiceMap[pomdpChoice] < data.simpleMatrix.getRowCount(), "Invalid entry in map for choice " << pomdpChoice);
174 actionRewards.value()[data.originalToSimpleChoiceMap[pomdpChoice]] = rewardModel.getStateActionReward(pomdpChoice);
175 }
176 }
177 STORM_LOG_THROW(!rewardModel.hasTransitionRewards(), storm::exceptions::NotSupportedException, "Transition rewards are currently not supported.");
178 return storm::models::sparse::StandardRewardModel<ValueType>(std::move(stateRewards), std::move(actionRewards));
179}
180
181template class BinaryPomdpTransformer<storm::RationalNumber>;
182
183template class BinaryPomdpTransformer<double>;
184template class BinaryPomdpTransformer<storm::RationalFunction>;
185} // namespace transformer
186} // namespace storm
storm::storage::SparseMatrix< ValueType > const & getTransitionMatrix() const
Retrieves the matrix representing the transitions of the model.
Definition Model.cpp:197
std::unordered_map< std::string, RewardModelType > const & getRewardModels() const
Retrieves the reward models.
Definition Model.cpp:699
bool hasStateValuations() const
Retrieves whether this model was build with state valuations.
Definition Model.cpp:349
storm::storage::sparse::StateValuations const & getStateValuations() const
Retrieves the valuations of the states of the model.
Definition Model.cpp:354
storm::models::sparse::StateLabeling const & getStateLabeling() const
Returns the state labeling associated with this model.
Definition Model.cpp:319
virtual uint_fast64_t getNumberOfStates() const override
Returns the number of states of the model.
Definition Model.cpp:162
uint_fast64_t getNumberOfChoices(uint_fast64_t state) const
This class represents a partially observable Markov decision process.
Definition Pomdp.h:15
uint32_t getObservation(uint64_t state) const
Definition Pomdp.cpp:63
bool hasTransitionRewards() const
Retrieves whether the reward model has transition rewards.
std::vector< ValueType > const & getStateRewardVector() const
Retrieves the state rewards of the reward model.
ValueType const & getStateActionReward(uint_fast64_t choiceIndex) const
Retrieves the state-action reward for the given choice.
bool hasStateRewards() const
Retrieves whether the reward model has state rewards.
bool hasStateActionRewards() const
Retrieves whether the reward model has state-action rewards.
This class manages the labeling of the state space with a number of (atomic) labels.
storm::storage::BitVector const & getStates(std::string const &label) const
Returns the labeling of states associated with the given label.
A bit vector that is internally represented as a vector of 64-bit values.
Definition BitVector.h:18
void resize(uint_fast64_t newLength, bool init=false)
Resizes the bit vector to hold the given new number of bits.
void set(uint_fast64_t index, bool value=true)
Sets the given truth value at the given index.
size_t size() const
Retrieves the number of bits this bit vector can store.
A class that can be used to build a sparse matrix by adding value by value.
StateValuations blowup(std::vector< uint64_t > const &mapNewToOld) const
PomdpTransformationResult< ValueType > transform(storm::models::sparse::Pomdp< ValueType > const &pomdp, bool transformSimple, bool keepStateValuations=false) const
#define STORM_LOG_ASSERT(cond, message)
Definition macros.h:11
#define STORM_LOG_THROW(cond, exception, message)
Definition macros.h:30
RewardModelType transformRewardModel(RewardModelType const &originalRewardModel, storm::storage::BitVector const &subsystem, storm::storage::BitVector const &subsystemActions, bool makeRowGroupingTrivial)
LabParser.cpp.
Definition cli.cpp:18
std::optional< storm::storage::sparse::StateValuations > stateValuations
std::unordered_map< std::string, RewardModelType > rewardModels
storm::storage::SparseMatrix< ValueType > transitionMatrix
storm::models::sparse::StateLabeling stateLabeling
std::optional< std::vector< uint32_t > > observabilityClasses
bool operator()(BinaryPomdpTransformerRowGroup const &lhs, BinaryPomdpTransformerRowGroup const &rhs) const
BinaryPomdpTransformerRowGroup(uint64_t origState, uint64_t firstRow, uint64_t endRow, uint32_t origStateObservation)
std::vector< BinaryPomdpTransformerRowGroup > split() const
std::shared_ptr< storm::models::sparse::Pomdp< ValueType > > transformedPomdp