Browse Source

Added reward parsing:

* Transition-based rewards are parsed using the existing (Deterministic)SparseTransitionsParser.
* State-based rewards are parsed using a new SparseStateRewardParser that parses lines consisting of a state and an associated reward.
* The Dtmc class now stores the two reward models.
* The DtmcParser class now parses up to one transition-based and one state-based reward file. They may, however, be omitted in which case the respective reward model is set to null.
tempestpy_adaptions
dehnert 12 years ago
parent
commit
b26a731383
  1. 10
      src/modelChecker/GmmxxDtmcPrctlModelChecker.h
  2. 38
      src/models/Dtmc.h
  3. 24
      src/parser/DtmcParser.cpp
  4. 3
      src/parser/DtmcParser.h
  5. 2
      src/utility/Settings.cpp

10
src/modelChecker/GmmxxDtmcPrctlModelChecker.h

@ -14,7 +14,7 @@
#include "src/modelChecker/DtmcPrctlModelChecker.h"
#include "src/solver/GraphAnalyzer.h"
#include "src/utility/Vector.h"
#include "src/utility/ConstTemplates.h"
#include "src/utility/Settings.h"
#include "gmm/gmm_matrix.h"
@ -56,7 +56,7 @@ public:
// Create the vector with which to multiply.
std::vector<Type>* result = new std::vector<Type>(this->getModel().getNumberOfStates());
mrmc::utility::setVectorValues(result, *rightStates, static_cast<Type>(1.0));
mrmc::utility::setVectorValues(result, *rightStates, mrmc::utility::constGetOne<Type>());
// Now perform matrix-vector multiplication as long as we meet the bound of the formula.
std::vector<Type>* swap = nullptr;
@ -84,7 +84,7 @@ public:
// Create the vector with which to multiply and initialize it correctly.
std::vector<Type> x(this->getModel().getNumberOfStates());
mrmc::utility::setVectorValues(&x, *nextStates, static_cast<Type>(1.0));
mrmc::utility::setVectorValues(&x, *nextStates, mrmc::utility::constGetOne<Type>());
// Delete obsolete sub-result.
delete nextStates;
@ -221,8 +221,8 @@ public:
}
// Set values of resulting vector that are known exactly.
mrmc::utility::setVectorValues<Type>(result, notExistsPhiUntilPsiStates, static_cast<Type>(0));
mrmc::utility::setVectorValues<Type>(result, alwaysPhiUntilPsiStates, static_cast<Type>(1.0));
mrmc::utility::setVectorValues<Type>(result, notExistsPhiUntilPsiStates, mrmc::utility::constGetZero<Type>());
mrmc::utility::setVectorValues<Type>(result, alwaysPhiUntilPsiStates, mrmc::utility::constGetOne<Type>());
return result;
}

38
src/models/Dtmc.h

