Browse Source

The PRISM parser can now parse DTMC models that do not use synchronization.

tempestpy_adaptions
dehnert 12 years ago
parent
commit
a17c99902b
  1. 44
      src/adapters/IntermediateRepresentationAdapter.h
  2. 20
      src/parser/PrismParser.cpp
  3. 17
      src/storm.cpp
  4. 105
      test.input

44
src/adapters/IntermediateRepresentationAdapter.h

@ -9,6 +9,7 @@
#define STORM_IR_INTERMEDIATEREPRESENTATIONADAPTER_H_
#include "src/storage/SparseMatrix.h"
#include "src/utility/Settings.h"
#include <tuple>
#include <unordered_map>
@ -19,6 +20,10 @@
typedef std::pair<std::vector<bool>, std::vector<int_fast64_t>> StateType;
#include "log4cplus/logger.h"
#include "log4cplus/loggingmacros.h"
extern log4cplus::Logger logger;
namespace storm {
namespace adapters {
@ -48,7 +53,7 @@ class IntermediateRepresentationAdapter {
public:
template<class T>
static storm::storage::SparseMatrix<T>* toSparseMatrix(storm::ir::Program const& program) {
LOG4CPLUS_INFO(logger, "Creating sparse matrix for probabilistic program.");
uint_fast64_t numberOfIntegerVariables = 0;
uint_fast64_t numberOfBooleanVariables = 0;
for (uint_fast64_t i = 0; i < program.getNumberOfModules(); ++i) {
@ -108,6 +113,8 @@ public:
StateType* currentState = stateQueue.front();
stateQueue.pop();
bool hasTransition = false;
// Iterate over all modules.
for (uint_fast64_t i = 0; i < program.getNumberOfModules(); ++i) {
storm::ir::Module const& module = program.getModule(i);
@ -118,21 +125,27 @@ public:
// Check if this command is enabled in the current state.
if (command.getGuard()->getValueAsBool(*currentState)) {
hasTransition = true;
std::unordered_map<StateType*, double, StateHash, StateCompare> stateToProbabilityMap;
std::queue<StateType*> statesToDelete;
for (uint_fast64_t k = 0; k < command.getNumberOfUpdates(); ++k) {
storm::ir::Update const& update = command.getUpdate(k);
StateType* newState = new StateType(*currentState);
// std::cout << "took state: " << newState->first << "/" << newState->second << std::endl;
std::map<std::string, storm::ir::Assignment> const& booleanAssignmentMap = update.getBooleanAssignments();
for (auto assignedVariable : booleanAssignmentMap) {
setValue(newState, booleanVariableToIndexMap[assignedVariable.first], assignedVariable.second.getExpression()->getValueAsBool(*currentState));
}
std::map<std::string, storm::ir::Assignment> const& integerAssignmentMap = update.getIntegerAssignments();
for (auto assignedVariable : integerAssignmentMap) {
// std::cout << "evaluting " << assignedVariable.second.getExpression()->toString() << " as " << assignedVariable.second.getExpression()->getValueAsInt(*currentState) << std::endl;
setValue(newState, integerVariableToIndexMap[assignedVariable.first], assignedVariable.second.getExpression()->getValueAsInt(*currentState));
}
// std::cout << "applied: " << update.toString() << std::endl;
// std::cout << "got: " << newState->first << "/" << newState->second << std::endl;
auto probIt = stateToProbabilityMap.find(newState);
if (probIt != stateToProbabilityMap.end()) {
stateToProbabilityMap[newState] += update.getLikelihoodExpression()->getValueAsDouble(*currentState);
@ -144,7 +157,7 @@ public:
auto it = stateToIndexMap.find(newState);
if (it != stateToIndexMap.end()) {
// Delete the state object directly as we have already seen that state.
delete newState;
statesToDelete.push(newState);
} else {
// Add state to the queue of states that are still to be explored.
stateQueue.push(newState);
@ -157,18 +170,31 @@ public:
++nextIndex;
}
}
while (!statesToDelete.empty()) {
delete statesToDelete.front();
statesToDelete.pop();
}
}
}
}
}
std::cout << "Found " << allStates.size() << " reachable states and " << totalNumberOfTransitions << " transitions.";
if (!hasTransition) {
if (storm::settings::instance()->isSet("fix-deadlocks")) {
++totalNumberOfTransitions;
} else {
LOG4CPLUS_ERROR(logger, "Error while creating sparse matrix from probabilistic program: found deadlock state.");
throw storm::exceptions::WrongFileFormatException() << "Error while creating sparse matrix from probabilistic program: found deadlock state.";
}
}
}
storm::storage::SparseMatrix<T>* resultMatrix = new storm::storage::SparseMatrix<T>(allStates.size());
resultMatrix->initialize(totalNumberOfTransitions);
uint_fast64_t currentIndex = 0;
for (StateType* currentState : allStates) {
bool hasTransition = false;
// Iterate over all modules.
for (uint_fast64_t i = 0; i < program.getNumberOfModules(); ++i) {
storm::ir::Module const& module = program.getModule(i);
@ -179,6 +205,7 @@ public:
// Check if this command is enabled in the current state.
if (command.getGuard()->getValueAsBool(*currentState)) {
hasTransition = true;
std::map<uint_fast64_t, double> stateIndexToProbabilityMap;
for (uint_fast64_t k = 0; k < command.getNumberOfUpdates(); ++k) {
storm::ir::Update const& update = command.getUpdate(k);
@ -212,15 +239,22 @@ public:
}
}
}
if (!hasTransition) {
resultMatrix->addNextValue(currentIndex, currentIndex, 1);
}
++currentIndex;
}
resultMatrix->finalize();
LOG4CPLUS_INFO(logger, "Created sparse matrix with " << allStates.size() << " reachable states and " << totalNumberOfTransitions << " transitions.");
// Now free all the elements we allocated.
for (auto element : allStates) {
delete element;
}
return resultMatrix;
}

20
src/parser/PrismParser.cpp

@ -88,7 +88,7 @@ struct PrismParser::PrismGrammar : qi::grammar<Iterator, storm::ir::Program(), q
atomicIntegerExpression.name("integer expression");
integerMultExpression %= atomicIntegerExpression[qi::_val = qi::_1] >> *(qi::lit("*") >> atomicIntegerExpression)[qi::_val = phoenix::construct<std::shared_ptr<storm::ir::expressions::BaseExpression>>(phoenix::new_<storm::ir::expressions::BinaryNumericalFunctionExpression>(storm::ir::expressions::BaseExpression::int_, qi::_val, qi::_1, storm::ir::expressions::BinaryNumericalFunctionExpression::TIMES))];
integerMultExpression.name("integer expression");
integerPlusExpression = integerMultExpression[qi::_val = qi::_1] >> *(qi::lit("+") >> integerMultExpression)[qi::_val = phoenix::construct<std::shared_ptr<storm::ir::expressions::BaseExpression>>(phoenix::new_<storm::ir::expressions::BinaryNumericalFunctionExpression>(storm::ir::expressions::BaseExpression::int_, qi::_val, qi::_1, storm::ir::expressions::BinaryNumericalFunctionExpression::PLUS))];
integerPlusExpression = integerMultExpression[qi::_val = qi::_1] >> *((qi::lit("+")[qi::_a = true] | qi::lit("-")[qi::_a = false]) >> integerMultExpression)[phoenix::if_(qi::_a) [qi::_val = phoenix::construct<std::shared_ptr<storm::ir::expressions::BaseExpression>>(phoenix::new_<storm::ir::expressions::BinaryNumericalFunctionExpression>(storm::ir::expressions::BaseExpression::int_, qi::_val, qi::_1, storm::ir::expressions::BinaryNumericalFunctionExpression::PLUS)) ] .else_ [qi::_val = phoenix::construct<std::shared_ptr<storm::ir::expressions::BaseExpression>>(phoenix::new_<storm::ir::expressions::BinaryNumericalFunctionExpression>(storm::ir::expressions::BaseExpression::int_, qi::_val, qi::_1, storm::ir::expressions::BinaryNumericalFunctionExpression::MINUS))]];
integerPlusExpression.name("integer expression");
integerExpression %= integerPlusExpression;
integerExpression.name("integer expression");
@ -98,7 +98,7 @@ struct PrismParser::PrismGrammar : qi::grammar<Iterator, storm::ir::Program(), q
constantAtomicIntegerExpression.name("constant integer expression");
constantIntegerMultExpression %= constantAtomicIntegerExpression[qi::_val = qi::_1] >> *(qi::lit("*") >> constantAtomicIntegerExpression)[qi::_val = phoenix::construct<std::shared_ptr<storm::ir::expressions::BaseExpression>>(phoenix::new_<storm::ir::expressions::BinaryNumericalFunctionExpression>(storm::ir::expressions::BaseExpression::int_, qi::_val, qi::_1, storm::ir::expressions::BinaryNumericalFunctionExpression::TIMES))];
constantIntegerMultExpression.name("constant integer expression");
constantIntegerPlusExpression = constantIntegerMultExpression[qi::_val = qi::_1] >> *(qi::lit("+") >> constantIntegerMultExpression)[qi::_val = phoenix::construct<std::shared_ptr<storm::ir::expressions::BaseExpression>>(phoenix::new_<storm::ir::expressions::BinaryNumericalFunctionExpression>(storm::ir::expressions::BaseExpression::int_, qi::_val, qi::_1, storm::ir::expressions::BinaryNumericalFunctionExpression::PLUS))];
constantIntegerPlusExpression = constantIntegerMultExpression[qi::_val = qi::_1] >> *((qi::lit("+")[qi::_a = true] | qi::lit("-")[qi::_a = false]) >> constantIntegerMultExpression)[phoenix::if_(qi::_a) [qi::_val = phoenix::construct<std::shared_ptr<storm::ir::expressions::BaseExpression>>(phoenix::new_<storm::ir::expressions::BinaryNumericalFunctionExpression>(storm::ir::expressions::BaseExpression::int_, qi::_val, qi::_1, storm::ir::expressions::BinaryNumericalFunctionExpression::PLUS)) ] .else_ [qi::_val = phoenix::construct<std::shared_ptr<storm::ir::expressions::BaseExpression>>(phoenix::new_<storm::ir::expressions::BinaryNumericalFunctionExpression>(storm::ir::expressions::BaseExpression::int_, qi::_val, qi::_1, storm::ir::expressions::BinaryNumericalFunctionExpression::MINUS))]];
constantIntegerPlusExpression.name("constant integer expression");
constantIntegerExpression %= constantIntegerPlusExpression;
constantIntegerExpression.name("constant integer expression");
@ -106,9 +106,9 @@ struct PrismParser::PrismGrammar : qi::grammar<Iterator, storm::ir::Program(), q
// This block defines all expressions of type double that are by syntax constant. That is, they are evaluable given the values for all constants.
constantAtomicDoubleExpression %= (qi::lit("(") >> constantDoubleExpression >> qi::lit(")") | doubleConstantExpression);
constantAtomicDoubleExpression.name("constant double expression");
constantDoubleMultExpression %= constantAtomicDoubleExpression[qi::_val = qi::_1] >> *(qi::lit("*") >> constantAtomicDoubleExpression)[qi::_val = phoenix::construct<std::shared_ptr<storm::ir::expressions::BaseExpression>>(phoenix::new_<storm::ir::expressions::BinaryNumericalFunctionExpression>(storm::ir::expressions::BaseExpression::double_, qi::_val, qi::_1, storm::ir::expressions::BinaryNumericalFunctionExpression::TIMES))];
constantDoubleMultExpression %= constantAtomicDoubleExpression[qi::_val = qi::_1] >> *((qi::lit("*")[qi::_a = true] | qi::lit("/")[qi::_a = false]) >> constantAtomicDoubleExpression)[phoenix::if_(qi::_a) [qi::_val = phoenix::construct<std::shared_ptr<storm::ir::expressions::BaseExpression>>(phoenix::new_<storm::ir::expressions::BinaryNumericalFunctionExpression>(storm::ir::expressions::BaseExpression::double_, qi::_val, qi::_1, storm::ir::expressions::BinaryNumericalFunctionExpression::TIMES)) ] .else_ [qi::_val = phoenix::construct<std::shared_ptr<storm::ir::expressions::BaseExpression>>(phoenix::new_<storm::ir::expressions::BinaryNumericalFunctionExpression>(storm::ir::expressions::BaseExpression::double_, qi::_val, qi::_1, storm::ir::expressions::BinaryNumericalFunctionExpression::DIVIDE))]];
constantDoubleMultExpression.name("constant double expression");
constantDoublePlusExpression %= constantDoubleMultExpression[qi::_val = qi::_1] >> *(qi::lit("+") >> constantDoubleMultExpression)[qi::_val = phoenix::construct<std::shared_ptr<storm::ir::expressions::BaseExpression>>(phoenix::new_<storm::ir::expressions::BinaryNumericalFunctionExpression>(storm::ir::expressions::BaseExpression::double_, qi::_val, qi::_1, storm::ir::expressions::BinaryNumericalFunctionExpression::PLUS))];
constantDoublePlusExpression %= constantDoubleMultExpression[qi::_val = qi::_1] >> *((qi::lit("+")[qi::_a = true] | qi::lit("-")[qi::_a = false]) >> constantDoubleMultExpression)[phoenix::if_(qi::_a) [qi::_val = phoenix::construct<std::shared_ptr<storm::ir::expressions::BaseExpression>>(phoenix::new_<storm::ir::expressions::BinaryNumericalFunctionExpression>(storm::ir::expressions::BaseExpression::double_, qi::_val, qi::_1, storm::ir::expressions::BinaryNumericalFunctionExpression::PLUS)) ] .else_ [qi::_val = phoenix::construct<std::shared_ptr<storm::ir::expressions::BaseExpression>>(phoenix::new_<storm::ir::expressions::BinaryNumericalFunctionExpression>(storm::ir::expressions::BaseExpression::double_, qi::_val, qi::_1, storm::ir::expressions::BinaryNumericalFunctionExpression::MINUS))]];
constantDoublePlusExpression.name("constant double expression");
constantDoubleExpression %= constantDoublePlusExpression;
constantDoubleExpression.name("constant double expression");
@ -186,9 +186,9 @@ struct PrismParser::PrismGrammar : qi::grammar<Iterator, storm::ir::Program(), q
commandDefinition.name("command");
// This block defines all entities that are neede for parsing variable definitions.
booleanVariableDefinition = (freeIdentifierName >> qi::lit(":") >> qi::lit("bool") > -(qi::lit("init") > constantBooleanExpression[qi::_b = phoenix::construct<std::shared_ptr<storm::ir::expressions::BaseExpression>>(qi::_1)]) > qi::lit(";"))[phoenix::push_back(qi::_r1, phoenix::construct<storm::ir::BooleanVariable>(phoenix::val(nextBooleanVariableIndex), phoenix::val(qi::_1), qi::_b)), phoenix::insert(qi::_r2, phoenix::construct<std::pair<std::string, uint_fast64_t>>(qi::_1, phoenix::val(nextBooleanVariableIndex))), qi::_a = phoenix::construct<std::shared_ptr<storm::ir::expressions::VariableExpression>>(phoenix::new_<storm::ir::expressions::VariableExpression>(storm::ir::expressions::BaseExpression::bool_, phoenix::val(nextBooleanVariableIndex), qi::_1)), phoenix::bind(booleanVariables_.add, qi::_1, qi::_a), phoenix::bind(booleanVariableNames_.add, qi::_1, qi::_1), phoenix::bind(localBooleanVariables_.add, qi::_1, qi::_1), phoenix::ref(nextBooleanVariableIndex)++];
booleanVariableDefinition = (freeIdentifierName >> qi::lit(":") >> qi::lit("bool") > -(qi::lit("init") > constantBooleanExpression[qi::_b = phoenix::construct<std::shared_ptr<storm::ir::expressions::BaseExpression>>(qi::_1)]) > qi::lit(";"))[phoenix::push_back(qi::_r1, phoenix::construct<storm::ir::BooleanVariable>(phoenix::ref(nextBooleanVariableIndex), phoenix::val(qi::_1), qi::_b)), phoenix::insert(qi::_r2, phoenix::construct<std::pair<std::string, uint_fast64_t>>(qi::_1, phoenix::ref(nextBooleanVariableIndex))), qi::_a = phoenix::construct<std::shared_ptr<storm::ir::expressions::VariableExpression>>(phoenix::new_<storm::ir::expressions::VariableExpression>(storm::ir::expressions::BaseExpression::bool_, phoenix::ref(nextBooleanVariableIndex), qi::_1)), phoenix::bind(booleanVariables_.add, qi::_1, qi::_a), phoenix::bind(booleanVariableNames_.add, qi::_1, qi::_1), phoenix::bind(localBooleanVariables_.add, qi::_1, qi::_1), phoenix::ref(nextBooleanVariableIndex) = phoenix::ref(nextBooleanVariableIndex) + 1];
booleanVariableDefinition.name("boolean variable declaration");
integerVariableDefinition = (freeIdentifierName > qi::lit(":") > qi::lit("[") > constantIntegerExpression > qi::lit("..") > constantIntegerExpression > qi::lit("]") > -(qi::lit("init") > constantIntegerExpression[qi::_b = phoenix::construct<std::shared_ptr<storm::ir::expressions::BaseExpression>>(qi::_1)]) > qi::lit(";"))[phoenix::push_back(qi::_r1, phoenix::construct<storm::ir::IntegerVariable>(phoenix::val(nextIntegerVariableIndex), qi::_1, qi::_2, qi::_3, qi::_b)), phoenix::insert(qi::_r2, phoenix::construct<std::pair<std::string, uint_fast64_t>>(qi::_1, phoenix::val(nextIntegerVariableIndex))), qi::_a = phoenix::construct<std::shared_ptr<storm::ir::expressions::VariableExpression>>(phoenix::new_<storm::ir::expressions::VariableExpression>(storm::ir::expressions::BaseExpression::int_, phoenix::val(nextIntegerVariableIndex), qi::_1)), phoenix::bind(integerVariables_.add, qi::_1, qi::_a), phoenix::bind(integerVariableNames_.add, qi::_1, qi::_1), phoenix::bind(localIntegerVariables_.add, qi::_1, qi::_1), phoenix::ref(nextIntegerVariableIndex)++];
integerVariableDefinition = (freeIdentifierName > qi::lit(":") > qi::lit("[") > constantIntegerExpression > qi::lit("..") > constantIntegerExpression > qi::lit("]") > -(qi::lit("init") > constantIntegerExpression[qi::_b = phoenix::construct<std::shared_ptr<storm::ir::expressions::BaseExpression>>(qi::_1)]) > qi::lit(";"))[phoenix::push_back(qi::_r1, phoenix::construct<storm::ir::IntegerVariable>(phoenix::ref(nextIntegerVariableIndex), qi::_1, qi::_2, qi::_3, qi::_b)), phoenix::insert(qi::_r2, phoenix::construct<std::pair<std::string, uint_fast64_t>>(qi::_1, phoenix::ref(nextIntegerVariableIndex))), qi::_a = phoenix::construct<std::shared_ptr<storm::ir::expressions::VariableExpression>>(phoenix::new_<storm::ir::expressions::VariableExpression>(storm::ir::expressions::BaseExpression::int_, phoenix::ref(nextIntegerVariableIndex), qi::_1)), phoenix::bind(integerVariables_.add, qi::_1, qi::_a), phoenix::bind(integerVariableNames_.add, qi::_1, qi::_1), phoenix::bind(localIntegerVariables_.add, qi::_1, qi::_1), phoenix::ref(nextIntegerVariableIndex) = phoenix::ref(nextIntegerVariableIndex) + 1];
integerVariableDefinition.name("integer variable declaration");
variableDefinition = (booleanVariableDefinition(qi::_r1, qi::_r3) | integerVariableDefinition(qi::_r2, qi::_r4));
variableDefinition.name("variable declaration");
@ -297,18 +297,18 @@ struct PrismParser::PrismGrammar : qi::grammar<Iterator, storm::ir::Program(), q
// Rules with integer result type.
qi::rule<Iterator, std::shared_ptr<storm::ir::expressions::BaseExpression>(), Skipper> integerExpression;
qi::rule<Iterator, std::shared_ptr<storm::ir::expressions::BaseExpression>(), Skipper> integerPlusExpression;
qi::rule<Iterator, std::shared_ptr<storm::ir::expressions::BaseExpression>(), qi::locals<bool>, Skipper> integerPlusExpression;
qi::rule<Iterator, std::shared_ptr<storm::ir::expressions::BaseExpression>(), Skipper> integerMultExpression;
qi::rule<Iterator, std::shared_ptr<storm::ir::expressions::BaseExpression>(), Skipper> atomicIntegerExpression;
qi::rule<Iterator, std::shared_ptr<storm::ir::expressions::BaseExpression>(), Skipper> constantIntegerExpression;
qi::rule<Iterator, std::shared_ptr<storm::ir::expressions::BaseExpression>(), Skipper> constantIntegerPlusExpression;
qi::rule<Iterator, std::shared_ptr<storm::ir::expressions::BaseExpression>(), qi::locals<bool>, Skipper> constantIntegerPlusExpression;
qi::rule<Iterator, std::shared_ptr<storm::ir::expressions::BaseExpression>(), Skipper> constantIntegerMultExpression;
qi::rule<Iterator, std::shared_ptr<storm::ir::expressions::BaseExpression>(), Skipper> constantAtomicIntegerExpression;
// Rules with double result type.
qi::rule<Iterator, std::shared_ptr<storm::ir::expressions::BaseExpression>(), Skipper> constantDoubleExpression;
qi::rule<Iterator, std::shared_ptr<storm::ir::expressions::BaseExpression>(), Skipper> constantDoublePlusExpression;
qi::rule<Iterator, std::shared_ptr<storm::ir::expressions::BaseExpression>(), Skipper> constantDoubleMultExpression;
qi::rule<Iterator, std::shared_ptr<storm::ir::expressions::BaseExpression>(), qi::locals<bool>, Skipper> constantDoublePlusExpression;
qi::rule<Iterator, std::shared_ptr<storm::ir::expressions::BaseExpression>(), qi::locals<bool>, Skipper> constantDoubleMultExpression;
qi::rule<Iterator, std::shared_ptr<storm::ir::expressions::BaseExpression>(), Skipper> constantAtomicDoubleExpression;
// Rules for variable recognition.

17
src/storm.cpp

@ -198,7 +198,7 @@ void testCheckingSynchronousLeader(storm::models::Dtmc<double>& dtmc, uint_fast6
storm::formula::BoundedUntil<double>* boundedUntilFormula = new storm::formula::BoundedUntil<double>(new storm::formula::Ap<double>("true"), electedFormula, 1);
probFormula = new storm::formula::ProbabilisticNoBoundOperator<double>(boundedUntilFormula);
for (uint_fast64_t L = 1; L < 5; ++L) {
for (uint_fast64_t L = 9; L < 10; ++L) {
boundedUntilFormula->setBound(L*(n + 1));
mc->check(*probFormula);
}
@ -224,12 +224,11 @@ void testChecking() {
if (parser.getType() == storm::models::DTMC) {
std::shared_ptr<storm::models::Dtmc<double>> dtmc = parser.getModel<storm::models::Dtmc<double>>();
dtmc->printModelInformationToStream(std::cout);
// testCheckingDie(*dtmc);
// testCheckingCrowds(*dtmc);
testCheckingSynchronousLeader(*dtmc, 5);
}
else std::cout << "Input is not DTMC" << std::endl;
// testCheckingDie(*dtmc);
// testCheckingCrowds(*dtmc);
// testCheckingSynchronousLeader(*dtmc, 4);
}
/*!
@ -237,9 +236,9 @@ void testChecking() {
*/
int main(const int argc, const char* argv[]) {
initializeLogger();
// if (!parseOptions(argc, argv)) {
// return 0;
//}
if (!parseOptions(argc, argv)) {
return 0;
}
// printHeader(argc, argv);
// testChecking();
@ -247,7 +246,7 @@ int main(const int argc, const char* argv[]) {
storm::parser::PrismParser parser;
std::shared_ptr<storm::ir::Program> program = parser.parseFile("test.input");
storm::storage::SparseMatrix<double>* result = storm::adapters::IntermediateRepresentationAdapter::toSparseMatrix<double>(*program);
result->print();
// result->print();
delete result;
cleanUp();

105
test.input

@ -1,24 +1,89 @@
// Knuth's model of a fair die using only fair coins
dtmc
module die
// local state
s : [0..7] init 0;
// value of the dice
d : [0..6] init 0;
[] s=0 -> 0.5 : (s'=1) + 0.5 : (s'=2);
[] s=1 -> 0.5 : (s'=3) + 0.5 : (s'=4);
[] s=2 -> 0.5 : (s'=5) + 0.5 : (s'=6);
[] s=3 -> 0.5 : (s'=1) + 0.5 : (s'=7) & (d'=1);
[] s=4 -> 0.5 : (s'=7) & (d'=2) + 0.5 : (s'=7) & (d'=3);
[] s=5 -> 0.5 : (s'=7) & (d'=4) + 0.5 : (s'=7) & (d'=5);
[] s=6 -> 0.5 : (s'=2) + 0.5 : (s'=7) & (d'=6);
[] s=7 -> 1: (s'=7);
// probability of forwarding
const double PF = 0.8;
const double notPF = .2; // must be 1-PF
// probability that a crowd member is bad
const double badC = .167;
// probability that a crowd member is good
const double goodC = 0.833;
// Total number of protocol runs to analyze
const int TotalRuns = 5;
// size of the crowd
const int CrowdSize = 10;
module crowds
// protocol phase
phase: [0..4] init 0;
// crowd member good (or bad)
good: bool init false;
// number of protocol runs
runCount: [0..TotalRuns] init 0;
// observe_i is the number of times the attacker observed crowd member i
observe0: [0..TotalRuns] init 0;
observe1: [0..TotalRuns] init 0;
observe2: [0..TotalRuns] init 0;
observe3: [0..TotalRuns] init 0;
observe4: [0..TotalRuns] init 0;
observe5: [0..TotalRuns] init 0;
observe6: [0..TotalRuns] init 0;
observe7: [0..TotalRuns] init 0;
observe8: [0..TotalRuns] init 0;
observe9: [0..TotalRuns] init 0;
// the last seen crowd member
lastSeen: [0..CrowdSize - 1] init 0;
// get the protocol started
[] phase=0 & runCount<TotalRuns -> 1: (phase'=1) & (runCount'=runCount+1) & (lastSeen'=0);
// decide whether crowd member is good or bad according to given probabilities
[] phase=1 -> goodC : (phase'=2) & (good'=true) + badC : (phase'=2) & (good'=false);
// if the current member is a good member, update the last seen index (chosen uniformly)
[] phase=2 & good -> 1/10 : (lastSeen'=0) & (phase'=3) + 1/10 : (lastSeen'=1) & (phase'=3) + 1/10 : (lastSeen'=2) & (phase'=3) + 1/10 : (lastSeen'=3) & (phase'=3) + 1/10 : (lastSeen'=4) & (phase'=3) + 1/10 : (lastSeen'=5) & (phase'=3) + 1/10 : (lastSeen'=6) & (phase'=3) + 1/10 : (lastSeen'=7) & (phase'=3) + 1/10 : (lastSeen'=8) & (phase'=3) + 1/10 : (lastSeen'=9) & (phase'=3);
// if the current member is a bad member, record the most recently seen index
[] phase=2 & !good & lastSeen=0 & observe0 < TotalRuns -> 1: (observe0'=observe0+1) & (phase'=4);
[] phase=2 & !good & lastSeen=1 & observe1 < TotalRuns -> 1: (observe1'=observe1+1) & (phase'=4);
[] phase=2 & !good & lastSeen=2 & observe2 < TotalRuns -> 1: (observe2'=observe2+1) & (phase'=4);
[] phase=2 & !good & lastSeen=3 & observe3 < TotalRuns -> 1: (observe3'=observe3+1) & (phase'=4);
[] phase=2 & !good & lastSeen=4 & observe4 < TotalRuns -> 1: (observe4'=observe4+1) & (phase'=4);
[] phase=2 & !good & lastSeen=5 & observe5 < TotalRuns -> 1: (observe5'=observe5+1) & (phase'=4);
[] phase=2 & !good & lastSeen=6 & observe6 < TotalRuns -> 1: (observe6'=observe6+1) & (phase'=4);
[] phase=2 & !good & lastSeen=7 & observe7 < TotalRuns -> 1: (observe7'=observe7+1) & (phase'=4);
[] phase=2 & !good & lastSeen=8 & observe8 < TotalRuns -> 1: (observe8'=observe8+1) & (phase'=4);
[] phase=2 & !good & lastSeen=9 & observe9 < TotalRuns -> 1: (observe9'=observe9+1) & (phase'=4);
// good crowd members forward with probability PF and deliver otherwise
[] phase=3 -> PF : (phase'=1) + notPF : (phase'=4);
// deliver the message and start over
[] phase=4 -> 1: (phase'=0);
endmodule
rewards "coin_flips"
[] s<7 : 1;
endrewards
label "observe0Greater1" = observe0>1;
label "observe1Greater1" = observe1>1;
label "observe2Greater1" = observe2>1;
label "observe3Greater1" = observe3>1;
label "observe4Greater1" = observe4>1;
label "observe5Greater1" = observe5>1;
label "observe6Greater1" = observe6>1;
label "observe7Greater1" = observe7>1;
label "observe8Greater1" = observe8>1;
label "observe9Greater1" = observe9>1;
label "observeIGreater1" = observe1>1|observe2>1|observe3>1|observe4>1|observe5>1|observe6>1|observe7>1|observe8>1|observe9>1;
label "observeOnlyTrueSender" = observe0>1&observe1<=1&observe2<=1&observe3<=1&observe4<=1&observe5<=1&observe6<=1&observe7<=1&observe8<=1&observe9<=1;
Loading…
Cancel
Save