TimQu
9 years ago
6 changed files with 788 additions and 1 deletions
-
2src/models/sparse/StochasticTwoPlayerGame.cpp
-
161src/utility/ModelInstantiator.cpp
-
158src/utility/ModelInstantiator.h
-
266test/functional/utility/ModelInstantiatorTest.cpp
-
146test/functional/utility/brp16_2.pm
-
56test/functional/utility/coin2_2.pm
@ -0,0 +1,161 @@ |
|||||
|
/*
|
||||
|
* File: ModelInstantiator.cpp |
||||
|
* Author: Tim Quatmann |
||||
|
* |
||||
|
* Created on February 23, 2016 |
||||
|
*/ |
||||
|
|
||||
|
#include "src/utility/ModelInstantiator.h"
|
||||
|
#include "src/models/sparse/StandardRewardModel.h"
|
||||
|
|
||||
|
namespace storm { |
||||
|
namespace utility { |
||||
|
|
||||
|
template<typename ParametricSparseModelType, typename ConstantSparseModelType> |
||||
|
ModelInstantiator<ParametricSparseModelType, ConstantSparseModelType>::ModelInstantiator(ParametricSparseModelType const& parametricModel){ |
||||
|
//Now pre-compute the information for the equation system.
|
||||
|
initializeModelSpecificData(parametricModel); |
||||
|
initializeMatrixMapping(this->instantiatedModel->getTransitionMatrix(), this->functions, this->matrixMapping, parametricModel.getTransitionMatrix()); |
||||
|
|
||||
|
for(auto& rewModel : this->instantiatedModel->getRewardModels()) { |
||||
|
if(rewModel.second.hasStateRewards()){ |
||||
|
initializeVectorMapping(rewModel.second.getStateRewardVector(), this->functions, this->vectorMapping, parametricModel.getRewardModel(rewModel.first).getStateRewardVector()); |
||||
|
} |
||||
|
if(rewModel.second.hasStateActionRewards()){ |
||||
|
initializeVectorMapping(rewModel.second.getStateActionRewardVector(), this->functions, this->vectorMapping, parametricModel.getRewardModel(rewModel.first).getStateActionRewardVector()); |
||||
|
} |
||||
|
if(rewModel.second.hasTransitionRewards()){ |
||||
|
initializeMatrixMapping(rewModel.second.getTransitionRewardMatrix(), this->functions, this->matrixMapping, parametricModel.getRewardModel(rewModel.first).getTransitionRewardMatrix()); |
||||
|
} |
||||
|
} |
||||
|
} |
||||
|
|
||||
|
template<typename ParametricSparseModelType, typename ConstantType> |
||||
|
ModelInstantiator<ParametricSparseModelType, ConstantType>::~ModelInstantiator() { |
||||
|
//Intentionally left empty
|
||||
|
} |
||||
|
|
||||
|
template<typename ParametricSparseModelType, typename ConstantSparseModelType> |
||||
|
storm::storage::SparseMatrix<typename ConstantSparseModelType::ValueType> ModelInstantiator<ParametricSparseModelType, ConstantSparseModelType>::buildDummyMatrix(storm::storage::SparseMatrix<ParametricType> const& parametricMatrix) const{ |
||||
|
storm::storage::SparseMatrixBuilder<ConstantType> matrixBuilder(parametricMatrix.getRowCount(), |
||||
|
parametricMatrix.getColumnCount(), |
||||
|
parametricMatrix.getEntryCount(), |
||||
|
true, // no force dimensions
|
||||
|
true, //Custom row grouping
|
||||
|
parametricMatrix.getRowGroupCount()); |
||||
|
for(std::size_t rowGroup = 0; rowGroup < parametricMatrix.getRowGroupCount(); ++rowGroup){ |
||||
|
matrixBuilder.newRowGroup(parametricMatrix.getRowGroupIndices()[rowGroup]); |
||||
|
for(std::size_t row = parametricMatrix.getRowGroupIndices()[rowGroup]; row < parametricMatrix.getRowGroupIndices()[rowGroup+1]; ++row){ |
||||
|
ConstantType dummyValue = storm::utility::one<ConstantType>(); |
||||
|
for(auto const& paramEntry : parametricMatrix.getRow(row)){ |
||||
|
matrixBuilder.addNextValue(row, paramEntry.getColumn(), dummyValue); |
||||
|
dummyValue = storm::utility::zero<ConstantType>(); |
||||
|
} |
||||
|
} |
||||
|
} |
||||
|
return matrixBuilder.build(); |
||||
|
} |
||||
|
|
||||
|
template<typename ParametricSparseModelType, typename ConstantSparseModelType> |
||||
|
std::unordered_map<std::string, typename ConstantSparseModelType::RewardModelType> ModelInstantiator<ParametricSparseModelType, ConstantSparseModelType>::buildDummyRewardModels(std::unordered_map<std::string, typename ParametricSparseModelType::RewardModelType> const& parametricRewardModel) const { |
||||
|
std::unordered_map<std::string, typename ConstantSparseModelType::RewardModelType> result; |
||||
|
for(auto const& paramRewardModel : parametricRewardModel){ |
||||
|
auto const& rewModel = paramRewardModel.second; |
||||
|
boost::optional<std::vector<ConstantType>> optionalStateRewardVector; |
||||
|
if(rewModel.hasStateRewards()) { |
||||
|
optionalStateRewardVector = std::vector<ConstantType>(rewModel.getStateRewardVector().size()); |
||||
|
} |
||||
|
boost::optional<std::vector<ConstantType>> optionalStateActionRewardVector; |
||||
|
if(rewModel.hasStateActionRewards()) { |
||||
|
optionalStateActionRewardVector = std::vector<ConstantType>(rewModel.getStateActionRewardVector().size()); |
||||
|
} |
||||
|
boost::optional<storm::storage::SparseMatrix<ConstantType>> optionalTransitionRewardMatrix; |
||||
|
if(rewModel.hasTransitionRewards()) { |
||||
|
optionalTransitionRewardMatrix = buildDummyMatrix(rewModel.getTransitionRewardMatrix()); |
||||
|
} |
||||
|
result.insert(std::make_pair(paramRewardModel.first, |
||||
|
storm::models::sparse::StandardRewardModel<ConstantType>(std::move(optionalStateRewardVector), std::move(optionalStateActionRewardVector), std::move(optionalTransitionRewardMatrix)) |
||||
|
)); |
||||
|
} |
||||
|
return result; |
||||
|
} |
||||
|
|
||||
|
template<typename ParametricSparseModelType, typename ConstantSparseModelType> |
||||
|
void ModelInstantiator<ParametricSparseModelType, ConstantSparseModelType>::initializeMatrixMapping(storm::storage::SparseMatrix<ConstantType>& constantMatrix, |
||||
|
std::unordered_map<ParametricType, ConstantType>& functions, |
||||
|
std::vector<std::pair<typename storm::storage::SparseMatrix<ConstantType>::iterator, ConstantType*>>& mapping, |
||||
|
storm::storage::SparseMatrix<ParametricType> const& parametricMatrix) const{ |
||||
|
ConstantType dummyValue = storm::utility::one<ConstantType>(); |
||||
|
auto constantEntryIt = constantMatrix.begin(); |
||||
|
auto parametricEntryIt = parametricMatrix.begin(); |
||||
|
while(parametricEntryIt != parametricMatrix.end()){ |
||||
|
STORM_LOG_ASSERT(parametricEntryIt->getColumn() == constantEntryIt->getColumn(), "Entries of parametric and constant matrix are not at the same position"); |
||||
|
if(storm::utility::isConstant(parametricEntryIt->getValue())){ |
||||
|
//Constant entries can be inserted directly
|
||||
|
constantEntryIt->setValue(storm::utility::convertNumber<ConstantType>(storm::utility::parametric::getConstantPart(parametricEntryIt->getValue()))); |
||||
|
} else { |
||||
|
//insert the new function and store that the current constantMatrix entry needs to be set to the value of this function
|
||||
|
auto functionsIt = functions.insert(std::make_pair(parametricEntryIt->getValue(), dummyValue)).first; |
||||
|
mapping.emplace_back(std::make_pair(constantEntryIt, &(functionsIt->second))); |
||||
|
//Note that references to elements of an unordered map remain valid after calling unordered_map::insert.
|
||||
|
} |
||||
|
++constantEntryIt; |
||||
|
++parametricEntryIt; |
||||
|
} |
||||
|
STORM_LOG_ASSERT(constantEntryIt == constantMatrix.end(), "Parametric matrix seems to have more or less entries then the constant matrix"); |
||||
|
//TODO: is this necessary?
|
||||
|
constantMatrix.updateNonzeroEntryCount(); |
||||
|
} |
||||
|
|
||||
|
template<typename ParametricSparseModelType, typename ConstantSparseModelType> |
||||
|
void ModelInstantiator<ParametricSparseModelType, ConstantSparseModelType>::initializeVectorMapping(std::vector<ConstantType>& constantVector, |
||||
|
std::unordered_map<ParametricType, ConstantType>& functions, |
||||
|
std::vector<std::pair<typename std::vector<ConstantType>::iterator, ConstantType*>>& mapping, |
||||
|
std::vector<ParametricType> const& parametricVector) const{ |
||||
|
ConstantType dummyValue = storm::utility::one<ConstantType>(); |
||||
|
auto constantEntryIt = constantVector.begin(); |
||||
|
auto parametricEntryIt = parametricVector.begin(); |
||||
|
while(parametricEntryIt != parametricVector.end()){ |
||||
|
if(storm::utility::isConstant(storm::utility::simplify(*parametricEntryIt))){ |
||||
|
//Constant entries can be inserted directly
|
||||
|
*constantEntryIt = storm::utility::convertNumber<ConstantType>(storm::utility::parametric::getConstantPart(*parametricEntryIt)); |
||||
|
} else { |
||||
|
//insert the new function and store that the current constantVector entry needs to be set to the value of this function
|
||||
|
auto functionsIt = functions.insert(std::make_pair(*parametricEntryIt, dummyValue)).first; |
||||
|
mapping.emplace_back(std::make_pair(constantEntryIt, &(functionsIt->second))); |
||||
|
//Note that references to elements of an unordered map remain valid after calling unordered_map::insert.
|
||||
|
} |
||||
|
++constantEntryIt; |
||||
|
++parametricEntryIt; |
||||
|
} |
||||
|
STORM_LOG_ASSERT(constantEntryIt == constantVector.end(), "Parametric vector seems to have more or less entries then the constant vector"); |
||||
|
} |
||||
|
|
||||
|
template<typename ParametricSparseModelType, typename ConstantSparseModelType> |
||||
|
ConstantSparseModelType const& ModelInstantiator<ParametricSparseModelType, ConstantSparseModelType>::instantiate(std::map<VariableType, CoefficientType>const& valuation){ |
||||
|
//Write results into the placeholders
|
||||
|
for(auto& functionResult : this->functions){ |
||||
|
functionResult.second=storm::utility::convertNumber<ConstantType>( |
||||
|
storm::utility::parametric::evaluate(functionResult.first, valuation)); |
||||
|
} |
||||
|
|
||||
|
//Write the instantiated values to the matrices and vectors according to the stored mappings
|
||||
|
for(auto& entryValuePair : this->matrixMapping){ |
||||
|
entryValuePair.first->setValue(*(entryValuePair.second)); |
||||
|
} |
||||
|
for(auto& entryValuePair : this->vectorMapping){ |
||||
|
*(entryValuePair.first)=*(entryValuePair.second); |
||||
|
} |
||||
|
|
||||
|
return *this->instantiatedModel; |
||||
|
} |
||||
|
|
||||
|
#ifdef STORM_HAVE_CARL
|
||||
|
template class ModelInstantiator<storm::models::sparse::Dtmc<storm::RationalFunction>, storm::models::sparse::Dtmc<double>>; |
||||
|
template class ModelInstantiator<storm::models::sparse::Mdp<storm::RationalFunction>, storm::models::sparse::Mdp<double>>; |
||||
|
template class ModelInstantiator<storm::models::sparse::Ctmc<storm::RationalFunction>, storm::models::sparse::Ctmc<double>>; |
||||
|
template class ModelInstantiator<storm::models::sparse::MarkovAutomaton<storm::RationalFunction>, storm::models::sparse::MarkovAutomaton<double>>; |
||||
|
template class ModelInstantiator<storm::models::sparse::StochasticTwoPlayerGame<storm::RationalFunction>, storm::models::sparse::StochasticTwoPlayerGame<double>>; |
||||
|
#endif
|
||||
|
} //namespace utility
|
||||
|
} //namespace storm
|
@ -0,0 +1,158 @@ |
|||||
|
/* |
||||
|
* File: ModelInstantiator.h |
||||
|
* Author: Tim Quatmann |
||||
|
* |
||||
|
* Created on February 23, 2016 |
||||
|
*/ |
||||
|
|
||||
|
#ifndef STORM_UTILITY_MODELINSTANTIATOR_H |
||||
|
#define STORM_UTILITY_MODELINSTANTIATOR_H |
||||
|
|
||||
|
#include <unordered_map> |
||||
|
#include <memory> |
||||
|
#include <type_traits> |
||||
|
|
||||
|
#include "src/models/sparse/Dtmc.h" |
||||
|
#include "src/models/sparse/Mdp.h" |
||||
|
#include "src/models/sparse/Ctmc.h" |
||||
|
#include "src/models/sparse/MarkovAutomaton.h" |
||||
|
#include "src/models/sparse/StochasticTwoPlayerGame.h" |
||||
|
#include "src/utility/parametric.h" |
||||
|
#include "src/utility/constants.h" |
||||
|
|
||||
|
namespace storm { |
||||
|
namespace utility{ |
||||
|
|
||||
|
|
||||
|
/*! |
||||
|
* This class allows efficient instantiation of the given parametric model. |
||||
|
* The key to efficiency is to evaluate every distinct transition- (or reward-) function only once |
||||
|
* instead of evaluating the same function for each occurrence in the model. |
||||
|
*/ |
||||
|
template<typename ParametricSparseModelType, typename ConstantSparseModelType> |
||||
|
class ModelInstantiator { |
||||
|
public: |
||||
|
typedef typename ParametricSparseModelType::ValueType ParametricType; |
||||
|
typedef typename storm::utility::parametric::VariableType<ParametricType>::type VariableType; |
||||
|
typedef typename storm::utility::parametric::CoefficientType<ParametricType>::type CoefficientType; |
||||
|
typedef typename ConstantSparseModelType::ValueType ConstantType; |
||||
|
|
||||
|
/*! |
||||
|
* Constructs a ModelInstantiator |
||||
|
* @param parametricModel The model that is to be instantiated |
||||
|
*/ |
||||
|
ModelInstantiator(ParametricSparseModelType const& parametricModel); |
||||
|
|
||||
|
/*! |
||||
|
* Destructs the ModelInstantiator |
||||
|
*/ |
||||
|
virtual ~ModelInstantiator(); |
||||
|
|
||||
|
/*! |
||||
|
* Evaluates the occurring parametric functions and retrieves the instantiated model |
||||
|
* @param valuation Maps each occurring variables to the value with which it should be substituted |
||||
|
* @return The instantiated model |
||||
|
*/ |
||||
|
ConstantSparseModelType const& instantiate(std::map<VariableType, CoefficientType>const& valuation); |
||||
|
|
||||
|
|
||||
|
private: |
||||
|
/*! |
||||
|
* Initializes the instantiatedModel with dummy data by considering the model-specific ingredients. |
||||
|
* Also initializes other model-specific data, e.g., the exitRate vector of a markov automaton |
||||
|
*/ |
||||
|
template<typename PMT = ParametricSparseModelType> |
||||
|
typename std::enable_if< |
||||
|
std::is_same<PMT,storm::models::sparse::Dtmc<typename ParametricSparseModelType::ValueType>>::value || |
||||
|
std::is_same<PMT,storm::models::sparse::Mdp<typename ParametricSparseModelType::ValueType>>::value || |
||||
|
std::is_same<PMT,storm::models::sparse::Ctmc<typename ParametricSparseModelType::ValueType>>::value |
||||
|
>::type |
||||
|
initializeModelSpecificData(PMT const& parametricModel) { |
||||
|
auto stateLabelingCopy = parametricModel.getStateLabeling(); |
||||
|
auto choiceLabelingCopy = parametricModel.getOptionalChoiceLabeling(); |
||||
|
this->instantiatedModel = std::make_shared<ConstantSparseModelType>(buildDummyMatrix(parametricModel.getTransitionMatrix()), std::move(stateLabelingCopy), buildDummyRewardModels(parametricModel.getRewardModels()), std::move(choiceLabelingCopy)); |
||||
|
} |
||||
|
|
||||
|
template<typename PMT = ParametricSparseModelType> |
||||
|
typename std::enable_if< |
||||
|
std::is_same<PMT,storm::models::sparse::MarkovAutomaton<typename ParametricSparseModelType::ValueType>>::value |
||||
|
>::type |
||||
|
initializeModelSpecificData(PMT const& parametricModel) { |
||||
|
auto stateLabelingCopy = parametricModel.getStateLabeling(); |
||||
|
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)); |
||||
|
|
||||
|
initializeVectorMapping(this->instantiatedModel->getExitRates(), this->functions, this->vectorMapping, parametricModel.getExitRates()); |
||||
|
} |
||||
|
|
||||
|
template<typename PMT = ParametricSparseModelType> |
||||
|
typename std::enable_if< |
||||
|
std::is_same<PMT,storm::models::sparse::StochasticTwoPlayerGame<typename ParametricSparseModelType::ValueType>>::value |
||||
|
>::type |
||||
|
initializeModelSpecificData(PMT const& parametricModel) { |
||||
|
auto player1MatrixCopy = parametricModel.getPlayer1Matrix(); |
||||
|
auto stateLabelingCopy = parametricModel.getStateLabeling(); |
||||
|
boost::optional<std::vector<typename storm::models::sparse::LabelSet>> player1ChoiceLabeling, player2ChoiceLabeling; |
||||
|
if(parametricModel.hasPlayer1ChoiceLabeling()) player1ChoiceLabeling = parametricModel.getPlayer1ChoiceLabeling(); |
||||
|
if(parametricModel.hasPlayer2ChoiceLabeling()) player2ChoiceLabeling = parametricModel.getPlayer2ChoiceLabeling(); |
||||
|
this->instantiatedModel = std::make_shared<ConstantSparseModelType>(std::move(player1MatrixCopy), buildDummyMatrix(parametricModel.getTransitionMatrix()), std::move(stateLabelingCopy), buildDummyRewardModels(parametricModel.getRewardModels()), std::move(player1ChoiceLabeling), std::move(player2ChoiceLabeling)); |
||||
|
} |
||||
|
|
||||
|
/*! |
||||
|
* Creates a matrix that has entries at the same position as the given matrix. |
||||
|
* The returned matrix is a stochastic matrix, i.e., the rows sum up to one. |
||||
|
*/ |
||||
|
storm::storage::SparseMatrix<ConstantType> buildDummyMatrix(storm::storage::SparseMatrix<ParametricType> const& parametricMatrix) const; |
||||
|
|
||||
|
/*! |
||||
|
* Creates a copy of the given reward models with the same names and with state(action)rewards / transitionrewards having the same entry-count and entry-positions. |
||||
|
*/ |
||||
|
std::unordered_map<std::string, typename ConstantSparseModelType::RewardModelType> buildDummyRewardModels(std::unordered_map<std::string, typename ParametricSparseModelType::RewardModelType> const& parametricRewardModel) const; |
||||
|
|
||||
|
/*! |
||||
|
* Connects the occurring functions with the corresponding matrix entries |
||||
|
* |
||||
|
* @note constantMatrix and parametricMatrix should have entries at the same positions |
||||
|
* |
||||
|
* @param constantMatrix The matrix to which the evaluation results are written |
||||
|
* @param functions Occurring functions are inserted in this map |
||||
|
* @param mapping The connections of functions to matrix entries are push_backed into this |
||||
|
* @param parametricMatrix the source matrix with the functions to consider. |
||||
|
*/ |
||||
|
void initializeMatrixMapping(storm::storage::SparseMatrix<ConstantType>& constantMatrix, |
||||
|
std::unordered_map<ParametricType, ConstantType>& functions, |
||||
|
std::vector<std::pair<typename storm::storage::SparseMatrix<ConstantType>::iterator, ConstantType*>>& mapping, |
||||
|
storm::storage::SparseMatrix<ParametricType> const& parametricMatrix) const; |
||||
|
|
||||
|
/*! |
||||
|
* Connects the occurring functions with the corresponding vector entries |
||||
|
* |
||||
|
* @note constantVector and parametricVector should have the same size |
||||
|
* |
||||
|
* @param constantVector The vector to which the evaluation results are written |
||||
|
* @param functions Occurring functions with their placeholders are inserted in this map |
||||
|
* @param mapping The connections of functions to vector entries are push_backed into this |
||||
|
* @param parametricVector the source vector with the functions to consider. |
||||
|
*/ |
||||
|
void initializeVectorMapping(std::vector<ConstantType>& constantVector, |
||||
|
std::unordered_map<ParametricType, ConstantType>& functions, |
||||
|
std::vector<std::pair<typename std::vector<ConstantType>::iterator, ConstantType*>>& mapping, |
||||
|
std::vector<ParametricType> const& parametricVector) const; |
||||
|
|
||||
|
/// The resulting model |
||||
|
std::shared_ptr<ConstantSparseModelType> instantiatedModel; |
||||
|
/// the occurring functions together with the corresponding placeholders for their evaluated result |
||||
|
std::unordered_map<ParametricType, ConstantType> functions; |
||||
|
/// Connection of matrix entries with placeholders |
||||
|
std::vector<std::pair<typename storm::storage::SparseMatrix<ConstantType>::iterator, ConstantType*>> matrixMapping; |
||||
|
/// Connection of Vector entries with placeholders |
||||
|
std::vector<std::pair<typename std::vector<ConstantType>::iterator, ConstantType*>> vectorMapping; |
||||
|
|
||||
|
|
||||
|
}; |
||||
|
}//Namespace utility |
||||
|
} //namespace storm |
||||
|
#endif /* STORM_UTILITY_MODELINSTANTIATOR_H */ |
||||
|
|
@ -0,0 +1,266 @@ |
|||||
|
#include "gtest/gtest.h"
|
||||
|
#include "storm-config.h"
|
||||
|
|
||||
|
#ifdef STORM_HAVE_CARL
|
||||
|
|
||||
|
#include "src/adapters/CarlAdapter.h"
|
||||
|
#include<carl/numbers/numbers.h>
|
||||
|
#include<carl/core/VariablePool.h>
|
||||
|
|
||||
|
|
||||
|
#include "src/settings/SettingsManager.h"
|
||||
|
#include "src/settings/modules/GeneralSettings.h"
|
||||
|
|
||||
|
#include "utility/storm.h"
|
||||
|
#include "utility/ModelInstantiator.h"
|
||||
|
#include "src/models/sparse/Model.h"
|
||||
|
#include "src/models/sparse/Dtmc.h"
|
||||
|
#include "src/models/sparse/Mdp.h"
|
||||
|
|
||||
|
TEST(ModelInstantiatorTest, Brp_Prob) { |
||||
|
carl::VariablePool::getInstance().clear(); |
||||
|
|
||||
|
std::string const& programFile = STORM_CPP_TESTS_BASE_PATH "/functional/utility/brp16_2.pm"; |
||||
|
std::string const& formulaAsString = "P=? [F s=5 ]"; |
||||
|
std::string const& constantsAsString = ""; //e.g. pL=0.9,TOACK=0.5
|
||||
|
|
||||
|
// Program and formula
|
||||
|
storm::prism::Program program = storm::parseProgram(programFile); |
||||
|
program.checkValidity(); |
||||
|
std::vector<std::shared_ptr<storm::logic::Formula>> formulas = storm::parseFormulasForProgram(formulaAsString, program); |
||||
|
ASSERT_TRUE(formulas.size()==1); |
||||
|
// Parametric model
|
||||
|
typename storm::builder::ExplicitPrismModelBuilder<storm::RationalFunction>::Options options = storm::builder::ExplicitPrismModelBuilder<storm::RationalFunction>::Options(*formulas[0]); |
||||
|
options.addConstantDefinitionsFromString(program, constantsAsString); |
||||
|
options.preserveFormula(*formulas[0]); |
||||
|
std::shared_ptr<storm::models::sparse::Dtmc<storm::RationalFunction>> dtmc = storm::builder::ExplicitPrismModelBuilder<storm::RationalFunction>().translateProgram(program, options)->as<storm::models::sparse::Dtmc<storm::RationalFunction>>(); |
||||
|
|
||||
|
storm::utility::ModelInstantiator<storm::models::sparse::Dtmc<storm::RationalFunction>, storm::models::sparse::Dtmc<double>> modelInstantiator(*dtmc); |
||||
|
EXPECT_FALSE(dtmc->hasRewardModel()); |
||||
|
|
||||
|
{ |
||||
|
std::map<storm::Variable, storm::RationalNumber> valuation; |
||||
|
storm::Variable const& pL = carl::VariablePool::getInstance().findVariableWithName("pL"); |
||||
|
ASSERT_NE(pL, carl::Variable::NO_VARIABLE); |
||||
|
storm::Variable const& pK = carl::VariablePool::getInstance().findVariableWithName("pK"); |
||||
|
ASSERT_NE(pK, carl::Variable::NO_VARIABLE); |
||||
|
valuation.insert(std::make_pair(pL,carl::rationalize<storm::RationalNumber>(0.8))); |
||||
|
valuation.insert(std::make_pair(pK,carl::rationalize<storm::RationalNumber>(0.9))); |
||||
|
|
||||
|
storm::models::sparse::Dtmc<double> const& instantiated(modelInstantiator.instantiate(valuation)); |
||||
|
|
||||
|
ASSERT_EQ(dtmc->getTransitionMatrix().getRowGroupIndices(), instantiated.getTransitionMatrix().getRowGroupIndices()); |
||||
|
for(std::size_t rowGroup = 0; rowGroup < dtmc->getTransitionMatrix().getRowGroupCount(); ++rowGroup){ |
||||
|
for(std::size_t row = dtmc->getTransitionMatrix().getRowGroupIndices()[rowGroup]; row < dtmc->getTransitionMatrix().getRowGroupIndices()[rowGroup+1]; ++row){ |
||||
|
auto instantiatedEntry = instantiated.getTransitionMatrix().getRow(row).begin(); |
||||
|
for(auto const& paramEntry : dtmc->getTransitionMatrix().getRow(row)){ |
||||
|
EXPECT_EQ(paramEntry.getColumn(), instantiatedEntry->getColumn()); |
||||
|
double evaluatedValue = carl::toDouble(paramEntry.getValue().evaluate(valuation)); |
||||
|
EXPECT_EQ(evaluatedValue, instantiatedEntry->getValue()); |
||||
|
++instantiatedEntry; |
||||
|
} |
||||
|
EXPECT_EQ(instantiated.getTransitionMatrix().getRow(row).end(),instantiatedEntry); |
||||
|
} |
||||
|
} |
||||
|
EXPECT_EQ(dtmc->getStateLabeling(), instantiated.getStateLabeling()); |
||||
|
EXPECT_EQ(dtmc->getOptionalChoiceLabeling(), instantiated.getOptionalChoiceLabeling()); |
||||
|
|
||||
|
storm::modelchecker::SparseDtmcPrctlModelChecker<storm::models::sparse::Dtmc<double>> modelchecker(instantiated); |
||||
|
std::unique_ptr<storm::modelchecker::CheckResult> chkResult = modelchecker.check(*formulas[0]); |
||||
|
storm::modelchecker::ExplicitQuantitativeCheckResult<double>& quantitativeChkResult = chkResult->asExplicitQuantitativeCheckResult<double>(); |
||||
|
EXPECT_NEAR(0.2989278941, quantitativeChkResult[*instantiated.getInitialStates().begin()], storm::settings::generalSettings().getPrecision()); |
||||
|
} |
||||
|
|
||||
|
{ |
||||
|
std::map<storm::Variable, storm::RationalNumber> valuation; |
||||
|
storm::Variable const& pL = carl::VariablePool::getInstance().findVariableWithName("pL"); |
||||
|
ASSERT_NE(pL, carl::Variable::NO_VARIABLE); |
||||
|
storm::Variable const& pK = carl::VariablePool::getInstance().findVariableWithName("pK"); |
||||
|
ASSERT_NE(pK, carl::Variable::NO_VARIABLE); |
||||
|
valuation.insert(std::make_pair(pL,carl::rationalize<storm::RationalNumber>(1))); |
||||
|
valuation.insert(std::make_pair(pK,carl::rationalize<storm::RationalNumber>(1))); |
||||
|
|
||||
|
storm::models::sparse::Dtmc<double> const& instantiated(modelInstantiator.instantiate(valuation)); |
||||
|
|
||||
|
ASSERT_EQ(dtmc->getTransitionMatrix().getRowGroupIndices(), instantiated.getTransitionMatrix().getRowGroupIndices()); |
||||
|
for(std::size_t rowGroup = 0; rowGroup < dtmc->getTransitionMatrix().getRowGroupCount(); ++rowGroup){ |
||||
|
for(std::size_t row = dtmc->getTransitionMatrix().getRowGroupIndices()[rowGroup]; row < dtmc->getTransitionMatrix().getRowGroupIndices()[rowGroup+1]; ++row){ |
||||
|
auto instantiatedEntry = instantiated.getTransitionMatrix().getRow(row).begin(); |
||||
|
for(auto const& paramEntry : dtmc->getTransitionMatrix().getRow(row)){ |
||||
|
EXPECT_EQ(paramEntry.getColumn(), instantiatedEntry->getColumn()); |
||||
|
double evaluatedValue = carl::toDouble(paramEntry.getValue().evaluate(valuation)); |
||||
|
EXPECT_EQ(evaluatedValue, instantiatedEntry->getValue()); |
||||
|
++instantiatedEntry; |
||||
|
} |
||||
|
EXPECT_EQ(instantiated.getTransitionMatrix().getRow(row).end(),instantiatedEntry); |
||||
|
} |
||||
|
} |
||||
|
EXPECT_EQ(dtmc->getStateLabeling(), instantiated.getStateLabeling()); |
||||
|
EXPECT_EQ(dtmc->getOptionalChoiceLabeling(), instantiated.getOptionalChoiceLabeling()); |
||||
|
|
||||
|
storm::modelchecker::SparseDtmcPrctlModelChecker<storm::models::sparse::Dtmc<double>> modelchecker(instantiated); |
||||
|
std::unique_ptr<storm::modelchecker::CheckResult> chkResult = modelchecker.check(*formulas[0]); |
||||
|
storm::modelchecker::ExplicitQuantitativeCheckResult<double>& quantitativeChkResult = chkResult->asExplicitQuantitativeCheckResult<double>(); |
||||
|
EXPECT_EQ(0.0 , quantitativeChkResult[*instantiated.getInitialStates().begin()]); |
||||
|
} |
||||
|
|
||||
|
{ |
||||
|
std::map<storm::Variable, storm::RationalNumber> valuation; |
||||
|
storm::Variable const& pL = carl::VariablePool::getInstance().findVariableWithName("pL"); |
||||
|
ASSERT_NE(pL, carl::Variable::NO_VARIABLE); |
||||
|
storm::Variable const& pK = carl::VariablePool::getInstance().findVariableWithName("pK"); |
||||
|
ASSERT_NE(pK, carl::Variable::NO_VARIABLE); |
||||
|
valuation.insert(std::make_pair(pL,carl::rationalize<storm::RationalNumber>(1))); |
||||
|
valuation.insert(std::make_pair(pK,carl::rationalize<storm::RationalNumber>(0.9))); |
||||
|
|
||||
|
storm::models::sparse::Dtmc<double> const& instantiated(modelInstantiator.instantiate(valuation)); |
||||
|
|
||||
|
ASSERT_EQ(dtmc->getTransitionMatrix().getRowGroupIndices(), instantiated.getTransitionMatrix().getRowGroupIndices()); |
||||
|
for(std::size_t rowGroup = 0; rowGroup < dtmc->getTransitionMatrix().getRowGroupCount(); ++rowGroup){ |
||||
|
for(std::size_t row = dtmc->getTransitionMatrix().getRowGroupIndices()[rowGroup]; row < dtmc->getTransitionMatrix().getRowGroupIndices()[rowGroup+1]; ++row){ |
||||
|
auto instantiatedEntry = instantiated.getTransitionMatrix().getRow(row).begin(); |
||||
|
for(auto const& paramEntry : dtmc->getTransitionMatrix().getRow(row)){ |
||||
|
EXPECT_EQ(paramEntry.getColumn(), instantiatedEntry->getColumn()); |
||||
|
double evaluatedValue = carl::toDouble(paramEntry.getValue().evaluate(valuation)); |
||||
|
EXPECT_EQ(evaluatedValue, instantiatedEntry->getValue()); |
||||
|
++instantiatedEntry; |
||||
|
} |
||||
|
EXPECT_EQ(instantiated.getTransitionMatrix().getRow(row).end(),instantiatedEntry); |
||||
|
} |
||||
|
} |
||||
|
EXPECT_EQ(dtmc->getStateLabeling(), instantiated.getStateLabeling()); |
||||
|
EXPECT_EQ(dtmc->getOptionalChoiceLabeling(), instantiated.getOptionalChoiceLabeling()); |
||||
|
|
||||
|
storm::modelchecker::SparseDtmcPrctlModelChecker<storm::models::sparse::Dtmc<double>> modelchecker(instantiated); |
||||
|
std::unique_ptr<storm::modelchecker::CheckResult> chkResult = modelchecker.check(*formulas[0]); |
||||
|
storm::modelchecker::ExplicitQuantitativeCheckResult<double>& quantitativeChkResult = chkResult->asExplicitQuantitativeCheckResult<double>(); |
||||
|
EXPECT_NEAR(0.01588055832, quantitativeChkResult[*instantiated.getInitialStates().begin()], storm::settings::generalSettings().getPrecision()); |
||||
|
} |
||||
|
} |
||||
|
|
||||
|
TEST(ModelInstantiatorTest, Brp_Rew) { |
||||
|
carl::VariablePool::getInstance().clear(); |
||||
|
|
||||
|
std::string const& programFile = STORM_CPP_TESTS_BASE_PATH "/functional/utility/brp16_2.pm"; |
||||
|
std::string const& formulaAsString = "R=? [F ((s=5) | (s=0&srep=3)) ]"; |
||||
|
std::string const& constantsAsString = ""; //e.g. pL=0.9,TOACK=0.5
|
||||
|
|
||||
|
// Program and formula
|
||||
|
storm::prism::Program program = storm::parseProgram(programFile); |
||||
|
program.checkValidity(); |
||||
|
std::vector<std::shared_ptr<storm::logic::Formula>> formulas = storm::parseFormulasForProgram(formulaAsString, program); |
||||
|
ASSERT_TRUE(formulas.size()==1); |
||||
|
// Parametric model
|
||||
|
typename storm::builder::ExplicitPrismModelBuilder<storm::RationalFunction>::Options options = storm::builder::ExplicitPrismModelBuilder<storm::RationalFunction>::Options(*formulas[0]); |
||||
|
options.addConstantDefinitionsFromString(program, constantsAsString); |
||||
|
options.preserveFormula(*formulas[0]); |
||||
|
std::shared_ptr<storm::models::sparse::Dtmc<storm::RationalFunction>> dtmc = storm::builder::ExplicitPrismModelBuilder<storm::RationalFunction>().translateProgram(program, options)->as<storm::models::sparse::Dtmc<storm::RationalFunction>>(); |
||||
|
|
||||
|
storm::utility::ModelInstantiator<storm::models::sparse::Dtmc<storm::RationalFunction>, storm::models::sparse::Dtmc<double>> modelInstantiator(*dtmc); |
||||
|
|
||||
|
{ |
||||
|
std::map<storm::Variable, storm::RationalNumber> valuation; |
||||
|
storm::Variable const& pL = carl::VariablePool::getInstance().findVariableWithName("pL"); |
||||
|
ASSERT_NE(pL, carl::Variable::NO_VARIABLE); |
||||
|
storm::Variable const& pK = carl::VariablePool::getInstance().findVariableWithName("pK"); |
||||
|
ASSERT_NE(pK, carl::Variable::NO_VARIABLE); |
||||
|
storm::Variable const& TOMsg = carl::VariablePool::getInstance().findVariableWithName("TOMsg"); |
||||
|
ASSERT_NE(pK, carl::Variable::NO_VARIABLE); |
||||
|
storm::Variable const& TOAck = carl::VariablePool::getInstance().findVariableWithName("TOAck"); |
||||
|
ASSERT_NE(pK, carl::Variable::NO_VARIABLE); |
||||
|
valuation.insert(std::make_pair(pL,carl::rationalize<storm::RationalNumber>(0.9))); |
||||
|
valuation.insert(std::make_pair(pK,carl::rationalize<storm::RationalNumber>(0.3))); |
||||
|
valuation.insert(std::make_pair(TOMsg,carl::rationalize<storm::RationalNumber>(0.3))); |
||||
|
valuation.insert(std::make_pair(TOAck,carl::rationalize<storm::RationalNumber>(0.5))); |
||||
|
|
||||
|
storm::models::sparse::Dtmc<double> const& instantiated(modelInstantiator.instantiate(valuation)); |
||||
|
|
||||
|
ASSERT_EQ(dtmc->getTransitionMatrix().getRowGroupIndices(), instantiated.getTransitionMatrix().getRowGroupIndices()); |
||||
|
for(std::size_t rowGroup = 0; rowGroup < dtmc->getTransitionMatrix().getRowGroupCount(); ++rowGroup){ |
||||
|
for(std::size_t row = dtmc->getTransitionMatrix().getRowGroupIndices()[rowGroup]; row < dtmc->getTransitionMatrix().getRowGroupIndices()[rowGroup+1]; ++row){ |
||||
|
auto instantiatedEntry = instantiated.getTransitionMatrix().getRow(row).begin(); |
||||
|
for(auto const& paramEntry : dtmc->getTransitionMatrix().getRow(row)){ |
||||
|
EXPECT_EQ(paramEntry.getColumn(), instantiatedEntry->getColumn()); |
||||
|
double evaluatedValue = carl::toDouble(paramEntry.getValue().evaluate(valuation)); |
||||
|
EXPECT_EQ(evaluatedValue, instantiatedEntry->getValue()); |
||||
|
++instantiatedEntry; |
||||
|
} |
||||
|
EXPECT_EQ(instantiated.getTransitionMatrix().getRow(row).end(),instantiatedEntry); |
||||
|
} |
||||
|
} |
||||
|
ASSERT_TRUE(instantiated.hasUniqueRewardModel()); |
||||
|
EXPECT_FALSE(instantiated.getUniqueRewardModel()->second.hasStateRewards()); |
||||
|
EXPECT_FALSE(instantiated.getUniqueRewardModel()->second.hasTransitionRewards()); |
||||
|
EXPECT_TRUE(instantiated.getUniqueRewardModel()->second.hasStateActionRewards()); |
||||
|
ASSERT_TRUE(dtmc->getUniqueRewardModel()->second.hasStateActionRewards()); |
||||
|
std::size_t stateActionEntries = dtmc->getUniqueRewardModel()->second.getStateActionRewardVector().size(); |
||||
|
ASSERT_EQ(stateActionEntries, instantiated.getUniqueRewardModel()->second.getStateActionRewardVector().size()); |
||||
|
for(std::size_t i =0; i<stateActionEntries; ++i){ |
||||
|
double evaluatedValue = carl::toDouble(dtmc->getUniqueRewardModel()->second.getStateActionRewardVector()[i].evaluate(valuation)); |
||||
|
EXPECT_EQ(evaluatedValue, instantiated.getUniqueRewardModel()->second.getStateActionRewardVector()[i]); |
||||
|
} |
||||
|
EXPECT_EQ(dtmc->getStateLabeling(), instantiated.getStateLabeling()); |
||||
|
EXPECT_EQ(dtmc->getOptionalChoiceLabeling(), instantiated.getOptionalChoiceLabeling()); |
||||
|
|
||||
|
storm::modelchecker::SparseDtmcPrctlModelChecker<storm::models::sparse::Dtmc<double>> modelchecker(instantiated); |
||||
|
std::unique_ptr<storm::modelchecker::CheckResult> chkResult = modelchecker.check(*formulas[0]); |
||||
|
storm::modelchecker::ExplicitQuantitativeCheckResult<double>& quantitativeChkResult = chkResult->asExplicitQuantitativeCheckResult<double>(); |
||||
|
EXPECT_NEAR(1.308324495, quantitativeChkResult[*instantiated.getInitialStates().begin()], storm::settings::generalSettings().getPrecision()); |
||||
|
} |
||||
|
|
||||
|
} |
||||
|
|
||||
|
|
||||
|
TEST(ModelInstantiatorTest, consensus) { |
||||
|
carl::VariablePool::getInstance().clear(); |
||||
|
|
||||
|
std::string const& programFile = STORM_CPP_TESTS_BASE_PATH "/functional/utility/coin2_2.pm"; |
||||
|
std::string const& formulaAsString = "Pmin=? [F \"finished\"&\"all_coins_equal_1\" ]"; |
||||
|
std::string const& constantsAsString = ""; //e.g. pL=0.9,TOACK=0.5
|
||||
|
|
||||
|
// Program and formula
|
||||
|
storm::prism::Program program = storm::parseProgram(programFile); |
||||
|
program.checkValidity(); |
||||
|
std::vector<std::shared_ptr<storm::logic::Formula>> formulas = storm::parseFormulasForProgram(formulaAsString, program); |
||||
|
ASSERT_TRUE(formulas.size()==1); |
||||
|
// Parametric model
|
||||
|
typename storm::builder::ExplicitPrismModelBuilder<storm::RationalFunction>::Options options = storm::builder::ExplicitPrismModelBuilder<storm::RationalFunction>::Options(*formulas[0]); |
||||
|
options.addConstantDefinitionsFromString(program, constantsAsString); |
||||
|
options.preserveFormula(*formulas[0]); |
||||
|
std::shared_ptr<storm::models::sparse::Mdp<storm::RationalFunction>> mdp = storm::builder::ExplicitPrismModelBuilder<storm::RationalFunction>().translateProgram(program, options)->as<storm::models::sparse::Mdp<storm::RationalFunction>>(); |
||||
|
|
||||
|
storm::utility::ModelInstantiator<storm::models::sparse::Mdp<storm::RationalFunction>, storm::models::sparse::Mdp<double>> modelInstantiator(*mdp); |
||||
|
|
||||
|
std::map<storm::Variable, storm::RationalNumber> valuation; |
||||
|
storm::Variable const& p1 = carl::VariablePool::getInstance().findVariableWithName("p1"); |
||||
|
ASSERT_NE(p1, carl::Variable::NO_VARIABLE); |
||||
|
storm::Variable const& p2 = carl::VariablePool::getInstance().findVariableWithName("p2"); |
||||
|
ASSERT_NE(p2, carl::Variable::NO_VARIABLE); |
||||
|
valuation.insert(std::make_pair(p1,carl::rationalize<storm::RationalNumber>(0.51))); |
||||
|
valuation.insert(std::make_pair(p2,carl::rationalize<storm::RationalNumber>(0.49))); |
||||
|
storm::models::sparse::Mdp<double> const& instantiated(modelInstantiator.instantiate(valuation)); |
||||
|
|
||||
|
ASSERT_EQ(mdp->getTransitionMatrix().getRowGroupIndices(), instantiated.getTransitionMatrix().getRowGroupIndices()); |
||||
|
for(std::size_t rowGroup = 0; rowGroup < mdp->getTransitionMatrix().getRowGroupCount(); ++rowGroup){ |
||||
|
for(std::size_t row = mdp->getTransitionMatrix().getRowGroupIndices()[rowGroup]; row < mdp->getTransitionMatrix().getRowGroupIndices()[rowGroup+1]; ++row){ |
||||
|
auto instantiatedEntry = instantiated.getTransitionMatrix().getRow(row).begin(); |
||||
|
for(auto const& paramEntry : mdp->getTransitionMatrix().getRow(row)){ |
||||
|
EXPECT_EQ(paramEntry.getColumn(), instantiatedEntry->getColumn()); |
||||
|
double evaluatedValue = carl::toDouble(paramEntry.getValue().evaluate(valuation)); |
||||
|
EXPECT_EQ(evaluatedValue, instantiatedEntry->getValue()); |
||||
|
++instantiatedEntry; |
||||
|
} |
||||
|
EXPECT_EQ(instantiated.getTransitionMatrix().getRow(row).end(),instantiatedEntry); |
||||
|
} |
||||
|
} |
||||
|
EXPECT_EQ(mdp->getStateLabeling(), instantiated.getStateLabeling()); |
||||
|
EXPECT_EQ(mdp->getOptionalChoiceLabeling(), instantiated.getOptionalChoiceLabeling()); |
||||
|
|
||||
|
storm::modelchecker::SparseMdpPrctlModelChecker<storm::models::sparse::Mdp<double>> modelchecker(instantiated); |
||||
|
std::unique_ptr<storm::modelchecker::CheckResult> chkResult = modelchecker.check(*formulas[0]); |
||||
|
storm::modelchecker::ExplicitQuantitativeCheckResult<double>& quantitativeChkResult = chkResult->asExplicitQuantitativeCheckResult<double>(); |
||||
|
EXPECT_NEAR(0.3526577219, quantitativeChkResult[*instantiated.getInitialStates().begin()], storm::settings::generalSettings().getPrecision()); |
||||
|
|
||||
|
} |
||||
|
|
||||
|
#endif
|
@ -0,0 +1,146 @@ |
|||||
|
// bounded retransmission protocol [D'AJJL01] |
||||
|
// gxn/dxp 23/05/2001 |
||||
|
|
||||
|
dtmc |
||||
|
|
||||
|
// number of chunks |
||||
|
const int N = 16; |
||||
|
// maximum number of retransmissions |
||||
|
const int MAX = 2; |
||||
|
|
||||
|
// reliability of channels |
||||
|
const double pL; |
||||
|
const double pK; |
||||
|
|
||||
|
// timeouts |
||||
|
const double TOMsg; |
||||
|
const double TOAck; |
||||
|
|
||||
|
module sender |
||||
|
|
||||
|
s : [0..6]; |
||||
|
// 0 idle |
||||
|
// 1 next_frame |
||||
|
// 2 wait_ack |
||||
|
// 3 retransmit |
||||
|
// 4 success |
||||
|
// 5 error |
||||
|
// 6 wait sync |
||||
|
srep : [0..3]; |
||||
|
// 0 bottom |
||||
|
// 1 not ok (nok) |
||||
|
// 2 do not know (dk) |
||||
|
// 3 ok (ok) |
||||
|
nrtr : [0..MAX]; |
||||
|
i : [0..N]; |
||||
|
bs : bool; |
||||
|
s_ab : bool; |
||||
|
fs : bool; |
||||
|
ls : bool; |
||||
|
|
||||
|
// idle |
||||
|
[NewFile] (s=0) -> (s'=1) & (i'=1) & (srep'=0); |
||||
|
// next_frame |
||||
|
[aF] (s=1) -> (s'=2) & (fs'=(i=1)) & (ls'=(i=N)) & (bs'=s_ab) & (nrtr'=0); |
||||
|
// wait_ack |
||||
|
[aB] (s=2) -> (s'=4) & (s_ab'=!s_ab); |
||||
|
[TO_Msg] (s=2) -> (s'=3); |
||||
|
[TO_Ack] (s=2) -> (s'=3); |
||||
|
// retransmit |
||||
|
[aF] (s=3) & (nrtr<MAX) -> (s'=2) & (fs'=(i=1)) & (ls'=(i=N)) & (bs'=s_ab) & (nrtr'=nrtr+1); |
||||
|
[] (s=3) & (nrtr=MAX) & (i<N) -> (s'=5) & (srep'=1); |
||||
|
[] (s=3) & (nrtr=MAX) & (i=N) -> (s'=5) & (srep'=2); |
||||
|
// success |
||||
|
[] (s=4) & (i<N) -> (s'=1) & (i'=i+1); |
||||
|
[] (s=4) & (i=N) -> (s'=0) & (srep'=3); |
||||
|
// error |
||||
|
[SyncWait] (s=5) -> (s'=6); |
||||
|
// wait sync |
||||
|
[SyncWait] (s=6) -> (s'=0) & (s_ab'=false); |
||||
|
|
||||
|
endmodule |
||||
|
|
||||
|
module receiver |
||||
|
|
||||
|
r : [0..5]; |
||||
|
// 0 new_file |
||||
|
// 1 fst_safe |
||||
|
// 2 frame_received |
||||
|
// 3 frame_reported |
||||
|
// 4 idle |
||||
|
// 5 resync |
||||
|
rrep : [0..4]; |
||||
|
// 0 bottom |
||||
|
// 1 fst |
||||
|
// 2 inc |
||||
|
// 3 ok |
||||
|
// 4 nok |
||||
|
fr : bool; |
||||
|
lr : bool; |
||||
|
br : bool; |
||||
|
r_ab : bool; |
||||
|
recv : bool; |
||||
|
|
||||
|
|
||||
|
// new_file |
||||
|
[SyncWait] (r=0) -> (r'=0); |
||||
|
[aG] (r=0) -> (r'=1) & (fr'=fs) & (lr'=ls) & (br'=bs) & (recv'=T); |
||||
|
// fst_safe_frame |
||||
|
[] (r=1) -> (r'=2) & (r_ab'=br); |
||||
|
// frame_received |
||||
|
[] (r=2) & (r_ab=br) & (fr=true) & (lr=false) -> (r'=3) & (rrep'=1); |
||||
|
[] (r=2) & (r_ab=br) & (fr=false) & (lr=false) -> (r'=3) & (rrep'=2); |
||||
|
[] (r=2) & (r_ab=br) & (fr=false) & (lr=true) -> (r'=3) & (rrep'=3); |
||||
|
[aA] (r=2) & !(r_ab=br) -> (r'=4); |
||||
|
// frame_reported |
||||
|
[aA] (r=3) -> (r'=4) & (r_ab'=!r_ab); |
||||
|
// idle |
||||
|
[aG] (r=4) -> (r'=2) & (fr'=fs) & (lr'=ls) & (br'=bs) & (recv'=T); |
||||
|
[SyncWait] (r=4) & (ls=true) -> (r'=5); |
||||
|
[SyncWait] (r=4) & (ls=false) -> (r'=5) & (rrep'=4); |
||||
|
// resync |
||||
|
[SyncWait] (r=5) -> (r'=0) & (rrep'=0); |
||||
|
|
||||
|
endmodule |
||||
|
|
||||
|
// prevents more than one file being sent |
||||
|
module tester |
||||
|
|
||||
|
T : bool; |
||||
|
|
||||
|
[NewFile] (T=false) -> (T'=true); |
||||
|
|
||||
|
endmodule |
||||
|
|
||||
|
module channelK |
||||
|
|
||||
|
k : [0..2]; |
||||
|
|
||||
|
// idle |
||||
|
[aF] (k=0) -> pK : (k'=1) + 1-pK : (k'=2); |
||||
|
// sending |
||||
|
[aG] (k=1) -> (k'=0); |
||||
|
// lost |
||||
|
[TO_Msg] (k=2) -> (k'=0); |
||||
|
|
||||
|
endmodule |
||||
|
|
||||
|
module channelL |
||||
|
|
||||
|
l : [0..2]; |
||||
|
|
||||
|
// idle |
||||
|
[aA] (l=0) -> pL : (l'=1) + 1-pL : (l'=2); |
||||
|
// sending |
||||
|
[aB] (l=1) -> (l'=0); |
||||
|
// lost |
||||
|
[TO_Ack] (l=2) -> (l'=0); |
||||
|
|
||||
|
endmodule |
||||
|
|
||||
|
rewards |
||||
|
[TO_Msg] true : TOMsg; |
||||
|
[TO_Ack] true : TOAck; |
||||
|
endrewards |
||||
|
|
||||
|
|
@ -0,0 +1,56 @@ |
|||||
|
//Randomised Consensus Protocol |
||||
|
|
||||
|
mdp |
||||
|
const double p1; // in [0.2 , 0.8] |
||||
|
const double p2; // in [0.2 , 0.8] |
||||
|
|
||||
|
|
||||
|
const int N=2; |
||||
|
const int K=2; |
||||
|
const int range = 2*(K+1)*N; |
||||
|
const int counter_init = (K+1)*N; |
||||
|
const int left = N; |
||||
|
const int right = 2*(K+1)*N - N; |
||||
|
|
||||
|
// shared coin |
||||
|
global counter : [0..range] init counter_init; |
||||
|
|
||||
|
module process1 |
||||
|
|
||||
|
// program counter |
||||
|
pc1 : [0..3]; |
||||
|
// 0 - flip |
||||
|
// 1 - write |
||||
|
// 2 - check |
||||
|
// 3 - finished |
||||
|
|
||||
|
// local coin |
||||
|
coin1 : [0..1]; |
||||
|
|
||||
|
// flip coin |
||||
|
[] (pc1=0) -> p1 : (coin1'=0) & (pc1'=1) + 1 - p1 : (coin1'=1) & (pc1'=1); |
||||
|
// write tails -1 (reset coin to add regularity) |
||||
|
[] (pc1=1) & (coin1=0) & (counter>0) -> (counter'=counter-1) & (pc1'=2) & (coin1'=0); |
||||
|
// write heads +1 (reset coin to add regularity) |
||||
|
[] (pc1=1) & (coin1=1) & (counter<range) -> (counter'=counter+1) & (pc1'=2) & (coin1'=0); |
||||
|
// check |
||||
|
// decide tails |
||||
|
[] (pc1=2) & (counter<=left) -> (pc1'=3) & (coin1'=0); |
||||
|
// decide heads |
||||
|
[] (pc1=2) & (counter>=right) -> (pc1'=3) & (coin1'=1); |
||||
|
// flip again |
||||
|
[] (pc1=2) & (counter>left) & (counter<right) -> (pc1'=0); |
||||
|
// loop (all loop together when done) |
||||
|
[done] (pc1=3) -> (pc1'=3); |
||||
|
|
||||
|
endmodule |
||||
|
|
||||
|
module process2 = process1[pc1=pc2,coin1=coin2,p1=p2] endmodule |
||||
|
label "finished" = pc1=3 &pc2=3 ; |
||||
|
label "all_coins_equal_1" = coin1=1 &coin2=1 ; |
||||
|
rewards "steps" |
||||
|
true : 1; |
||||
|
endrewards |
||||
|
|
||||
|
|
||||
|
|
Write
Preview
Loading…
Cancel
Save
Reference in new issue