Browse Source

Avoided conversion of memory states. They are now directly represented as 64 bit integers

tempestpy_adaptions
TimQu 7 years ago
parent
commit
36c3a4d9ef
  1. 6
      src/storm/logic/BoundedUntilFormula.cpp
  2. 74
      src/storm/modelchecker/multiobjective/rewardbounded/MemoryStateManager.cpp
  3. 40
      src/storm/modelchecker/multiobjective/rewardbounded/MemoryStateManager.h
  4. 145
      src/storm/modelchecker/multiobjective/rewardbounded/MultiDimensionalRewardUnfolding.cpp
  5. 2
      src/storm/modelchecker/multiobjective/rewardbounded/MultiDimensionalRewardUnfolding.h
  6. 236
      src/storm/modelchecker/multiobjective/rewardbounded/ProductModel.cpp
  7. 29
      src/storm/modelchecker/multiobjective/rewardbounded/ProductModel.h
  8. 10
      src/test/storm/modelchecker/SparseMdpMultiDimensionalRewardUnfoldingTest.cpp

6
src/storm/logic/BoundedUntilFormula.cpp

@ -121,7 +121,7 @@ namespace storm {
}
Formula const& BoundedUntilFormula::getRightSubformula() const {
STORM_LOG_ASSERT(leftSubformula.size() == 1, "The right subformula is not unique.");
STORM_LOG_ASSERT(rightSubformula.size() == 1, "The right subformula is not unique.");
return *rightSubformula.at(0);
}
@ -323,10 +323,10 @@ namespace storm {
if (this->isMultiDimensional()) {
out << "]";
}
this->getRightSubformula().writeToStream(out);
}
this->getRightSubformula().writeToStream(out);
return out;
}
}

74
src/storm/modelchecker/multiobjective/rewardbounded/MemoryStateManager.cpp

@ -0,0 +1,74 @@
#include "storm/modelchecker/multiobjective/rewardbounded/MemoryStateManager.h"
#include "storm/utility/macros.h"
#include "storm/exceptions/IllegalArgumentException.h"
namespace storm {
namespace modelchecker {
namespace multiobjective {
MemoryStateManager::MemoryStateManager(uint64_t dimensionCount) : dimensionCount(dimensionCount), dimensionBitMask(1ull), relevantBitsMask((1ull << dimensionCount) - 1), stateCount(dimensionBitMask << dimensionCount) {
STORM_LOG_THROW(dimensionCount > 0, storm::exceptions::IllegalArgumentException, "Invoked MemoryStateManager with zero dimension count.");
STORM_LOG_THROW(dimensionCount <= 64, storm::exceptions::IllegalArgumentException, "Invoked MemoryStateManager with too many dimensions.");
}
uint64_t const& MemoryStateManager::getDimensionCount() const {
return dimensionCount;
}
uint64_t const& MemoryStateManager::getMemoryStateCount() const {
return stateCount;
}
MemoryStateManager::MemoryState MemoryStateManager::getInitialMemoryState() const {
STORM_LOG_ASSERT(dimensionCount > 0, "Invoked MemoryStateManager with zero dimension count.");
return relevantBitsMask;
}
bool MemoryStateManager::isRelevantDimension(MemoryState const& state, uint64_t dimension) const {
STORM_LOG_ASSERT(dimensionCount > 0, "Invoked MemoryStateManager with zero dimension count.");
return (state & (dimensionBitMask << dimension)) != 0;
}
void MemoryStateManager::setRelevantDimension(MemoryState& state, uint64_t dimension, bool value) const {
STORM_LOG_ASSERT(dimensionCount > 0, "Invoked MemoryStateManager with zero dimension count.");
STORM_LOG_ASSERT(dimension < dimensionCount, "Tried to set a dimension that is larger then the number of considered dimensions");
if (value) {
state |= (dimensionBitMask << dimension);
} else {
state &= ~(dimensionBitMask << dimension);
}
assert(state < getMemoryStateCount());
}
void MemoryStateManager::setRelevantDimensions(MemoryState& state, storm::storage::BitVector const& dimensions, bool value) const {
STORM_LOG_ASSERT(dimensionCount > 0, "Invoked MemoryStateManager with zero dimension count.");
STORM_LOG_ASSERT(dimensions.size() == dimensionCount, "Invalid size of given bitset.");
if (value) {
for (auto const& d : dimensions) {
state |= (dimensionBitMask << d);
}
} else {
for (auto const& d : dimensions) {
state &= ~(dimensionBitMask << d);
}
}
assert(state < getMemoryStateCount());
}
std::string MemoryStateManager::toString(MemoryState const& state) const {
STORM_LOG_ASSERT(dimensionCount > 0, "Invoked EpochManager with zero dimension count.");
std::string res = "[";
res += (isRelevantDimension(state, 0) ? "1" : "0");
for (uint64_t d = 1; d < dimensionCount; ++d) {
res += ", ";
res += (isRelevantDimension(state, d) ? "1" : "0");
}
res += "]";
return res;
}
}
}
}

40
src/storm/modelchecker/multiobjective/rewardbounded/MemoryStateManager.h

@ -0,0 +1,40 @@
#pragma once
#include <vector>
#include <set>
#include <string>
#include "storm/storage/BitVector.h"
namespace storm {
namespace modelchecker {
namespace multiobjective {
class MemoryStateManager {
public:
typedef uint64_t MemoryState;
MemoryStateManager(uint64_t dimensionCount);
uint64_t const& getDimensionCount() const;
uint64_t const& getMemoryStateCount() const;
MemoryState getInitialMemoryState() const;
bool isRelevantDimension(MemoryState const& state, uint64_t dimension) const;
void setRelevantDimension(MemoryState& state, uint64_t dimension, bool value = true) const;
void setRelevantDimensions(MemoryState& state, storm::storage::BitVector const& dimensions, bool value = true) const;
std::string toString(MemoryState const& state) const;
private:
uint64_t const dimensionCount;
uint64_t const dimensionBitMask;
uint64_t const relevantBitsMask;
uint64_t const stateCount;
};
}
}
}

145
src/storm/modelchecker/multiobjective/rewardbounded/MultiDimensionalRewardUnfolding.cpp

