Browse Source

statesOfCoalition are now part of LRAViHelper

tempestpy_adaptions
Stefan Pranger 4 years ago
parent
commit
7eb5cd98a5
  1. 2
      src/storm/modelchecker/helper/infinitehorizon/SparseNondeterministicGameInfiniteHorizonHelper.cpp
  2. 21
      src/storm/modelchecker/helper/infinitehorizon/internal/LraViHelper.cpp
  3. 41
      src/storm/modelchecker/helper/infinitehorizon/internal/LraViHelper.h

2
src/storm/modelchecker/helper/infinitehorizon/SparseNondeterministicGameInfiniteHorizonHelper.cpp

@ -116,7 +116,7 @@ namespace storm {
if (this->isContinuousTime()) { if (this->isContinuousTime()) {
STORM_LOG_THROW(false, storm::exceptions::InternalException, "We cannot handle continuous time games."); STORM_LOG_THROW(false, storm::exceptions::InternalException, "We cannot handle continuous time games.");
} else { } else {
storm::modelchecker::helper::internal::LraViHelper<ValueType, storm::storage::MaximalEndComponent, storm::modelchecker::helper::internal::LraViTransitionsType::GameNondetTsNoIs> viHelper(mec, this->_transitionMatrix, aperiodicFactor);
storm::modelchecker::helper::internal::LraViHelper<ValueType, storm::storage::MaximalEndComponent, storm::modelchecker::helper::internal::LraViTransitionsType::GameNondetTsNoIs> viHelper(mec, this->_transitionMatrix, aperiodicFactor, nullptr, nullptr, &statesOfCoalition);
return viHelper.performValueIteration(env, stateRewardsGetter, actionRewardsGetter, nullptr, &this->getOptimizationDirection(), optimalChoices); return viHelper.performValueIteration(env, stateRewardsGetter, actionRewardsGetter, nullptr, &this->getOptimizationDirection(), optimalChoices);
} }
} }

21
src/storm/modelchecker/helper/infinitehorizon/internal/LraViHelper.cpp