@ -39,9 +39,14 @@ public:
* @param stateLabeling The labeling that assigns a set of atomic
* propositions to each state.
*/
Dtmc(std::shared_ptr<mrmc::storage::SquareSparseMatrix<T>> probabilityMatrix, std::shared_ptr<mrmc::models::AtomicPropositionsLabeling> stateLabeling)
: probabilityMatrix(probabilityMatrix), stateLabeling(stateLabeling), backwardTransitions(nullptr) {
if (! this->checkValidityProbabilityMatrix()) {
Dtmc(std::shared_ptr<mrmc::storage::SquareSparseMatrix<T>> probabilityMatrix,
std::shared_ptr<mrmc::models::AtomicPropositionsLabeling> stateLabeling,
std::shared_ptr<std::vector<T>> stateRewards = nullptr,
std::shared_ptr<mrmc::storage::SquareSparseMatrix<T>> transitionRewardMatrix = nullptr)
: probabilityMatrix(probabilityMatrix), stateLabeling(stateLabeling),
stateRewards(stateRewards), transitionRewardMatrix(transitionRewardMatrix),
backwardTransitions(nullptr) {
if (!this->checkValidityProbabilityMatrix()) {
std::cerr << "Probability matrix is invalid" << std::endl;
}
}
@ -52,11 +57,12 @@ public:
* @param dtmc A reference to the DTMC that is to be copied.
*/
Dtmc(const Dtmc<T> &dtmc) : probabilityMatrix(dtmc.probabilityMatrix),
stateLabeling(dtmc.stateLabeling) {
stateLabeling(dtmc.stateLabeling), stateRewards(dtmc.stateRewards),
transitionRewardMatrix(dtmc.transitionRewardMatrix) {
if (dtmc.backardTransitions != nullptr) {
this->backwardTransitions = new mrmc::models::GraphTransitions<T>(*dtmc.backwardTransitions);
}
if (! this->checkValidityProbabilityMatrix()) {
if (!this->checkValidityProbabilityMatrix()) {
std::cerr << "Probability matrix is invalid" << std::endl;
}
}
@ -108,6 +114,22 @@ public:
return this->probabilityMatrix;
}
/*!
* Returns a pointer to the matrix representing the transition rewards.
* @return A pointer to the matrix representing the transition rewards.
*/
std::shared_ptr<mrmc::storage::SquareSparseMatrix<T>> getTransitionRewardMatrix() const {
return this->transitionRewardMatrix;
}
/*!
* Returns a pointer to the vector representing the state rewards.
* @return A pointer to the vector representing the state rewards.
*/
std::shared_ptr<std::vector<T>> getStateRewards() const {
return this->stateRewards;
}
/*!
*
*/
@ -170,6 +192,12 @@ private:
/*! The labeling of the states of the DTMC. */
std::shared_ptr<mrmc::models::AtomicPropositionsLabeling> stateLabeling;
/*! The state-based rewards of the DTMC. */
std::shared_ptr<std::vector<T>> stateRewards;
/*! The transition-based rewards of the DTMC. */
std::shared_ptr<mrmc::storage::SquareSparseMatrix<T>> transitionRewardMatrix;
/*!
* A data structure that stores the predecessors for all states. This is
* needed for backwards directed searches.

24
src/parser/DtmcParser.cpp

@ -8,6 +8,7 @@
#include "DtmcParser.h"
#include "DeterministicSparseTransitionParser.h"
#include "AtomicPropositionLabelingParser.h"
#include "SparseStateRewardParser.h"
namespace mrmc {
namespace parser {
@ -19,15 +20,30 @@ namespace parser {
*
* @param transitionSystemFile String containing the location of the transition file (....tra)
* @param labelingFile String containing the location of the labeling file (....lab)
* @param stateRewardFile String containing the location of the state reward file (...srew)
* @param transitionRewardFile String containing the location of the transition reward file (...trew)
*/
DtmcParser::DtmcParser(std::string const & transitionSystemFile, std::string const & labelingFile) {
DtmcParser::DtmcParser(std::string const & transitionSystemFile, std::string const & labelingFile,
std::string const & stateRewardFile, std::string const & transitionRewardFile) {
mrmc::parser::DeterministicSparseTransitionParser tp(transitionSystemFile);
uint_fast64_t nodeCount = tp.getMatrix()->getRowCount();
uint_fast64_t stateCount = tp.getMatrix()->getRowCount();
mrmc::parser::AtomicPropositionLabelingParser lp(nodeCount, labelingFile);
std::shared_ptr<std::vector<double>> stateRewards = nullptr;
std::shared_ptr<mrmc::storage::SquareSparseMatrix<double>> transitionRewards = nullptr;
dtmc = std::shared_ptr<mrmc::models::Dtmc<double> >(new mrmc::models::Dtmc<double>(tp.getMatrix(), lp.getLabeling()));
mrmc::parser::AtomicPropositionLabelingParser lp(stateCount, labelingFile);
if (stateRewardFile != "") {
mrmc::parser::SparseStateRewardParser srp(stateCount, stateRewardFile);
stateRewards = srp.getStateRewards();
}
if (transitionRewardFile != "") {
mrmc::parser::DeterministicSparseTransitionParser trp(transitionRewardFile);
transitionRewards = trp.getMatrix();
}
dtmc = std::shared_ptr<mrmc::models::Dtmc<double>>(new mrmc::models::Dtmc<double>(tp.getMatrix(), lp.getLabeling(), stateRewards, transitionRewards));
}
} /* namespace parser */
} /* namespace mrmc */

3
src/parser/DtmcParser.h

@ -24,7 +24,8 @@ namespace parser {
*/
class DtmcParser: public mrmc::parser::Parser {
public:
DtmcParser(std::string const & transitionSystemFile, std::string const & labelingFile);
DtmcParser(std::string const & transitionSystemFile, std::string const & labelingFile,
std::string const & stateRewardFile = "", std::string const & transitionRewardFile = "");
std::shared_ptr<mrmc::models::Dtmc<double>> getDtmc() {
return this->dtmc;

2
src/utility/Settings.cpp

@ -130,6 +130,8 @@ void Settings::initDescriptions() {
("test-prctl", bpo::value<std::string>(), "name of prctl file")
("trafile", bpo::value<std::string>()->required(), "name of the .tra file")
("labfile", bpo::value<std::string>()->required(), "name of the .lab file")
("transrew", bpo::value<std::string>()->default_value(""), "name of transition reward file")
("staterew", bpo::value<std::string>()->default_value(""), "name of state reward file")
;
}

Loading…
Cancel
Save