@ -6,7 +6,6 @@
#include "storm/utility/macros.h"
#include "storm/logic/Formulas.h"
#include "storm/storage/memorystructure/MemoryStructureBuilder.h"
#include "storm/modelchecker/propositional/SparsePropositionalModelChecker.h"
#include "storm/modelchecker/results/ExplicitQualitativeCheckResult.h"
@ -83,7 +82,10 @@ namespace storm {
}
dimension.memoryLabel = memLabel;
dimension.isUpperBounded = subformula.hasUpperBound(dim);
// for simplicity we do not allow intervals or unbounded formulas.
STORM_LOG_THROW(subformula.hasLowerBound(dim) != dimension.isUpperBounded, storm::exceptions::NotSupportedException, "Bounded until formulas are only supported by this method if they consider either an upper bound or a lower bound. Got " << subformula << " instead.");
// lower bounded until formulas with non-trivial left hand side are excluded as this would require some additional effort (in particular the ProductModel::transformMemoryState method).
STORM_LOG_THROW(dimension.isUpperBounded || subformula.getLeftSubformula(dim).isTrueFormula(), storm::exceptions::NotSupportedException, "Lower bounded until formulas are only supported by this method if the left subformula is 'true'. Got " << subformula << " instead.");
if (subformula.getTimeBoundReference(dim).isTimeBound() || subformula.getTimeBoundReference(dim).isStepBound()) {
dimensionWiseEpochSteps.push_back(std::vector<uint64_t>(model.getNumberOfChoices(), 1));
dimension.scalingFactor = storm::utility::one<ValueType>();
@ -124,12 +126,11 @@ namespace storm {
for (uint64_t currDim = dim; currDim < dim + boundedUntilFormula.getDimension(); ++currDim ) {
if (!boundedUntilFormula.hasMultiDimensionalSubformulas() || dimensions[currDim].isUpperBounded) {
dimensions[currDim].dependentDimensions = objDimensions;
std::cout << "dimension " << currDim << " has depDims: " << dimensions[currDim].dependentDimensions << std::endl;
} else {
dimensions[currDim].dependentDimensions = storm::storage::BitVector(dimensions.size(), false);
dimensions[currDim].dependentDimensions.set(currDim, true);
std::cout << "dimension " << currDim << " has depDims: " << dimensions[currDim].dependentDimensions << std::endl;
}
// std::cout << "dimension " << currDim << " has depDims: " << dimensions[currDim].dependentDimensions << std::endl;
}
dim += boundedUntilFormula.getDimension();
} else if (objectives[objIndex].formula->isRewardOperatorFormula() && objectives[objIndex].formula->getSubformula().isCumulativeRewardFormula()) {
@ -138,7 +139,6 @@ namespace storm {
++dim;
}
std::cout << "obj " << objIndex << " has dimensions " << objDimensions << std::endl;
objectiveDimensions.push_back(std::move(objDimensions));
}
assert(dim == dimensions.size());
@ -172,14 +172,7 @@ namespace storm {
template<typename ValueType, bool SingleObjectiveMode>
void MultiDimensionalRewardUnfolding<ValueType, SingleObjectiveMode>::initializeMemoryProduct(std::vector<Epoch> const& epochSteps) {
// build the memory structure
auto memoryStructure = computeMemoryStructure();
// build a mapping between the different representations of memory states
auto memoryStateMap = computeMemoryStateMap(memoryStructure);
productModel = std::make_unique<ProductModel<ValueType>>(model, memoryStructure, dimensions, objectiveDimensions, epochManager, std::move(memoryStateMap), epochSteps);
productModel = std::make_unique<ProductModel<ValueType>>(model, objectives, dimensions, objectiveDimensions, epochManager, epochSteps);
}
template<typename ValueType, bool SingleObjectiveMode>
@ -205,6 +198,7 @@ namespace storm {
}
STORM_LOG_THROW(!bound.containsVariables(), storm::exceptions::NotSupportedException, "The bound " << bound << " contains undefined constants.");
ValueType discretizedBound = storm::utility::convertNumber<ValueType>(bound.evaluateAsRational());
STORM_LOG_THROW(dimensions[dim].isUpperBounded || isStrict || !storm::utility::isZero(discretizedBound), storm::exceptions::NotSupportedException, "Lower bounds need to be either strict or greater than zero.");
discretizedBound /= dimensions[dim].scalingFactor;
if (storm::utility::isInteger(discretizedBound)) {
if (isStrict == dimensions[dim].isUpperBounded) {
@ -296,18 +290,17 @@ namespace storm {
uint64_t productChoice = epochModelToProductChoiceMap[reducedChoice];
uint64_t productState = productModel->getProductStateFromChoice(productChoice);
auto const& memoryState = productModel->getMemoryState(productState);
auto const& memoryStateBv = productModel->convertMemoryState(memoryState);
Epoch successorEpoch = epochManager.getSuccessorEpoch(epoch, productModel->getSteps()[productChoice]);
EpochClass successorEpochClass = epochManager.getEpochClass(successorEpoch);
// Find out whether objective reward is earned for the current choice
// Objective reward is not earned if
// a) there is an upper bounded subObjective that is still relevant but the corresponding reward bound is passed after taking the choice
// a) there is an upper bounded subObjective that is __still_relevant__ but the corresponding reward bound is passed after taking the choice
// b) there is a lower bounded subObjective and the corresponding reward bound is not passed yet.
for (uint64_t objIndex = 0; objIndex < this->objectives.size(); ++objIndex) {
bool rewardEarned = !storm::utility::isZero(epochModel.objectiveRewards[objIndex][reducedChoice]);
if (rewardEarned) {
for (auto const& dim : objectiveDimensions[objIndex]) {
if (dimensions[dim].isUpperBounded == epochManager.isBottomDimension(successorEpoch, dim) && memoryStateBv.get(dim)) {
if (dimensions[dim].isUpperBounded == epochManager.isBottomDimension(successorEpoch, dim) && productModel->getMemoryStateManager().isRelevantDimension(memoryState, dim)) {
rewardEarned = false;
break;
}
@ -583,129 +576,15 @@ namespace storm {
template<typename ValueType, bool SingleObjectiveMode>
typename MultiDimensionalRewardUnfolding<ValueType, SingleObjectiveMode>::SolutionType const& MultiDimensionalRewardUnfolding<ValueType, SingleObjectiveMode>::getInitialStateResult(Epoch const& epoch) {
STORM_LOG_ASSERT(model.getInitialStates().getNumberOfSetBits() == 1, "The model has multiple initial states.");
STORM_LOG_ASSERT(productModel->getProduct().getInitialStates().getNumberOfSetBits() == 1, "The product has multiple initial states.");
return getStateSolution(epoch, *productModel->getProduct().getInitialStates().begin());
STORM_LOG_ASSERT(!epochManager.hasBottomDimension(epoch), "Tried to get the initial state result in an epoch that still contains at least one bottom dimension.");
return getStateSolution(epoch, productModel->getInitialProductState(*model.getInitialStates().begin(), model.getInitialStates()));
}
template<typename ValueType, bool SingleObjectiveMode>
typename MultiDimensionalRewardUnfolding<ValueType, SingleObjectiveMode>::SolutionType const& MultiDimensionalRewardUnfolding<ValueType, SingleObjectiveMode>::getInitialStateResult(Epoch const& epoch, uint64_t initialStateIndex) {
STORM_LOG_ASSERT(model.getInitialStates().get(initialStateIndex), "The given model state is not an initial state.");
for (uint64_t memState = 0; memState < productModel->getNumberOfMemoryState(); ++memState) {
uint64_t productState = productModel->getProductState(initialStateIndex, memState);
if (productModel->getProduct().getInitialStates().get(productState)) {
return getStateSolution(epoch, productState);
}
}
STORM_LOG_THROW(false, storm::exceptions::UnexpectedException, "Could not find the initial product state corresponding to the given initial model state.");
return getStateSolution(epoch, -1ull);
}
template<typename ValueType, bool SingleObjectiveMode>
storm::storage::MemoryStructure MultiDimensionalRewardUnfolding<ValueType, SingleObjectiveMode>::computeMemoryStructure() const {
storm::modelchecker::SparsePropositionalModelChecker<storm::models::sparse::Mdp<ValueType>> mc(model);
// Create a memory structure that remembers whether (sub)objectives are satisfied
storm::storage::MemoryStructure memory = storm::storage::MemoryStructureBuilder<ValueType>::buildTrivialMemoryStructure(model);
for (uint64_t objIndex = 0; objIndex < objectives.size(); ++objIndex) {
if (!objectives[objIndex].formula->isProbabilityOperatorFormula()) {
continue;
}
std::vector<uint64_t> dimensionIndexMap;
for (auto const& globalDimensionIndex : objectiveDimensions[objIndex]) {
dimensionIndexMap.push_back(globalDimensionIndex);
}
// collect the memory states for this objective
std::vector<storm::storage::BitVector> objMemStates;
storm::storage::BitVector m(dimensionIndexMap.size(), false);
for (; !m.full(); m.increment()) {
objMemStates.push_back(~m);
}
objMemStates.push_back(~m);
assert(objMemStates.size() == 1ull << dimensionIndexMap.size());
// build objective memory
auto objMemoryBuilder = storm::storage::MemoryStructureBuilder<ValueType>(objMemStates.size(), model);
// Get the set of states that for all subobjectives satisfy either the left or the right subformula
storm::storage::BitVector constraintStates(model.getNumberOfStates(), true);
for (auto const& dim : objectiveDimensions[objIndex]) {
auto const& dimension = dimensions[dim];
STORM_LOG_ASSERT(dimension.formula->isBoundedUntilFormula(), "Unexpected Formula type");
constraintStates &=
(mc.check(dimension.formula->asBoundedUntilFormula().getLeftSubformula())->asExplicitQualitativeCheckResult().getTruthValuesVector() |
mc.check(dimension.formula->asBoundedUntilFormula().getRightSubformula())->asExplicitQualitativeCheckResult().getTruthValuesVector());
}
// Build the transitions between the memory states
for (uint64_t memState = 0; memState < objMemStates.size(); ++memState) {
auto const& memStateBV = objMemStates[memState];
for (uint64_t memStatePrime = 0; memStatePrime < objMemStates.size(); ++memStatePrime) {
auto const& memStatePrimeBV = objMemStates[memStatePrime];
if (memStatePrimeBV.isSubsetOf(memStateBV)) {
std::shared_ptr<storm::logic::Formula const> transitionFormula = storm::logic::Formula::getTrueFormula();
for (auto const& subObjIndex : memStateBV) {
std::shared_ptr<storm::logic::Formula const> subObjFormula = dimensions[dimensionIndexMap[subObjIndex]].formula->asBoundedUntilFormula().getRightSubformula().asSharedPointer();
if (memStatePrimeBV.get(subObjIndex)) {
subObjFormula = std::make_shared<storm::logic::UnaryBooleanStateFormula>(storm::logic::UnaryBooleanStateFormula::OperatorType::Not, subObjFormula);
}
transitionFormula = std::make_shared<storm::logic::BinaryBooleanStateFormula>(storm::logic::BinaryBooleanStateFormula::OperatorType::And, transitionFormula, subObjFormula);
}
storm::storage::BitVector transitionStates = mc.check(*transitionFormula)->asExplicitQualitativeCheckResult().getTruthValuesVector();
if (memStatePrimeBV.empty()) {
transitionStates |= ~constraintStates;
} else {
transitionStates &= constraintStates;
}
objMemoryBuilder.setTransition(memState, memStatePrime, transitionStates);
// Set the initial states
if (memStateBV.full()) {
storm::storage::BitVector initialTransitionStates = model.getInitialStates() & transitionStates;
// At this point we can check whether there is an initial state that already satisfies all subObjectives.
// Such a situation is not supported as we can not reduce this (easily) to an expected reward computation.
STORM_LOG_THROW(!memStatePrimeBV.empty() || initialTransitionStates.empty() || initialTransitionStates.isDisjointFrom(constraintStates), storm::exceptions::NotSupportedException, "The objective " << *objectives[objIndex].formula << " is already satisfied in an initial state. This special case is not supported.");
for (auto const& initState : initialTransitionStates) {
objMemoryBuilder.setInitialMemoryState(initState, memStatePrime);
}
}
}
}
}
// Build the memory labels
for (uint64_t memState = 0; memState < objMemStates.size(); ++memState) {
auto const& memStateBV = objMemStates[memState];
for (auto const& subObjIndex : memStateBV) {
objMemoryBuilder.setLabel(memState, dimensions[dimensionIndexMap[subObjIndex]].memoryLabel.get());
}
}
auto objMemory = objMemoryBuilder.build();
memory = memory.product(objMemory);
}
return memory;
}
template<typename ValueType, bool SingleObjectiveMode>
std::vector<storm::storage::BitVector> MultiDimensionalRewardUnfolding<ValueType, SingleObjectiveMode>::computeMemoryStateMap(storm::storage::MemoryStructure const& memory) const {
// Compute a mapping between the different representations of memory states
std::vector<storm::storage::BitVector> result;
result.reserve(memory.getNumberOfStates());
for (uint64_t memState = 0; memState < memory.getNumberOfStates(); ++memState) {
storm::storage::BitVector relevantSubObjectives(epochManager.getDimensionCount(), false);
std::set<std::string> stateLabels = memory.getStateLabeling().getLabelsOfState(memState);
for (uint64_t dim = 0; dim < epochManager.getDimensionCount(); ++dim) {
if (dimensions[dim].memoryLabel && stateLabels.find(dimensions[dim].memoryLabel.get()) != stateLabels.end()) {
relevantSubObjectives.set(dim, true);
}
}
result.push_back(std::move(relevantSubObjectives));
}
return result;
STORM_LOG_ASSERT(!epochManager.hasBottomDimension(epoch), "Tried to get the initial state result in an epoch that still contains at least one bottom dimension.");
return getStateSolution(epoch, productModel->getInitialProductState(initialStateIndex, model.getInitialStates()));
}
template class MultiDimensionalRewardUnfolding<double, true>;

2
src/storm/modelchecker/multiobjective/rewardbounded/MultiDimensionalRewardUnfolding.h

@ -90,8 +90,6 @@ namespace storm {
void initializeMemoryProduct(std::vector<Epoch> const& epochSteps);
storm::storage::MemoryStructure computeMemoryStructure() const;
std::vector<storm::storage::BitVector> computeMemoryStateMap(storm::storage::MemoryStructure const& memory) const;
template<bool SO = SingleObjectiveMode, typename std::enable_if<SO, int>::type = 0>
SolutionType getScaledSolution(SolutionType const& solution, ValueType const& scalingFactor) const;

236
src/storm/modelchecker/multiobjective/rewardbounded/ProductModel.cpp

@ -5,6 +5,7 @@
#include "storm/logic/Formulas.h"
#include "storm/logic/CloneVisitor.h"
#include "storm/storage/memorystructure/SparseModelMemoryProduct.h"
#include "storm/storage/memorystructure/MemoryStructureBuilder.h"
#include "storm/modelchecker/propositional/SparsePropositionalModelChecker.h"
#include "storm/modelchecker/results/ExplicitQualitativeCheckResult.h"
@ -17,15 +18,19 @@ namespace storm {
namespace multiobjective {
template<typename ValueType>
ProductModel<ValueType>::ProductModel(storm::models::sparse::Mdp<ValueType> const& model, storm::storage::MemoryStructure const& memory, std::vector<Dimension<ValueType>> const& dimensions, std::vector<storm::storage::BitVector> const& objectiveDimensions, EpochManager const& epochManager, std::vector<storm::storage::BitVector>&& memoryStateMap, std::vector<Epoch> const& originalModelSteps) : dimensions(dimensions), objectiveDimensions(objectiveDimensions), epochManager(epochManager), memoryStateMap(std::move(memoryStateMap)) {
ProductModel<ValueType>::ProductModel(storm::models::sparse::Mdp<ValueType> const& model, std::vector<storm::modelchecker::multiobjective::Objective<ValueType>> const& objectives, std::vector<Dimension<ValueType>> const& dimensions, std::vector<storm::storage::BitVector> const& objectiveDimensions, EpochManager const& epochManager, std::vector<Epoch> const& originalModelSteps) : dimensions(dimensions), objectiveDimensions(objectiveDimensions), epochManager(epochManager), memoryStateManager(dimensions.size()) {
storm::storage::MemoryStructure memory = computeMemoryStructure(model, objectives);
assert(memoryStateManager.getMemoryStateCount() == memory.getNumberOfStates());
std::vector<MemoryState> memoryStateMap = computeMemoryStateMap(memory);
storm::storage::SparseModelMemoryProduct<ValueType> productBuilder(memory.product(model));
setReachableProductStates(productBuilder, originalModelSteps);
setReachableProductStates(productBuilder, originalModelSteps, memoryStateMap);
product = productBuilder.build()->template as<storm::models::sparse::Mdp<ValueType>>();
uint64_t numModelStates = productBuilder.getOriginalModel().getNumberOfStates();
uint64_t numMemoryStates = productBuilder.getMemory().getNumberOfStates();
uint64_t numMemoryStates = memoryStateManager.getMemoryStateCount();
uint64_t numProductStates = getProduct().getNumberOfStates();
// Compute a mappings from product states to model/memory states and back
@ -33,12 +38,12 @@ namespace storm {
productToModelStateMap.resize(numProductStates, std::numeric_limits<uint64_t>::max());
productToMemoryStateMap.resize(numProductStates, std::numeric_limits<uint64_t>::max());
for (uint64_t modelState = 0; modelState < numModelStates; ++modelState) {
for (uint64_t memoryState = 0; memoryState < numMemoryStates; ++memoryState) {
if (productBuilder.isStateReachable(modelState, memoryState)) {
uint64_t productState = productBuilder.getResultState(modelState, memoryState);
modelMemoryToProductStateMap[modelState * numMemoryStates + memoryState] = productState;
for (uint64_t memoryStateIndex = 0; memoryStateIndex < numMemoryStates; ++memoryStateIndex) {
if (productBuilder.isStateReachable(modelState, memoryStateIndex)) {
uint64_t productState = productBuilder.getResultState(modelState, memoryStateIndex);
modelMemoryToProductStateMap[modelState * numMemoryStates + memoryStateMap[memoryStateIndex]] = productState;
productToModelStateMap[productState] = modelState;
productToMemoryStateMap[productState] = memoryState;
productToMemoryStateMap[productState] = memoryStateMap[memoryStateIndex];
}
}
}
@ -60,9 +65,9 @@ namespace storm {
for (uint64_t choiceOffset = 0; choiceOffset < numChoices; ++choiceOffset) {
Epoch const& step = originalModelSteps[firstChoice + choiceOffset];
if (step != 0) {
for (uint64_t memState = 0; memState < numMemoryStates; ++memState) {
if (productStateExists(modelState, memState)) {
uint64_t productState = getProductState(modelState, memState);
for (MemoryState const& memoryState : memoryStateMap) {
if (productStateExists(modelState, memoryState)) {
uint64_t productState = getProductState(modelState, memoryState);
uint64_t productChoice = getProduct().getTransitionMatrix().getRowGroupIndices()[productState] + choiceOffset;
assert(productChoice < getProduct().getTransitionMatrix().getRowGroupIndices()[productState + 1]);
steps[productChoice] = step;
@ -71,24 +76,150 @@ namespace storm {
}
}
}
// getProduct().writeDotToStream(std::cout);
computeReachableStatesInEpochClasses();
}
template<typename ValueType>
storm::storage::MemoryStructure ProductModel<ValueType>::computeMemoryStructure(storm::models::sparse::Mdp<ValueType> const& model, std::vector<storm::modelchecker::multiobjective::Objective<ValueType>> const& objectives) const {
storm::modelchecker::SparsePropositionalModelChecker<storm::models::sparse::Mdp<ValueType>> mc(model);
// Create a memory structure that remembers whether (sub)objectives are satisfied
storm::storage::MemoryStructure memory = storm::storage::MemoryStructureBuilder<ValueType>::buildTrivialMemoryStructure(model);
for (uint64_t objIndex = 0; objIndex < objectives.size(); ++objIndex) {
if (!objectives[objIndex].formula->isProbabilityOperatorFormula()) {
continue;
}
std::vector<uint64_t> dimensionIndexMap;
for (auto const& globalDimensionIndex : objectiveDimensions[objIndex]) {
dimensionIndexMap.push_back(globalDimensionIndex);
}
bool objectiveContainsLowerBound = false;
for (auto const& globalDimensionIndex : objectiveDimensions[objIndex]) {
if (!dimensions[globalDimensionIndex].isUpperBounded) {
objectiveContainsLowerBound = true;
break;
}
}
// collect the memory states for this objective
std::vector<storm::storage::BitVector> objMemStates;
storm::storage::BitVector m(dimensionIndexMap.size(), false);
for (; !m.full(); m.increment()) {
objMemStates.push_back(~m);
}
objMemStates.push_back(~m);
assert(objMemStates.size() == 1ull << dimensionIndexMap.size());
// build objective memory
auto objMemoryBuilder = storm::storage::MemoryStructureBuilder<ValueType>(objMemStates.size(), model);
// Get the set of states that for all subobjectives satisfy either the left or the right subformula
storm::storage::BitVector constraintStates(model.getNumberOfStates(), true);
for (auto const& dim : objectiveDimensions[objIndex]) {
auto const& dimension = dimensions[dim];
STORM_LOG_ASSERT(dimension.formula->isBoundedUntilFormula(), "Unexpected Formula type");
constraintStates &=
(mc.check(dimension.formula->asBoundedUntilFormula().getLeftSubformula())->asExplicitQualitativeCheckResult().getTruthValuesVector() |
mc.check(dimension.formula->asBoundedUntilFormula().getRightSubformula())->asExplicitQualitativeCheckResult().getTruthValuesVector());
}
// Build the transitions between the memory states
for (uint64_t memState = 0; memState < objMemStates.size(); ++memState) {
auto const& memStateBV = objMemStates[memState];
for (uint64_t memStatePrime = 0; memStatePrime < objMemStates.size(); ++memStatePrime) {
auto const& memStatePrimeBV = objMemStates[memStatePrime];
if (memStatePrimeBV.isSubsetOf(memStateBV)) {
std::shared_ptr<storm::logic::Formula const> transitionFormula = storm::logic::Formula::getTrueFormula();
for (auto const& subObjIndex : memStateBV) {
std::shared_ptr<storm::logic::Formula const> subObjFormula = dimensions[dimensionIndexMap[subObjIndex]].formula->asBoundedUntilFormula().getRightSubformula().asSharedPointer();
if (memStatePrimeBV.get(subObjIndex)) {
subObjFormula = std::make_shared<storm::logic::UnaryBooleanStateFormula>(storm::logic::UnaryBooleanStateFormula::OperatorType::Not, subObjFormula);
}
transitionFormula = std::make_shared<storm::logic::BinaryBooleanStateFormula>(storm::logic::BinaryBooleanStateFormula::OperatorType::And, transitionFormula, subObjFormula);
}
storm::storage::BitVector transitionStates = mc.check(*transitionFormula)->asExplicitQualitativeCheckResult().getTruthValuesVector();
if (memStatePrimeBV.empty()) {
transitionStates |= ~constraintStates;
} else {
transitionStates &= constraintStates;
}
objMemoryBuilder.setTransition(memState, memStatePrime, transitionStates);
// Set the initial states
if (memStateBV.full()) {
storm::storage::BitVector initialTransitionStates = model.getInitialStates() & transitionStates;
// At this point we can check whether there is an initial state that already satisfies all subObjectives.
// Such a situation is not supported as we can not reduce this (easily) to an expected reward computation.
STORM_LOG_THROW(!memStatePrimeBV.empty() || initialTransitionStates.empty() || objectiveContainsLowerBound || initialTransitionStates.isDisjointFrom(constraintStates), storm::exceptions::NotSupportedException, "The objective " << *objectives[objIndex].formula << " is already satisfied in an initial state. This special case is not supported.");
for (auto const& initState : initialTransitionStates) {
objMemoryBuilder.setInitialMemoryState(initState, memStatePrime);
}
}
}
}
}
// Build the memory labels
for (uint64_t memState = 0; memState < objMemStates.size(); ++memState) {
auto const& memStateBV = objMemStates[memState];
for (auto const& subObjIndex : memStateBV) {
objMemoryBuilder.setLabel(memState, dimensions[dimensionIndexMap[subObjIndex]].memoryLabel.get());
}
}
auto objMemory = objMemoryBuilder.build();
memory = memory.product(objMemory);
}
return memory;
}
template<typename ValueType>
std::vector<typename ProductModel<ValueType>::MemoryState> ProductModel<ValueType>::computeMemoryStateMap(storm::storage::MemoryStructure const& memory) const {
// Compute a mapping between the different representations of memory states
std::vector<MemoryState> result;
result.reserve(memory.getNumberOfStates());
for (uint64_t memStateIndex = 0; memStateIndex < memory.getNumberOfStates(); ++memStateIndex) {
MemoryState memState = memoryStateManager.getInitialMemoryState();
std::set<std::string> stateLabels = memory.getStateLabeling().getLabelsOfState(memStateIndex);
for (uint64_t dim = 0; dim < epochManager.getDimensionCount(); ++dim) {
if (dimensions[dim].memoryLabel && stateLabels.find(dimensions[dim].memoryLabel.get()) != stateLabels.end()) {
memoryStateManager.setRelevantDimension(memState, dim, true);
} else {
memoryStateManager.setRelevantDimension(memState, dim, false);
}
}
result.push_back(std::move(memState));
}
return result;
}
template<typename ValueType>
void ProductModel<ValueType>::setReachableProductStates(storm::storage::SparseModelMemoryProduct<ValueType>& productBuilder, std::vector<Epoch> const& originalModelSteps) const {
void ProductModel<ValueType>::setReachableProductStates(storm::storage::SparseModelMemoryProduct<ValueType>& productBuilder, std::vector<Epoch> const& originalModelSteps, std::vector<MemoryState> const& memoryStateMap) const {
std::vector<uint64_t> inverseMemoryStateMap(memoryStateMap.size(), std::numeric_limits<uint64_t>::max());
for (uint64_t memStateIndex = 0; memStateIndex < memoryStateMap.size(); ++memStateIndex) {
inverseMemoryStateMap[memoryStateMap[memStateIndex]] = memStateIndex;
}
auto const& memory = productBuilder.getMemory();
auto const& model = productBuilder.getOriginalModel();
auto const& modelTransitions = model.getTransitionMatrix();
std::vector<storm::storage::BitVector> reachableStates(memory.getNumberOfStates(), storm::storage::BitVector(model.getNumberOfStates(), false));
std::vector<storm::storage::BitVector> reachableStates(memoryStateManager.getMemoryStateCount(), storm::storage::BitVector(model.getNumberOfStates(), false));
// Initialize the reachable states with the initial states
EpochClass initEpochClass = epochManager.getEpochClass(epochManager.getZeroEpoch());
storm::storage::BitVector initMemState(epochManager.getDimensionCount(), true);
auto memStateIt = memory.getInitialMemoryStates().begin();
for (auto const& initState : model.getInitialStates()) {
uint64_t transformedMemoryState = transformMemoryState(*memStateIt, initEpochClass, convertMemoryState(initMemState));
uint64_t transformedMemoryState = transformMemoryState(memoryStateMap[*memStateIt], initEpochClass, memoryStateManager.getInitialMemoryState());
reachableStates[transformedMemoryState].set(initState, true);
++memStateIt;
}
@ -105,8 +236,8 @@ namespace storm {
auto const& epochClass = *epochClassIt;
// Find the remaining set of reachable states via DFS.
std::vector<std::pair<uint64_t, uint64_t>> dfsStack;
for (uint64_t memState = 0; memState < memory.getNumberOfStates(); ++memState) {
std::vector<std::pair<uint64_t, MemoryState>> dfsStack;
for (MemoryState const& memState : memoryStateMap) {
for (auto const& modelState : reachableStates[memState]) {
dfsStack.emplace_back(modelState, memState);
}
@ -114,14 +245,15 @@ namespace storm {
while (!dfsStack.empty()) {
uint64_t currentModelState = dfsStack.back().first;
uint64_t currentMemoryState = dfsStack.back().second;
MemoryState currentMemoryState = dfsStack.back().second;
uint64_t currentMemoryStateIndex = inverseMemoryStateMap[currentMemoryState];
dfsStack.pop_back();
for (uint64_t choice = modelTransitions.getRowGroupIndices()[currentModelState]; choice != modelTransitions.getRowGroupIndices()[currentModelState + 1]; ++choice) {
for (auto transitionIt = modelTransitions.getRow(choice).begin(); transitionIt < modelTransitions.getRow(choice).end(); ++transitionIt) {
uint64_t successorMemoryState = memory.getSuccessorMemoryState(currentMemoryState, transitionIt - modelTransitions.begin());
MemoryState successorMemoryState = memoryStateMap[memory.getSuccessorMemoryState(currentMemoryStateIndex, transitionIt - modelTransitions.begin())];
successorMemoryState = transformMemoryState(successorMemoryState, epochClass, currentMemoryState);
if (!reachableStates[successorMemoryState].get(transitionIt->getColumn())) {
reachableStates[successorMemoryState].set(transitionIt->getColumn(), true);
@ -131,9 +263,10 @@ namespace storm {
}
}
}
for (uint64_t memState = 0; memState < memoryStateMap.size(); ++memState) {
for (auto const& modelState : reachableStates[memState]) {
productBuilder.addReachableState(modelState, memState);
for (uint64_t memStateIndex = 0; memStateIndex < memoryStateManager.getMemoryStateCount(); ++memStateIndex) {
for (auto const& modelState : reachableStates[memoryStateMap[memStateIndex]]) {
productBuilder.addReachableState(modelState, memStateIndex);
}
}
}
@ -149,44 +282,37 @@ namespace storm {
}
template<typename ValueType>
bool ProductModel<ValueType>::productStateExists(uint64_t const& modelState, uint64_t const& memoryState) const {
STORM_LOG_ASSERT(!memoryStateMap.empty(), "Tried to retrieve whether a product state exists but the memoryStateMap is not yet initialized.");
return modelMemoryToProductStateMap[modelState * memoryStateMap.size() + memoryState] < getProduct().getNumberOfStates();
bool ProductModel<ValueType>::productStateExists(uint64_t const& modelState, MemoryState const& memoryState) const {
return modelMemoryToProductStateMap[modelState * memoryStateManager.getMemoryStateCount() + memoryState] < getProduct().getNumberOfStates();
}
template<typename ValueType>
uint64_t ProductModel<ValueType>::getProductState(uint64_t const& modelState, uint64_t const& memoryState) const {
uint64_t ProductModel<ValueType>::getProductState(uint64_t const& modelState, MemoryState const& memoryState) const {
STORM_LOG_ASSERT(productStateExists(modelState, memoryState), "Tried to obtain a state in the model-memory-product that does not exist");
return modelMemoryToProductStateMap[modelState * memoryStateMap.size() + memoryState];
return modelMemoryToProductStateMap[modelState * memoryStateManager.getMemoryStateCount() + memoryState];
}
template<typename ValueType>
uint64_t ProductModel<ValueType>::getModelState(uint64_t const& productState) const {
return productToModelStateMap[productState];
uint64_t ProductModel<ValueType>::getInitialProductState(uint64_t const& initialModelState, storm::storage::BitVector const& initialModelStates) const {
auto productInitStateIt = getProduct().getInitialStates().begin();
productInitStateIt += initialModelStates.getNumberOfSetBitsBeforeIndex(initialModelState);
STORM_LOG_ASSERT(getModelState(*productInitStateIt) == initialModelState, "Could not find the corresponding initial state in the product model.");
return transformProductState(*productInitStateIt, epochManager.getEpochClass(epochManager.getZeroEpoch()), memoryStateManager.getInitialMemoryState());
}
template<typename ValueType>
uint64_t ProductModel<ValueType>::getMemoryState(uint64_t const& productState) const {
return productToMemoryStateMap[productState];
}
template<typename ValueType>
storm::storage::BitVector const& ProductModel<ValueType>::convertMemoryState(uint64_t const& memoryState) const {
STORM_LOG_ASSERT(!memoryStateMap.empty(), "Tried to convert a memory state, but the memoryStateMap is not yet initialized.");
return memoryStateMap[memoryState];
uint64_t ProductModel<ValueType>::getModelState(uint64_t const& productState) const {
return productToModelStateMap[productState];
}
template<typename ValueType>
uint64_t ProductModel<ValueType>::convertMemoryState(storm::storage::BitVector const& memoryState) const {
STORM_LOG_ASSERT(!memoryStateMap.empty(), "Tried to convert a memory state, but the memoryStateMap is not yet initialized.");
auto memStateIt = std::find(memoryStateMap.begin(), memoryStateMap.end(), memoryState);
return memStateIt - memoryStateMap.begin();
typename ProductModel<ValueType>::MemoryState ProductModel<ValueType>::getMemoryState(uint64_t const& productState) const {
return productToMemoryStateMap[productState];
}
template<typename ValueType>
uint64_t ProductModel<ValueType>::getNumberOfMemoryState() const {
STORM_LOG_ASSERT(!memoryStateMap.empty(), "Tried to retrieve the number of memory states but the memoryStateMap is not yet initialized.");
return memoryStateMap.size();
MemoryStateManager const& ProductModel<ValueType>::getMemoryStateManager() const {
return memoryStateManager;
}
template<typename ValueType>
@ -352,9 +478,8 @@ namespace storm {
storm::storage::BitVector ecInStates(getProduct().getNumberOfStates(), false);
if (!epochManager.hasBottomDimensionEpochClass(epochClass)) {
uint64_t initMemState = convertMemoryState(storm::storage::BitVector(epochManager.getDimensionCount(), true));
for (auto const& initState : getProduct().getInitialStates()) {
uint64_t transformedInitState = transformProductState(initState, epochClass, initMemState);
uint64_t transformedInitState = transformProductState(initState, epochClass, memoryStateManager.getInitialMemoryState());
ecInStates.set(transformedInitState, true);
}
}
@ -429,30 +554,31 @@ namespace storm {
}
template<typename ValueType>
uint64_t ProductModel<ValueType>::transformMemoryState(uint64_t const& memoryState, EpochClass const& epochClass, uint64_t const& predecessorMemoryState) const {
storm::storage::BitVector memoryStateBv = convertMemoryState(memoryState);
storm::storage::BitVector const& predecessorMemoryStateBv = convertMemoryState(predecessorMemoryState);
typename ProductModel<ValueType>::MemoryState ProductModel<ValueType>::transformMemoryState(MemoryState const& memoryState, EpochClass const& epochClass, MemoryState const& predecessorMemoryState) const {
MemoryState memoryStatePrime = memoryState;
for (auto const& objDimensions : objectiveDimensions) {
for (auto const& dim : objDimensions) {
auto const& dimension = dimensions[dim];
bool dimUpperBounded = dimension.isUpperBounded;
bool dimBottom = epochManager.isBottomDimensionEpochClass(epochClass, dim);
if (dimUpperBounded && dimBottom && predecessorMemoryStateBv.get(dim)) {
if (dimUpperBounded && dimBottom && memoryStateManager.isRelevantDimension(predecessorMemoryState, dim)) {
STORM_LOG_ASSERT(objDimensions == dimension.dependentDimensions, "Unexpected set of dependent dimensions");
memoryStateBv &= ~objDimensions;
memoryStateManager.setRelevantDimensions(memoryStatePrime, objDimensions, false);
break;
} else if (!dimUpperBounded && !dimBottom && predecessorMemoryStateBv.get(dim)) {
memoryStateBv |= dimension.dependentDimensions;
} else if (!dimUpperBounded && !dimBottom && memoryStateManager.isRelevantDimension(predecessorMemoryState, dim)) {
memoryStateManager.setRelevantDimensions(memoryStatePrime, dimension.dependentDimensions, true);
}
}
}
// std::cout << "Transformed memory state " << memoryStateManager.toString(memoryState) << " at epoch class " << epochClass << " with predecessor " << memoryStateManager.toString(predecessorMemoryState) << " to " << memoryStateManager.toString(memoryStatePrime) << std::endl;
return convertMemoryState(memoryStateBv);
return memoryStatePrime;
}
template<typename ValueType>
uint64_t ProductModel<ValueType>::transformProductState(uint64_t const& productState, EpochClass const& epochClass, uint64_t const& predecessorMemoryState) const {
uint64_t ProductModel<ValueType>::transformProductState(uint64_t const& productState, EpochClass const& epochClass, MemoryState const& predecessorMemoryState) const {
return getProductState(getModelState(productState), transformMemoryState(getMemoryState(productState), epochClass, predecessorMemoryState));
}

29
src/storm/modelchecker/multiobjective/rewardbounded/ProductModel.h

@ -5,6 +5,7 @@
#include "storm/storage/BitVector.h"
#include "storm/modelchecker/multiobjective/Objective.h"
#include "storm/modelchecker/multiobjective/rewardbounded/EpochManager.h"
#include "storm/modelchecker/multiobjective/rewardbounded/MemoryStateManager.h"
#include "storm/modelchecker/multiobjective/rewardbounded/Dimension.h"
#include "storm/models/sparse/Mdp.h"
#include "storm/utility/vector.h"
@ -22,36 +23,38 @@ namespace storm {
typedef typename EpochManager::Epoch Epoch;
typedef typename EpochManager::EpochClass EpochClass;
typedef typename MemoryStateManager::MemoryState MemoryState;
ProductModel(storm::models::sparse::Mdp<ValueType> const& model, storm::storage::MemoryStructure const& memory, std::vector<Dimension<ValueType>> const& dimensions, std::vector<storm::storage::BitVector> const& objectiveDimensions, EpochManager const& epochManager, std::vector<storm::storage::BitVector>&& memoryStateMap, std::vector<Epoch> const& originalModelSteps);
ProductModel(storm::models::sparse::Mdp<ValueType> const& model, std::vector<storm::modelchecker::multiobjective::Objective<ValueType>> const& objectives, std::vector<Dimension<ValueType>> const& dimensions, std::vector<storm::storage::BitVector> const& objectiveDimensions, EpochManager const& epochManager, std::vector<Epoch> const& originalModelSteps);
storm::models::sparse::Mdp<ValueType> const& getProduct() const;
std::vector<Epoch> const& getSteps() const;
bool productStateExists(uint64_t const& modelState, uint64_t const& memoryState) const;
uint64_t getProductState(uint64_t const& modelState, uint64_t const& memoryState) const;
uint64_t getInitialProductState(uint64_t const& initialModelState, storm::storage::BitVector const& initialModelStates) const;
uint64_t getModelState(uint64_t const& productState) const;
uint64_t getMemoryState(uint64_t const& productState) const;
uint64_t convertMemoryState(storm::storage::BitVector const& memoryState) const;
storm::storage::BitVector const& convertMemoryState(uint64_t const& memoryState) const;
uint64_t getNumberOfMemoryState() const;
MemoryState getMemoryState(uint64_t const& productState) const;
MemoryStateManager const& getMemoryStateManager() const;
uint64_t getProductStateFromChoice(uint64_t const& productChoice) const;
std::vector<std::vector<ValueType>> computeObjectiveRewards(EpochClass const& epochClass, std::vector<storm::modelchecker::multiobjective::Objective<ValueType>> const& objectives) const;
storm::storage::BitVector const& getInStates(EpochClass const& epochClass) const;
uint64_t transformMemoryState(uint64_t const& memoryState, EpochClass const& epochClass, uint64_t const& predecessorMemoryState) const;
uint64_t transformProductState(uint64_t const& productState, EpochClass const& epochClass, uint64_t const& predecessorMemoryState) const;
MemoryState transformMemoryState(MemoryState const& memoryState, EpochClass const& epochClass, MemoryState const& predecessorMemoryState) const;
uint64_t transformProductState(uint64_t const& productState, EpochClass const& epochClass, MemoryState const& predecessorMemoryState) const;
private:
void setReachableProductStates(storm::storage::SparseModelMemoryProduct<ValueType>& productBuilder, std::vector<Epoch> const& originalModelSteps) const;
storm::storage::MemoryStructure computeMemoryStructure(storm::models::sparse::Mdp<ValueType> const& model, std::vector<storm::modelchecker::multiobjective::Objective<ValueType>> const& objectives) const;
std::vector<MemoryState> computeMemoryStateMap(storm::storage::MemoryStructure const& memory) const;
void setReachableProductStates(storm::storage::SparseModelMemoryProduct<ValueType>& productBuilder, std::vector<Epoch> const& originalModelSteps, std::vector<MemoryState> const& memoryStateMap) const;
void collectReachableEpochClasses(std::set<EpochClass, std::function<bool(EpochClass const&, EpochClass const&)>>& reachableEpochClasses, std::set<Epoch> const& possibleSteps) const;
@ -61,6 +64,7 @@ namespace storm {
std::vector<Dimension<ValueType>> const& dimensions;
std::vector<storm::storage::BitVector> const& objectiveDimensions;
EpochManager const& epochManager;
MemoryStateManager memoryStateManager;
std::shared_ptr<storm::models::sparse::Mdp<ValueType>> product;
std::vector<Epoch> steps;
@ -69,9 +73,8 @@ namespace storm {
std::vector<uint64_t> modelMemoryToProductStateMap;
std::vector<uint64_t> productToModelStateMap;
std::vector<uint64_t> productToMemoryStateMap;
std::vector<MemoryState> productToMemoryStateMap;
std::vector<uint64_t> choiceToStateMap;
std::vector<storm::storage::BitVector> memoryStateMap;
};
}

10
src/test/storm/modelchecker/SparseMdpMultiDimensionalRewardUnfoldingTest.cpp

@ -213,6 +213,7 @@ TEST(SparseMdpMultiDimensionalRewardUnfoldingTest, single_obj_lower_bounds) {
std::string formulasAsString = "Pmax=? [ F{\"a\"}>=1.1 x=1 ]";
formulasAsString += "; \n Pmin=? [ F{\"a\"}>=1.1 x=1 ]";
formulasAsString += "; \n Pmax=? [ F{\"b\"}>3 x=1 ]";
formulasAsString += "; \n Pmax=? [ F{\"b\"}>=2 x=0 ]";
formulasAsString += "; \n Pmax=? [ multi(F{\"a\"}<=0.2 x=1, F{\"b\"}>3 x=1) ]";
formulasAsString += "; \n Pmax=? [ multi(F{\"a\"}>0.2 x=1, F{\"b\"}>3 x=1) ]";
formulasAsString += "; \n Pmax=? [ multi(F{\"a\"}>=2/5 x=4, F{\"a\"}<=0 x=4) ]";
@ -245,16 +246,21 @@ TEST(SparseMdpMultiDimensionalRewardUnfoldingTest, single_obj_lower_bounds) {
result = storm::api::verifyWithSparseEngine(mdp, storm::api::createTask<storm::RationalNumber>(formulas[3], true));
ASSERT_TRUE(result->isExplicitQuantitativeCheckResult());
expectedResult = storm::utility::convertNumber<storm::RationalNumber, std::string>("27/64");
expectedResult = storm::utility::convertNumber<storm::RationalNumber, std::string>("9/16");
EXPECT_EQ(expectedResult, result->asExplicitQuantitativeCheckResult<storm::RationalNumber>()[initState]);
result = storm::api::verifyWithSparseEngine(mdp, storm::api::createTask<storm::RationalNumber>(formulas[4], true));
ASSERT_TRUE(result->isExplicitQuantitativeCheckResult());
expectedResult = storm::utility::convertNumber<storm::RationalNumber, std::string>("243/640");
expectedResult = storm::utility::convertNumber<storm::RationalNumber, std::string>("27/64");
EXPECT_EQ(expectedResult, result->asExplicitQuantitativeCheckResult<storm::RationalNumber>()[initState]);
result = storm::api::verifyWithSparseEngine(mdp, storm::api::createTask<storm::RationalNumber>(formulas[5], true));
ASSERT_TRUE(result->isExplicitQuantitativeCheckResult());
expectedResult = storm::utility::convertNumber<storm::RationalNumber, std::string>("243/640");
EXPECT_EQ(expectedResult, result->asExplicitQuantitativeCheckResult<storm::RationalNumber>()[initState]);
result = storm::api::verifyWithSparseEngine(mdp, storm::api::createTask<storm::RationalNumber>(formulas[6], true));
ASSERT_TRUE(result->isExplicitQuantitativeCheckResult());
expectedResult = storm::utility::convertNumber<storm::RationalNumber, std::string>("81/160");
EXPECT_EQ(expectedResult, result->asExplicitQuantitativeCheckResult<storm::RationalNumber>()[initState]);

Loading…
Cancel
Save