@ -26,7 +26,7 @@ namespace storm {
namespace internal { namespace internal {
template <typename ValueType, typename ComponentType, LraViTransitionsType TransitionsType> template <typename ValueType, typename ComponentType, LraViTransitionsType TransitionsType>
LraViHelper<ValueType, ComponentType, TransitionsType>::LraViHelper(ComponentType const& component, storm::storage::SparseMatrix<ValueType> const& transitionMatrix, ValueType const& aperiodicFactor, storm::storage::BitVector const* timedStates, std::vector<ValueType> const* exitRates) : _transitionMatrix(transitionMatrix), _timedStates(timedStates), _hasInstantStates(TransitionsType == LraViTransitionsType::DetTsNondetIs || TransitionsType == LraViTransitionsType::DetTsDetIs), _Tsx1IsCurrent(false) {
LraViHelper<ValueType, ComponentType, TransitionsType>::LraViHelper(ComponentType const& component, storm::storage::SparseMatrix<ValueType> const& transitionMatrix, ValueType const& aperiodicFactor, storm::storage::BitVector const* timedStates, std::vector<ValueType> const* exitRates, storm::storage::BitVector const* statesOfCoalition) : _transitionMatrix(transitionMatrix), _timedStates(timedStates), _hasInstantStates(TransitionsType == LraViTransitionsType::DetTsNondetIs || TransitionsType == LraViTransitionsType::DetTsDetIs), _Tsx1IsCurrent(false), _statesOfCoalition(statesOfCoalition) {
setComponent(component); setComponent(component);
// Run through the component and collect some data: // Run through the component and collect some data:
@ -167,6 +167,9 @@ namespace storm {
if (env.solver().lra().isMaximalIterationCountSet()) { if (env.solver().lra().isMaximalIterationCountSet()) {
maxIter = env.solver().lra().getMaximalIterationCount(); maxIter = env.solver().lra().getMaximalIterationCount();
} }
if(gameNondetTs()) {
STORM_LOG_ASSERT(_statesOfCoalition != nullptr, "Tried to solve LRA problem for a game, but coalition states have not been set.");
}
// start the iterations // start the iterations
ValueType result = storm::utility::zero<ValueType>(); ValueType result = storm::utility::zero<ValueType>();
@ -218,13 +221,13 @@ namespace storm {
} }
performIterationStep(env, dir, choices); performIterationStep(env, dir, choices);
} }
std::cout << "result (" << iter << " steps):" << std::endl;
//std::cout << "result (" << iter << " steps):" << std::endl;
storm::utility::vector::applyPointwise<ValueType, ValueType>(xNew(), xNew(), [&iter] (ValueType const& x_i) -> ValueType { return x_i / (double)iter; }); storm::utility::vector::applyPointwise<ValueType, ValueType>(xNew(), xNew(), [&iter] (ValueType const& x_i) -> ValueType { return x_i / (double)iter; });
for(int i = 0; i < xNew().size() ; i++ ) {
std::cout << std::setprecision(4) << i << "\t: " << xNew().at(i) << "\t" << xNew().at(i) * _uniformizationRate << "\t" << std::setprecision(16) << xOld().at(i) *_uniformizationRate << std::endl;
//if(i == 50) {std::cout << "only showing top 50 lines"; break; }
}
if(gameNondetTs()) result = xOld().at(0) * _uniformizationRate; // TODO is "init" always going to be .at(0) ?
//for(int i = 0; i < xNew().size() ; i++ ) {
// std::cout << std::setprecision(4) << i << "\t: " << xNew().at(i) << "\t" << xNew().at(i) * _uniformizationRate << "\t" << std::setprecision(16) << xOld().at(i) *_uniformizationRate << std::endl;
// //if(i == 50) {std::cout << "only showing top 50 lines"; break; }
//}
if(gameNondetTs()) result = (xOld().at(0) * _uniformizationRate)/(double)iter; // TODO is "init" always going to be .at(0) ?
return result; return result;
} }
@ -386,10 +389,10 @@ namespace storm {
} }
} else if(gameNondetTs()) { // TODO DRYness? exact same behaviour as case above? } else if(gameNondetTs()) { // TODO DRYness? exact same behaviour as case above?
if (choices == nullptr) { if (choices == nullptr) {
_TsMultiplier->multiplyAndReduce(env, *dir, xOld(), &_TsChoiceValues, xNew());
_TsMultiplier->multiplyAndReduce(env, *dir, xOld(), &_TsChoiceValues, xNew(), nullptr, _statesOfCoalition);
} else { } else {
std::vector<uint64_t> tsChoices(_TsTransitions.getRowGroupCount()); std::vector<uint64_t> tsChoices(_TsTransitions.getRowGroupCount());
_TsMultiplier->multiplyAndReduce(env, *dir, xOld(), &_TsChoiceValues, xNew(), &tsChoices);
_TsMultiplier->multiplyAndReduce(env, *dir, xOld(), &_TsChoiceValues, xNew(), &tsChoices, _statesOfCoalition);
setInputModelChoices(*choices, tsChoices); // no components -> no need for that call? setInputModelChoices(*choices, tsChoices); // no components -> no need for that call?
} }
} else { } else {

41
src/storm/modelchecker/helper/infinitehorizon/internal/LraViHelper.h

@ -8,12 +8,12 @@
namespace storm { namespace storm {
class Environment; class Environment;
namespace modelchecker { namespace modelchecker {
namespace helper { namespace helper {
namespace internal { namespace internal {
/*! /*!
* Specifies differnt kinds of transition types with which this helper can be used * Specifies differnt kinds of transition types with which this helper can be used
* Ts means timed states (cf. Markovian states in a Markov Automaton) and Is means instant states (cf. probabilistic states in a Markov automaton). * Ts means timed states (cf. Markovian states in a Markov Automaton) and Is means instant states (cf. probabilistic states in a Markov automaton).
@ -45,7 +45,7 @@ namespace storm {
public: public:
/// Function mapping from indices to values /// Function mapping from indices to values
typedef std::function<ValueType(uint64_t)> ValueGetter; typedef std::function<ValueType(uint64_t)> ValueGetter;
/*! /*!
* Initializes a new VI helper for the provided MEC or BSCC * Initializes a new VI helper for the provided MEC or BSCC
* @param component the MEC or BSCC * @param component the MEC or BSCC
@ -55,7 +55,7 @@ namespace storm {
* @param exitRates The exit rates of the timed states (relevant for continuous time models). If nullptr, all rates are assumed to be 1 (which corresponds to a discrete time model) * @param exitRates The exit rates of the timed states (relevant for continuous time models). If nullptr, all rates are assumed to be 1 (which corresponds to a discrete time model)
* @note All indices and vectors must be w.r.t. the states as described by the provided transition matrix * @note All indices and vectors must be w.r.t. the states as described by the provided transition matrix
*/ */
LraViHelper(ComponentType const& component, storm::storage::SparseMatrix<ValueType> const& transitionMatrix, ValueType const& aperiodicFactor, storm::storage::BitVector const* timedStates = nullptr, std::vector<ValueType> const* exitRates = nullptr);
LraViHelper(ComponentType const& component, storm::storage::SparseMatrix<ValueType> const& transitionMatrix, ValueType const& aperiodicFactor, storm::storage::BitVector const* timedStates = nullptr, std::vector<ValueType> const* exitRates = nullptr, storm::storage::BitVector const* statesOfCoalition = nullptr);
/*! /*!
* Performs value iteration with the given state- and action values. * Performs value iteration with the given state- and action values.
@ -69,9 +69,9 @@ namespace storm {
* @note it is possible to call this method multiple times with different values. However, other changes to the environment or the optimization direction might not have the expected effect due to caching. * @note it is possible to call this method multiple times with different values. However, other changes to the environment or the optimization direction might not have the expected effect due to caching.
*/ */
ValueType performValueIteration(Environment const& env, ValueGetter const& stateValueGetter, ValueGetter const& actionValueGetter, std::vector<ValueType> const* exitRates = nullptr, storm::solver::OptimizationDirection const* dir = nullptr, std::vector<uint64_t>* choices = nullptr); ValueType performValueIteration(Environment const& env, ValueGetter const& stateValueGetter, ValueGetter const& actionValueGetter, std::vector<ValueType> const* exitRates = nullptr, storm::solver::OptimizationDirection const* dir = nullptr, std::vector<uint64_t>* choices = nullptr);
private: private:
/*! /*!
* Initializes the value iterations with the provided values. * Initializes the value iterations with the provided values.
* Resets all information from potential previous calls. * Resets all information from potential previous calls.
@ -80,7 +80,7 @@ namespace storm {
* @param stateValueGetter Function that returns for each global choice index (w.r.t. the input transitions) the value (e.g. reward) for that choice * @param stateValueGetter Function that returns for each global choice index (w.r.t. the input transitions) the value (e.g. reward) for that choice
*/ */
void initializeNewValues(ValueGetter const& stateValueGetter, ValueGetter const& actionValueGetter, std::vector<ValueType> const* exitRates = nullptr); void initializeNewValues(ValueGetter const& stateValueGetter, ValueGetter const& actionValueGetter, std::vector<ValueType> const* exitRates = nullptr);
/*! /*!
* Performs a single iteration step. * Performs a single iteration step.
* @param env The environment. * @param env The environment.
@ -90,41 +90,41 @@ namespace storm {
* @pre when calling this the first time, initializeNewValues must have been called before. Moreover, prepareNextIteration must be called between two calls of this. * @pre when calling this the first time, initializeNewValues must have been called before. Moreover, prepareNextIteration must be called between two calls of this.
*/ */
void performIterationStep(Environment const& env, storm::solver::OptimizationDirection const* dir = nullptr, std::vector<uint64_t>* choices = nullptr); void performIterationStep(Environment const& env, storm::solver::OptimizationDirection const* dir = nullptr, std::vector<uint64_t>* choices = nullptr);
struct ConvergenceCheckResult { struct ConvergenceCheckResult {
bool isPrecisionAchieved; bool isPrecisionAchieved;
ValueType currentValue; ValueType currentValue;
}; };
/*! /*!
* Checks whether the curently computed value achieves the desired precision * Checks whether the curently computed value achieves the desired precision
*/ */
ConvergenceCheckResult checkConvergence(bool relative, ValueType precision) const; ConvergenceCheckResult checkConvergence(bool relative, ValueType precision) const;
/*! /*!
* Must be called between two calls of performIterationStep. * Must be called between two calls of performIterationStep.
*/ */
void prepareNextIteration(Environment const& env); void prepareNextIteration(Environment const& env);
/// Prepares the necessary solvers and multipliers for doing the iterations. /// Prepares the necessary solvers and multipliers for doing the iterations.
void prepareSolversAndMultipliers(Environment const& env, storm::solver::OptimizationDirection const* dir = nullptr); void prepareSolversAndMultipliers(Environment const& env, storm::solver::OptimizationDirection const* dir = nullptr);
void setInputModelChoices(std::vector<uint64_t>& choices, std::vector<uint64_t> const& localMecChoices, bool setChoiceZeroToMarkovianStates = false, bool setChoiceZeroToProbabilisticStates = false) const; void setInputModelChoices(std::vector<uint64_t>& choices, std::vector<uint64_t> const& localMecChoices, bool setChoiceZeroToMarkovianStates = false, bool setChoiceZeroToProbabilisticStates = false) const;
/// Returns true iff the given state is a timed state /// Returns true iff the given state is a timed state
bool isTimedState(uint64_t const& inputModelStateIndex) const; bool isTimedState(uint64_t const& inputModelStateIndex) const;
/// The result for timed states of the most recent iteration /// The result for timed states of the most recent iteration
std::vector<ValueType>& xNew(); std::vector<ValueType>& xNew();
std::vector<ValueType> const& xNew() const; std::vector<ValueType> const& xNew() const;
/// The result for timed states of the previous iteration /// The result for timed states of the previous iteration
std::vector<ValueType>& xOld(); std::vector<ValueType>& xOld();
std::vector<ValueType> const& xOld() const; std::vector<ValueType> const& xOld() const;
/// @return true iff there potentially is a nondeterministic choice at timed states /// @return true iff there potentially is a nondeterministic choice at timed states
bool nondetTs() const; bool nondetTs() const;
/// @return true iff there potentially is a nondeterministic choice at instant states. Returns false if there are no instant states. /// @return true iff there potentially is a nondeterministic choice at instant states. Returns false if there are no instant states.
bool nondetIs() const; bool nondetIs() const;
@ -132,14 +132,15 @@ namespace storm {
bool gameNondetTs() const; bool gameNondetTs() const;
void setComponent(ComponentType component); void setComponent(ComponentType component);
// We need to make sure that states/choices will be processed in ascending order // We need to make sure that states/choices will be processed in ascending order
typedef std::map<uint64_t, std::set<uint64_t>> InternalComponentType; typedef std::map<uint64_t, std::set<uint64_t>> InternalComponentType;
InternalComponentType _component; InternalComponentType _component;
storm::storage::SparseMatrix<ValueType> const& _transitionMatrix; storm::storage::SparseMatrix<ValueType> const& _transitionMatrix;
storm::storage::BitVector const* _timedStates; // e.g. Markovian states of a Markov automaton. storm::storage::BitVector const* _timedStates; // e.g. Markovian states of a Markov automaton.
storm::storage::BitVector const* _statesOfCoalition;
bool _hasInstantStates; bool _hasInstantStates;
ValueType _uniformizationRate; ValueType _uniformizationRate;
storm::storage::SparseMatrix<ValueType> _TsTransitions, _TsToIsTransitions, _IsTransitions, _IsToTsTransitions; storm::storage::SparseMatrix<ValueType> _TsTransitions, _TsToIsTransitions, _IsTransitions, _IsToTsTransitions;
@ -154,4 +155,4 @@ namespace storm {
} }
} }
} }
}
}
Loading…
Cancel
Save