Browse Source

Beautified the code a bit.

Former-commit-id: d4b4a738c1
main
dehnert 12 years ago
parent
commit
973e51bacb
  1. 7
      src/modelchecker/prctl/SparseMdpPrctlModelChecker.h
  2. 18
      src/utility/graph.h

7
src/modelchecker/prctl/SparseMdpPrctlModelChecker.h

@ -713,11 +713,11 @@ namespace storm {
std::vector<uint_fast64_t> scheduler(maybeStates.getNumberOfSetBits());
// For each of the states we need to resolve the nondeterministic choice based on the information we are given.
Type maxProbability = 0;
Type maxProbability = -storm::utility::constGetInfinity<Type>();
Type currentProbability = 0;
uint_fast64_t currentStateIndex = 0;
for (auto state : maybeStates) {
maxProbability = 0;
maxProbability = -storm::utility::constGetInfinity<Type>();
for (uint_fast64_t row = 0, rowEnd = this->getModel().getNondeterministicChoiceIndices()[state + 1] - this->getModel().getNondeterministicChoiceIndices()[state]; row < rowEnd; ++row) {
typename storm::storage::SparseMatrix<Type>::Rows currentRow = this->getModel().getTransitionMatrix().getRow(this->getModel().getNondeterministicChoiceIndices()[state] + row);
@ -725,9 +725,10 @@ namespace storm {
for (auto& transition : currentRow) {
currentProbability += transition.value() * distancesToTarget[transition.column()];
// currentProbability -= transition.value() * (1 - distancesToNonTarget[transition.column()]);
}
if (currentProbability < maxProbability) {
if (currentProbability > maxProbability) {
maxProbability = currentProbability;
scheduler[currentStateIndex] = row;
}

18
src/utility/graph.h

@ -22,10 +22,8 @@
extern log4cplus::Logger logger;
namespace storm {
namespace utility {
namespace graph {
namespace utility {
namespace graph {
/*!
* Performs a backwards breadt-first search trough the underlying graph structure
@ -767,7 +765,7 @@ namespace graph {
template<typename T>
struct DistanceCompare {
bool operator()(std::pair<T, uint_fast64_t> const& lhs, std::pair<T, uint_fast64_t> const& rhs) const {
return lhs.first < rhs.first || (lhs.first == rhs.first && lhs.second < rhs.second);
return lhs.first > rhs.first || (lhs.first == rhs.first && lhs.second > rhs.second);
}
};
@ -795,8 +793,8 @@ namespace graph {
// Set the probability to 1 for all starting states.
std::set<std::pair<T, uint_fast64_t>, DistanceCompare<T>> probabilityStateSet;
for (auto state : startingStates) {
probabilityStateSet.emplace(storm::utility::constGetZero<T>(), state);
probabilities[state] = 1;
probabilityStateSet.emplace(storm::utility::constGetOne<T>(), state);
probabilities[state] = storm::utility::constGetOne<T>();
}
// As long as there is one reachable state, we need to consider it.
@ -837,10 +835,8 @@ namespace graph {
return result;
}
} // namespace graph
} // namespace utility
} // namespace graph
} // namespace utility
} // namespace storm
#endif /* STORM_UTILITY_GRAPH_H_ */
Loading…
Cancel
Save