Browse Source

building markov automata from prism code

Former-commit-id: 791c49c7cf
tempestpy_adaptions
TimQu 8 years ago
parent
commit
f681206393
  1. 70
      src/builder/ExplicitModelBuilder.cpp
  2. 19
      src/builder/ExplicitModelBuilder.h
  3. 4
      src/generator/PrismNextStateGenerator.cpp
  4. 51
      src/models/sparse/MarkovAutomaton.cpp
  5. 41
      src/models/sparse/MarkovAutomaton.h
  6. 10
      src/storage/BitVector.cpp
  7. 11
      src/storage/BitVector.h
  8. 64
      src/storage/SparseMatrix.cpp
  9. 8
      src/storage/SparseMatrix.h
  10. 2
      src/utility/ModelInstantiator.h
  11. 26
      test/functional/builder/ExplicitPrismModelBuilderTest.cpp
  12. 19
      test/functional/builder/hybrid_states.ma
  13. 15
      test/functional/builder/simple.ma
  14. 50
      test/functional/builder/stream2.ma

70
src/builder/ExplicitModelBuilder.cpp

@ -5,6 +5,7 @@
#include "src/models/sparse/Dtmc.h"
#include "src/models/sparse/Ctmc.h"
#include "src/models/sparse/Mdp.h"
#include "src/models/sparse/MarkovAutomaton.h"
#include "src/models/sparse/StandardRewardModel.h"
#include "src/storage/expressions/ExpressionManager.h"
@ -133,6 +134,9 @@ namespace storm {
case storm::generator::ModelType::MDP:
result = std::shared_ptr<storm::models::sparse::Model<ValueType, RewardModelType>>(new storm::models::sparse::Mdp<ValueType, RewardModelType>(std::move(modelComponents.transitionMatrix), std::move(modelComponents.stateLabeling), std::move(modelComponents.rewardModels), std::move(modelComponents.choiceLabeling)));
break;
case storm::generator::ModelType::MA:
result = std::shared_ptr<storm::models::sparse::Model<ValueType, RewardModelType>>(new storm::models::sparse::MarkovAutomaton<ValueType, RewardModelType>(std::move(modelComponents.transitionMatrix), std::move(modelComponents.stateLabeling), *std::move(modelComponents.markovianStates), std::move(modelComponents.rewardModels), std::move(modelComponents.choiceLabeling)));
break;
default:
STORM_LOG_THROW(false, storm::exceptions::WrongFormatException, "Error while creating model: cannot handle this model type.");
break;
@ -165,12 +169,17 @@ namespace storm {
}
template <typename ValueType, typename RewardModelType, typename StateType>
boost::optional<std::vector<boost::container::flat_set<uint_fast64_t>>> ExplicitModelBuilder<ValueType, RewardModelType, StateType>::buildMatrices(storm::storage::SparseMatrixBuilder<ValueType>& transitionMatrixBuilder, std::vector<RewardModelBuilder<typename RewardModelType::ValueType>>& rewardModelBuilders) {
void ExplicitModelBuilder<ValueType, RewardModelType, StateType>::buildMatrices(storm::storage::SparseMatrixBuilder<ValueType>& transitionMatrixBuilder, std::vector<RewardModelBuilder<typename RewardModelType::ValueType>>& rewardModelBuilders, boost::optional<std::vector<boost::container::flat_set<uint_fast64_t>>>& choiceLabels , boost::optional<storm::storage::BitVector>& markovianChoices) {
// Create choice labels, if requested,
boost::optional<std::vector<boost::container::flat_set<uint_fast64_t>>> choiceLabels;
if (generator->getOptions().isBuildChoiceLabelsSet()) {
choiceLabels = std::vector<boost::container::flat_set<uint_fast64_t>>();
}
// Create markovian states bit vector, if required
if (generator->getModelType() == storm::generator::ModelType::MA) {
// The BitVector will be resized when the correct size is known
markovianChoices = storm::storage::BitVector(64, false);
}
// Create a callback for the next-state generator to enable it to request the index of states.
std::function<StateType (CompressedState const&)> stateToIdCallback = std::bind(&ExplicitModelBuilder<ValueType, RewardModelType, StateType>::getOrAddStateIndex, this, std::placeholders::_1);
@ -219,6 +228,12 @@ namespace storm {
// Insert empty choice labeling for added self-loop transitions.
choiceLabels.get().push_back(boost::container::flat_set<uint_fast64_t>());
}
if (generator->getModelType() == storm::generator::ModelType::MA) {
markovianChoices->enlargeLiberally(currentRow, false);
markovianChoices->set(currentRow);
}
if (!generator->isDeterministicModel()) {
transitionMatrixBuilder.newRowGroup(currentRow);
}
@ -262,6 +277,12 @@ namespace storm {
choiceLabels.get().push_back(choice.getChoiceLabels());
}
// If we keep track of the Markovian choices, store whether the current one is Markovian.
if( markovianChoices && choice.isMarkovian() ) {
markovianChoices->enlargeLiberally(currentRow, false);
markovianChoices->set(currentRow);
}
// Add the probabilistic behavior to the matrix.
for (auto const& stateProbabilityPair : choice) {
transitionMatrixBuilder.addNextValue(currentRow, stateProbabilityPair.first, stateProbabilityPair.second);
@ -280,6 +301,11 @@ namespace storm {
++currentRowGroup;
}
}
if (markovianChoices) {
// We now know the correct size
markovianChoices->resize(currentRow, false);
}
// If the exploration order was not breadth-first, we need to fix the entries in the matrix according to
// (reversed) mapping of row groups to indices.
@ -304,8 +330,6 @@ namespace storm {
// Fix (c).
this->stateStorage.stateToId.remap([&remapping] (StateType const& state) { return remapping[state]; } );
}
return choiceLabels;
}
template <typename ValueType, typename RewardModelType, typename StateType>
@ -323,7 +347,8 @@ namespace storm {
rewardModelBuilders.emplace_back(generator->getRewardModelInformation(i));
}
modelComponents.choiceLabeling = buildMatrices(transitionMatrixBuilder, rewardModelBuilders);
boost::optional<storm::storage::BitVector> markovianChoices;
buildMatrices(transitionMatrixBuilder, rewardModelBuilders, modelComponents.choiceLabeling, markovianChoices);
modelComponents.transitionMatrix = transitionMatrixBuilder.build();
// Now finalize all reward models.
@ -342,9 +367,44 @@ namespace storm {
}
}
if (generator->getModelType() == storm::generator::ModelType::MA) {
STORM_LOG_ASSERT(markovianChoices, "No information regarding markovian choices available.");
buildMarkovianStates(modelComponents, *markovianChoices);
}
return modelComponents;
}
template <typename ValueType, typename RewardModelType, typename StateType>
void ExplicitModelBuilder<ValueType, RewardModelType, StateType>::buildMarkovianStates(ModelComponents& modelComponents, storm::storage::BitVector const& markovianChoices) const {
modelComponents.markovianStates = storm::storage::BitVector(modelComponents.transitionMatrix.getRowGroupCount(), false);
// Check for each state whether it contains a markovian choice.
for (uint_fast64_t state = 0; state < modelComponents.transitionMatrix.getRowGroupCount(); ++state ) {
uint_fast64_t firstChoice = modelComponents.transitionMatrix.getRowGroupIndices()[state];
uint_fast64_t markovianChoice = markovianChoices.getNextSetIndex(firstChoice);
if (markovianChoice < modelComponents.transitionMatrix.getRowGroupIndices()[state+1]) {
// Found a markovian choice. Assert that there is not a second one.
STORM_LOG_THROW(markovianChoices.getNextSetIndex(markovianChoice+1) >= modelComponents.transitionMatrix.getRowGroupIndices()[state+1], storm::exceptions::WrongFormatException, "Multiple Markovian choices defined for some state");
modelComponents.markovianStates->set(state);
// Swap the first choice and the found markovian choice (if they are not equal)
if (firstChoice != markovianChoice) {
modelComponents.transitionMatrix.swapRows(firstChoice, markovianChoice);
for (auto& rewardModel : modelComponents.rewardModels) {
if (rewardModel.second.hasStateActionRewards()) {
std::swap(rewardModel.second.getStateActionRewardVector()[firstChoice], rewardModel.second.getStateActionRewardVector()[markovianChoice]);
}
if (rewardModel.second.hasTransitionRewards()) {
rewardModel.second.getTransitionRewardMatrix().swapRows(firstChoice, markovianChoice);
}
}
if (modelComponents.choiceLabeling) {
std::swap(modelComponents.choiceLabeling.get()[firstChoice], modelComponents.choiceLabeling.get()[markovianChoice]);
}
}
}
}
}
template <typename ValueType, typename RewardModelType, typename StateType>
storm::models::sparse::StateLabeling ExplicitModelBuilder<ValueType, RewardModelType, StateType>::buildStateLabeling() {
return generator->label(stateStorage.stateToId, stateStorage.initialStateIndices, stateStorage.deadlockStateIndices);

19
src/builder/ExplicitModelBuilder.h

@ -63,6 +63,9 @@ namespace storm {
// A vector that stores a labeling for each choice.
boost::optional<std::vector<boost::container::flat_set<uint_fast64_t>>> choiceLabeling;
// A vector that stores which states are markovian.
boost::optional<storm::storage::BitVector> markovianStates;
};
struct Options {
@ -138,10 +141,10 @@ namespace storm {
*
* @param transitionMatrixBuilder The builder of the transition matrix.
* @param rewardModelBuilders The builders for the selected reward models.
* @return A tuple containing a vector with all rows at which the nondeterministic choices of each state begin
* and a vector containing the labels associated with each choice.
* @param choiceLabels is set to a vector containing the labels associated with each choice (is only set if choice labels are requested).
* @param markovianChoices is set to a bit vector storing whether a choice is markovian (is only set if the model type requires this information).
*/
boost::optional<std::vector<boost::container::flat_set<uint_fast64_t>>> buildMatrices(storm::storage::SparseMatrixBuilder<ValueType>& transitionMatrixBuilder, std::vector<RewardModelBuilder<typename RewardModelType::ValueType>>& rewardModelBuilders);
void buildMatrices(storm::storage::SparseMatrixBuilder<ValueType>& transitionMatrixBuilder, std::vector<RewardModelBuilder<typename RewardModelType::ValueType>>& rewardModelBuilders, boost::optional<std::vector<boost::container::flat_set<uint_fast64_t>>>& choiceLabels , boost::optional<storm::storage::BitVector>& markovianChoices);
/*!
* Explores the state space of the given program and returns the components of the model as a result.
@ -150,6 +153,16 @@ namespace storm {
*/
ModelComponents buildModelComponents();
/*!
* Set the markovian states of the given modelComponents,
* makes sure that each state has at most one markovian choice,
* and makes this choice the first one of the corresponding state
*
* @param modelComponents The components of the model build so far
* @markovianChoices bit vector storing whether a choice is markovian
*/
void buildMarkovianStates(ModelComponents& modelComponents, storm::storage::BitVector const& markovianChoices) const;
/*!
* Builds the state labeling for the given program.
*

4
src/generator/PrismNextStateGenerator.cpp

@ -366,7 +366,7 @@ namespace storm {
continue;
}
result.push_back(Choice<ValueType>(command.getActionIndex()));
result.push_back(Choice<ValueType>(command.getActionIndex(), command.isMarkovian()));
Choice<ValueType>& choice = result.back();
// Remember the command labels only if we were asked to.
@ -564,4 +564,4 @@ namespace storm {
template class PrismNextStateGenerator<storm::RationalNumber>;
template class PrismNextStateGenerator<storm::RationalFunction>;
}
}
}

51
src/models/sparse/MarkovAutomaton.cpp

@ -1,12 +1,14 @@
#include "src/models/sparse/MarkovAutomaton.h"
#include "src/models/sparse/StandardRewardModel.h"
#include "src/exceptions/InvalidArgumentException.h"
#include "src/utility/constants.h"
#include "src/adapters/CarlAdapter.h"
#include "src/storage/FlexibleSparseMatrix.h"
#include "src/models/sparse/Dtmc.h"
#include "src/models/sparse/StandardRewardModel.h"
#include "src/solver/stateelimination/StateEliminator.h"
#include "src/storage/FlexibleSparseMatrix.h"
#include "src/utility/constants.h"
#include "src/utility/vector.h"
#include "src/utility/macros.h"
#include "src/exceptions/InvalidArgumentException.h"
namespace storm {
namespace models {
@ -33,6 +35,30 @@ namespace storm {
: NondeterministicModel<ValueType, RewardModelType>(storm::models::ModelType::MarkovAutomaton, std::move(transitionMatrix), std::move(stateLabeling), std::move(rewardModels), std::move(optionalChoiceLabeling)), markovianStates(markovianStates), exitRates(std::move(exitRates)), closed(this->checkIsClosed()) {
this->turnRatesToProbabilities();
}
template <typename ValueType, typename RewardModelType>
MarkovAutomaton<ValueType, RewardModelType>::MarkovAutomaton(storm::storage::SparseMatrix<ValueType> const& rateMatrix,
storm::models::sparse::StateLabeling const& stateLabeling,
storm::storage::BitVector const& markovianStates,
std::unordered_map<std::string, RewardModelType> const& rewardModels,
boost::optional<std::vector<LabelSet>> const& optionalChoiceLabeling)
: NondeterministicModel<ValueType, RewardModelType>(storm::models::ModelType::MarkovAutomaton, rateMatrix, stateLabeling, rewardModels, optionalChoiceLabeling), markovianStates(markovianStates) {
turnRatesToProbabilities();
this->closed = checkIsClosed();
}
template <typename ValueType, typename RewardModelType>
MarkovAutomaton<ValueType, RewardModelType>::MarkovAutomaton(storm::storage::SparseMatrix<ValueType>&& rateMatrix,
storm::models::sparse::StateLabeling&& stateLabeling,
storm::storage::BitVector&& markovianStates,
std::unordered_map<std::string, RewardModelType>&& rewardModels,
boost::optional<std::vector<LabelSet>>&& optionalChoiceLabeling)
: NondeterministicModel<ValueType, RewardModelType>(storm::models::ModelType::MarkovAutomaton, rateMatrix, stateLabeling, rewardModels, optionalChoiceLabeling), markovianStates(markovianStates) {
turnRatesToProbabilities();
this->closed = checkIsClosed();
}
template <typename ValueType, typename RewardModelType>
MarkovAutomaton<ValueType, RewardModelType>::MarkovAutomaton(storm::storage::SparseMatrix<ValueType>&& transitionMatrix,
@ -230,9 +256,18 @@ namespace storm {
template <typename ValueType, typename RewardModelType>
void MarkovAutomaton<ValueType, RewardModelType>::turnRatesToProbabilities() {
for (auto state : this->markovianStates) {
for (auto& transition : this->getTransitionMatrix().getRow(this->getTransitionMatrix().getRowGroupIndices()[state])) {
transition.setValue(transition.getValue() / this->exitRates[state]);
this->exitRates.resize(this->getNumberOfStates());
for (uint_fast64_t state = 0; state< this->getNumberOfStates(); ++state) {
uint_fast64_t row = this->getTransitionMatrix().getRowGroupIndices()[state];
if(this->markovianStates.get(state)) {
this->exitRates[state] = this->getTransitionMatrix().getRowSum(row);
for (auto& transition : this->getTransitionMatrix().getRow(row)) {
transition.setValue(transition.getValue() / this->exitRates[state]);
}
++row;
}
for(; row < this->getTransitionMatrix().getRowGroupIndices()[state+1]; ++row) {
STORM_LOG_THROW(storm::utility::isOne(this->getTransitionMatrix().getRowSum(row)), storm::exceptions::InvalidArgumentException, "Transitions of rateMatrix do not sum up to one for some non-Markovian choice.");
}
}
}

41
src/models/sparse/MarkovAutomaton.h

@ -31,7 +31,43 @@ namespace storm {
std::vector<ValueType> const& exitRates,
std::unordered_map<std::string, RewardModelType> const& rewardModels = std::unordered_map<std::string, RewardModelType>(),
boost::optional<std::vector<LabelSet>> const& optionalChoiceLabeling = boost::optional<std::vector<LabelSet>>());
/*!
* Constructs a model from the given data.
*
* For hybrid states (i.e., states with Markovian and probabilistic transitions), it is assumed that the first
* choice corresponds to the markovian transitions.
*
* @param rateMatrix The matrix representing the transitions in the model in terms of rates.
* @param stateLabeling The labeling of the states.
* @param markovianStates A bit vector indicating the Markovian states of the automaton (i.e., states with at least one markovian transition).
* @param rewardModels A mapping of reward model names to reward models.
* @param optionalChoiceLabeling A vector that represents the labels associated with the choices of each state.
*/
MarkovAutomaton(storm::storage::SparseMatrix<ValueType> const& rateMatrix,
storm::models::sparse::StateLabeling const& stateLabeling,
storm::storage::BitVector const& markovianStates,
std::unordered_map<std::string, RewardModelType> const& rewardModels = std::unordered_map<std::string, RewardModelType>(),
boost::optional<std::vector<LabelSet>> const& optionalChoiceLabeling = boost::optional<std::vector<LabelSet>>());
/*!
* Constructs a model from the given data.
*
* For hybrid states (i.e., states with Markovian and probabilistic transitions), it is assumed that the first
* choice corresponds to the markovian transitions.
*
* @param rateMatrix The matrix representing the transitions in the model in terms of rates.
* @param stateLabeling The labeling of the states.
* @param markovianStates A bit vector indicating the Markovian states of the automaton (i.e., states with at least one markovian transition).
* @param rewardModels A mapping of reward model names to reward models.
* @param optionalChoiceLabeling A vector that represents the labels associated with the choices of each state.
*/
MarkovAutomaton(storm::storage::SparseMatrix<ValueType>&& rateMatrix,
storm::models::sparse::StateLabeling&& stateLabeling,
storm::storage::BitVector&& markovianStates,
std::unordered_map<std::string, RewardModelType>&& rewardModels = std::unordered_map<std::string, RewardModelType>(),
boost::optional<std::vector<LabelSet>>&& optionalChoiceLabeling = boost::optional<std::vector<LabelSet>>());
/*!
* Constructs a model by moving the given data.
*
@ -48,7 +84,7 @@ namespace storm {
std::vector<ValueType> const& exitRates,
std::unordered_map<std::string, RewardModelType>&& rewardModels = std::unordered_map<std::string, RewardModelType>(),
boost::optional<std::vector<LabelSet>>&& optionalChoiceLabeling = boost::optional<std::vector<LabelSet>>());
/*!
* Constructs a model by moving the given data.
*
@ -160,7 +196,8 @@ namespace storm {
/*!
* Under the assumption that the Markovian choices of this Markov automaton are expressed in terms of
* rates in the transition matrix, this procedure turns the rates into the corresponding probabilities by
* dividing each entry by the exit rate of the state.
* dividing each entry by the sum of the rates for that choice.
* Also sets the exitRates accordingly and throws an exception if the values for a non-markovian choice do not sum up to one.
*/
void turnRatesToProbabilities();

10
src/storage/BitVector.cpp

@ -252,6 +252,16 @@ namespace storm {
truncateLastBucket();
}
}
void BitVector::enlargeLiberally(uint_fast64_t minimumLength, bool init) {
if(minimumLength > this->size()) {
uint_fast64_t newLength = this->bucketCount() << 6;
while(newLength < minimumLength) {
newLength = newLength << 1;
}
resize(newLength, init);
}
}
BitVector BitVector::operator&(BitVector const& other) const {
STORM_LOG_ASSERT(bitCount == other.bitCount, "Length of the bit vectors does not match.");

11
src/storage/BitVector.h

@ -234,6 +234,17 @@ namespace storm {
*/
void resize(uint_fast64_t newLength, bool init = false);
/*!
* Enlarges the bit vector such that it holds at least the given number of bits (but possibly more).
* This can be used to diminish reallocations when the final size of the bit vector is not known yet.
* The bit vector does not become smaller.
* New bits are initialized to the given value.
*
* @param minimumLength The minimum number of bits that the bit vector should hold.
* @param init The truth value to which to initialize newly created bits.
*/
void enlargeLiberally(uint_fast64_t minimumLength, bool init = false);
/*!
* Performs a logical "and" with the given bit vector. In case the sizes of the bit vectors do not match,
* only the matching portion is considered and the overlapping bits are set to 0.

64
src/storage/SparseMatrix.cpp

@ -658,6 +658,70 @@ namespace storm {
return bv;
}
template<typename ValueType>
void SparseMatrix<ValueType>::swapRows(index_type const& row1, index_type const& row2) {
if(row1==row2) {
return;
}
// Get the index of the row that has more / less entries than the other
index_type largerRow = getRow(row1).getNumberOfEntries() > getRow(row2).getNumberOfEntries() ? row1 : row2;
index_type smallerRow = largerRow == row1 ? row2 : row1;
index_type rowSizeDifference = getRow(largerRow).getNumberOfEntries() - getRow(smallerRow).getNumberOfEntries();
// Save contents of larger row
std::vector<MatrixEntry<index_type, value_type>> largerRowContents(getRow(largerRow).begin(), getRow(largerRow).end());
if(largerRow < smallerRow) {
auto writeIt = getRows(largerRow, smallerRow+1).begin();
// write smaller row in its new position
for(auto& smallerRowEntry : getRow(smallerRow)) {
*writeIt = std::move(smallerRowEntry);
++writeIt;
}
if(!storm::utility::isZero(rowSizeDifference)) {
// write the intermediate rows into their correct position
for(auto& intermediateRowEntry : getRows(largerRow+1, smallerRow)) {
*writeIt = std::move(intermediateRowEntry);
++writeIt;
}
}
// write the larger row
for(auto& largerRowEntry : largerRowContents) {
*writeIt = std::move(largerRowEntry);
++writeIt;
}
STORM_LOG_ASSERT(writeIt == getRow(smallerRow).end(), "Unexpected position of write iterator");
//Update row indications
for(index_type row = largerRow +1; row <= smallerRow; ++row) {
rowIndications[row] -= rowSizeDifference;
}
} else {
auto writeIt = getRows(smallerRow, largerRow+1).end() -1;
// write smaller row in its new position
for(auto smallerRowEntryIt = getRow(smallerRow).end() -1; smallerRowEntryIt != getRow(smallerRow).begin()-1; --smallerRowEntryIt) {
*writeIt = std::move(*smallerRowEntryIt);
--writeIt;
}
if(!storm::utility::isZero(rowSizeDifference)) {
// write the intermediate rows into their correct position
for(auto intermediateRowEntryIt = getRows(smallerRow+1, largerRow).end() -1; intermediateRowEntryIt != getRows(smallerRow+1, largerRow).begin()-1; --intermediateRowEntryIt) {
*writeIt = std::move(*intermediateRowEntryIt);
--writeIt;
}
}
// write the larger row
for(auto largerRowEntryIt = largerRowContents.rbegin(); largerRowEntryIt != largerRowContents.rend(); ++largerRowEntryIt) {
*writeIt = std::move(*largerRowEntryIt);
--writeIt;
}
STORM_LOG_ASSERT(writeIt == getRow(smallerRow).begin()-1, "Unexpected position of write iterator");
//Update row indications
for(index_type row = smallerRow +1; row <= largerRow; ++row) {
rowIndications[row] += rowSizeDifference;
}
}
}
template<typename ValueType>
ValueType SparseMatrix<ValueType>::getConstrainedRowSum(index_type row, storm::storage::BitVector const& constraint) const {
ValueType result = storm::utility::zero<ValueType>();

8
src/storage/SparseMatrix.h

@ -655,11 +655,19 @@ namespace storm {
* @return True if the rows have identical entries.
*/
bool compareRows(index_type i1, index_type i2) const;
/*!
* Finds duplicate rows in a rowgroup.
*/
BitVector duplicateRowsInRowgroups() const;
/**
* Swaps the two rows.
* @param row1 Index of first row
* @param row2 Index of second row
*/
void swapRows(index_type const& row1, index_type const& row2);
/*!
* Selects exactly one row from each row group of this matrix and returns the resulting matrix.
*

2
src/utility/ModelInstantiator.h

@ -82,7 +82,7 @@ namespace storm {
auto markovianStatesCopy = parametricModel.getMarkovianStates();
auto choiceLabelingCopy = parametricModel.getOptionalChoiceLabeling();
std::vector<ConstantType> exitRates(parametricModel.getExitRates().size(), storm::utility::one<ConstantType>());
this->instantiatedModel = std::make_shared<ConstantSparseModelType>(buildDummyMatrix(parametricModel.getTransitionMatrix()), std::move(stateLabelingCopy), std::move(markovianStatesCopy), std::move(exitRates), buildDummyRewardModels(parametricModel.getRewardModels()), std::move(choiceLabelingCopy));
this->instantiatedModel = std::make_shared<ConstantSparseModelType>(buildDummyMatrix(parametricModel.getTransitionMatrix()), std::move(stateLabelingCopy), std::move(markovianStatesCopy), std::move(exitRates), true, buildDummyRewardModels(parametricModel.getRewardModels()), std::move(choiceLabelingCopy));
initializeVectorMapping(this->instantiatedModel->getExitRates(), this->functions, this->vectorMapping, parametricModel.getExitRates());
}

26
test/functional/builder/ExplicitPrismModelBuilderTest.cpp

@ -1,6 +1,7 @@
#include "gtest/gtest.h"
#include "storm-config.h"
#include "src/models/sparse/StandardRewardModel.h"
#include "src/models/sparse/MarkovAutomaton.h"
#include "src/settings/SettingMemento.h"
#include "src/parser/PrismParser.h"
#include "src/builder/ExplicitModelBuilder.h"
@ -99,6 +100,31 @@ TEST(ExplicitPrismModelBuilderTest, Mdp) {
EXPECT_EQ(59ul, model->getNumberOfTransitions());
}
TEST(ExplicitPrismModelBuilderTest, Ma) {
storm::prism::Program program = storm::parser::PrismParser::parse(STORM_CPP_TESTS_BASE_PATH "/functional/builder/simple.ma");
std::shared_ptr<storm::models::sparse::Model<double>> model = storm::builder::ExplicitModelBuilder<double>(program).build();
EXPECT_EQ(5ul, model->getNumberOfStates());
EXPECT_EQ(8ul, model->getNumberOfTransitions());
ASSERT_TRUE(model->isOfType(storm::models::ModelType::MarkovAutomaton));
EXPECT_EQ(4ul, model->as<storm::models::sparse::MarkovAutomaton<double>>()->getMarkovianStates().getNumberOfSetBits());
program = storm::parser::PrismParser::parse(STORM_CPP_TESTS_BASE_PATH "/functional/builder/hybrid_states.ma");
model = storm::builder::ExplicitModelBuilder<double>(program).build();
EXPECT_EQ(5ul, model->getNumberOfStates());
EXPECT_EQ(13ul, model->getNumberOfTransitions());
ASSERT_TRUE(model->isOfType(storm::models::ModelType::MarkovAutomaton));
EXPECT_EQ(5ul, model->as<storm::models::sparse::MarkovAutomaton<double>>()->getMarkovianStates().getNumberOfSetBits());
program = storm::parser::PrismParser::parse(STORM_CPP_TESTS_BASE_PATH "/functional/builder/stream2.ma");
model = storm::builder::ExplicitModelBuilder<double>(program).build();
EXPECT_EQ(12ul, model->getNumberOfStates());
EXPECT_EQ(14ul, model->getNumberOfTransitions());
ASSERT_TRUE(model->isOfType(storm::models::ModelType::MarkovAutomaton));
EXPECT_EQ(7ul, model->as<storm::models::sparse::MarkovAutomaton<double>>()->getMarkovianStates().getNumberOfSetBits());
}
TEST(ExplicitPrismModelBuilderTest, FailComposition) {
storm::prism::Program program = storm::parser::PrismParser::parse(STORM_CPP_TESTS_BASE_PATH "/functional/builder/system_composition.nm");

19
test/functional/builder/hybrid_states.ma

@ -0,0 +1,19 @@
ma
module hybrid_states
s : [0..4];
[] (s=0) -> 0.8 : (s'=0) + 0.2 : (s'=2);
[] (s=0) -> 1 : (s' = 1);
<> (s=0) -> 3 : (s' = 1);
<> (s=1) -> 9 : (s'=0) + 1 : (s'=3);
<> (s=2) -> 12 : (s'=4);
[] (s=3) -> 1 : true;
[] (s=3) -> 1 : (s'=4);
<> (s=3) -> 2 : (s'=3) + 2 : (s'=4);
[] (s=4) -> 1 : true;
<> (s=4) -> 3 : (s'=3);
endmodule

15
test/functional/builder/simple.ma

@ -0,0 +1,15 @@
ma
module simple
s : [0..4];
[alpha] (s=0) -> 1 : (s' = 1);
[beta] (s=0) -> 0.8 : (s'=0) + 0.2 : (s'=2);
<> (s=1) -> 9 : (s'=0) + 1 : (s'=3);
<> (s=2) -> 12 : (s'=4);
<> (s>2) -> 1 : true;
endmodule

50
test/functional/builder/stream2.ma

@ -0,0 +1,50 @@
ma
const int N = 2; // num packages
const double inRate = 4;
const double processingRate = 5;
module streamingclient
s : [0..3]; // current state:
// 0: decide whether to start
// 1: buffering
// 2: running
// 3: done
n : [0..N]; // number of received packages
k : [0..N]; // number of processed packages
[buffer] s=0 & n<N -> 1 : (s'=1);
[start] s=0 & k<n -> 1 : (s'=2) & (k'=k+1);
<> s=1 -> inRate : (n'=n+1) & (s'=0);
<> s=2 & n<N & k<n -> inRate : (n'=n+1) + processingRate : (k'=k+1);
<> s=2 & n<N & k=n -> inRate : (n'=n+1) + processingRate : (s'=0);
<> s=2 & n=N & k<n -> processingRate : (k'=k+1);
<> s=2 & n=N & k=N -> processingRate : (s'=3);
<> s=3 -> 1 : true;
endmodule
// All packages received and buffer empty
label "done" = (s=3);
rewards "buffering"
s=1 : 1;
endrewards
rewards "initialbuffering"
s=1 & k = 0: 1;
endrewards
rewards "intermediatebuffering"
s=1 & k > 0: 1;
endrewards
rewards "numRestarts"
[start] k > 0 : 1;
endrewards
Loading…
Cancel
Save