24namespace modelchecker {
27template<
typename ValueType>
34template<
typename ValueType>
41template<
typename ValueType>
43 STORM_LOG_ASSERT(this->isProduceSchedulerSet(),
"Trying to get the produced optimal choices although no scheduler was requested.");
45 "Trying to get the produced optimal choices but none were available. Was there a computation call before?");
46 return this->_producedOptimalChoices.get();
49template<
typename ValueType>
51 STORM_LOG_ASSERT(this->isProduceSchedulerSet(),
"Trying to get the produced optimal choices although no scheduler was requested.");
53 "Trying to get the produced optimal choices but none were available. Was there a computation call before?");
54 return this->_producedOptimalChoices.get();
57template<
typename ValueType>
59 auto const& optimalChoices = getProducedOptimalChoices();
61 for (uint64_t state = 0; state < optimalChoices.size(); ++state) {
62 scheduler.
setChoice(optimalChoices[state], state);
67template<
typename ValueType>
69 if (this->_longRunComponentDecomposition ==
nullptr) {
71 this->createBackwardTransitions();
72 this->_computedLongRunComponentDecomposition =
73 std::make_unique<storm::storage::MaximalEndComponentDecomposition<ValueType>>(this->_transitionMatrix, *this->_backwardTransitions);
74 this->_longRunComponentDecomposition = this->_computedLongRunComponentDecomposition.get();
78template<
typename ValueType>
85 if (this->isProduceSchedulerSet()) {
86 if (!this->_producedOptimalChoices.is_initialized()) {
87 this->_producedOptimalChoices.emplace();
89 this->_producedOptimalChoices->resize(this->_transitionMatrix.getRowGroupCount());
92 auto trivialResult = this->computeLraForTrivialMec(env, stateRewardsGetter, actionRewardsGetter, component);
93 if (trivialResult.first) {
94 return trivialResult.second;
100 method != storm::solver::LraMethod::LinearProgramming) {
102 "Selecting 'LP' as the solution technique for long-run properties to guarantee exact results. If you want to override this, please explicitly "
103 "specify a different LRA method.");
104 method = storm::solver::LraMethod::LinearProgramming;
107 "Selecting 'VI' as the solution technique for long-run properties to guarantee sound results. If you want to override this, please explicitly "
108 "specify a different LRA method.");
109 method = storm::solver::LraMethod::ValueIteration;
111 STORM_LOG_ERROR_COND(!this->isProduceSchedulerSet() || method == storm::solver::LraMethod::ValueIteration,
112 "Scheduler generation not supported for the chosen LRA method. Try value-iteration.");
113 if (method == storm::solver::LraMethod::LinearProgramming) {
114 return computeLraForMecLp(env, stateRewardsGetter, actionRewardsGetter, component);
115 }
else if (method == storm::solver::LraMethod::ValueIteration) {
116 return computeLraForMecVi(env, stateRewardsGetter, actionRewardsGetter, component);
118 STORM_LOG_THROW(
false, storm::exceptions::InvalidSettingsException,
"Unsupported technique.");
122template<
typename ValueType>
127 if (component.
size() == 1) {
128 auto const& element = *component.
begin();
131 if (!this->isContinuousTime()) {
134 ValueType bestValue = actionRewardsGetter(*choiceIt);
135 uint64_t bestChoice = *choiceIt;
137 ValueType currentValue = actionRewardsGetter(*choiceIt);
138 if ((this->minimize() && currentValue < bestValue) || (this->maximize() && currentValue > bestValue)) {
139 bestValue = std::move(currentValue);
140 bestChoice = *choiceIt;
143 if (this->isProduceSchedulerSet()) {
144 this->_producedOptimalChoices.get()[state] = bestChoice - this->_transitionMatrix.getRowGroupIndices()[state];
146 bestValue += stateRewardsGetter(state);
147 return {
true, bestValue};
152 "Nondeterministic continuous time model without Markovian states... Is this a not a Markov Automaton?");
153 STORM_LOG_THROW(this->_markovianStates->get(state), storm::exceptions::InvalidOperationException,
154 "Markov Automaton has Zeno behavior. Computation of Long Run Average values not supported.");
156 if (this->isProduceSchedulerSet()) {
157 this->_producedOptimalChoices.get()[state] = 0;
159 ValueType result = stateRewardsGetter(state) +
160 (this->isContinuousTime() ? (*this->_exitRates)[state] * actionRewardsGetter(*choiceIt) : actionRewardsGetter(*choiceIt));
161 return {
true, result};
164 return {
false, storm::utility::zero<ValueType>()};
167template<
typename ValueType>
173 std::vector<uint64_t>* optimalChoices =
nullptr;
174 if (this->isProduceSchedulerSet()) {
175 optimalChoices = &this->_producedOptimalChoices.get();
179 if (this->isContinuousTime()) {
183 viHelper(mec, this->_transitionMatrix, aperiodicFactor, this->_markovianStates, this->_exitRates);
184 return viHelper.performValueIteration(env, stateRewardsGetter, actionRewardsGetter, this->_exitRates, &this->getOptimizationDirection(),
190 viHelper(mec, this->_transitionMatrix, aperiodicFactor);
191 return viHelper.performValueIteration(env, stateRewardsGetter, actionRewardsGetter,
nullptr, &this->getOptimizationDirection(), optimalChoices);
195template<
typename ValueType>
200 auto solver = storm::utility::solver::getLpSolver<ValueType>(
"LRA for MEC");
204 solver->setOptimizationDirection(invert(this->getOptimizationDirection()));
208 std::map<uint_fast64_t, storm::expressions::Variable> stateToVariableMap;
209 for (
auto const& stateChoicesPair : mec) {
210 std::string variableName =
"x" + std::to_string(stateChoicesPair.first);
211 stateToVariableMap[stateChoicesPair.first] = solver->addUnboundedContinuousVariable(variableName);
217 for (
auto const& stateChoicesPair : mec) {
218 uint_fast64_t state = stateChoicesPair.first;
219 bool stateIsMarkovian = this->_markovianStates && this->_markovianStates->get(state);
223 for (
auto choice : stateChoicesPair.second) {
224 std::vector<storm::expressions::Expression> summands;
225 auto matrixRow = this->_transitionMatrix.getRow(choice);
226 summands.reserve(matrixRow.getNumberOfEntries() + 2);
228 if (stateIsMarkovian) {
230 }
else if (!this->isContinuousTime()) {
231 summands.push_back(-k);
234 for (
auto const& element : matrixRow) {
235 summands.push_back(stateToVariableMap.at(element.getColumn()) * solver->getConstant(element.getValue()));
239 if (stateIsMarkovian) {
241 value = stateRewardsGetter(state) / (*this->_exitRates)[state] + actionRewardsGetter(choice);
242 }
else if (!this->isContinuousTime()) {
244 value = stateRewardsGetter(state) + actionRewardsGetter(choice);
247 value = actionRewardsGetter(choice);
249 summands.push_back(solver->getConstant(value));
251 if (this->minimize()) {
256 solver->addConstraint(
"s" + std::to_string(state) +
"," + std::to_string(choice), constraint);
261 STORM_LOG_THROW(!this->isProduceSchedulerSet(), storm::exceptions::NotImplementedException,
262 "Scheduler extraction is not yet implemented for LP based LRA method.");
263 return solver->getContinuousValue(k);
271template<
typename ValueType>
273 std::vector<uint64_t>
const& inputToSspStateMap, uint64_t
const& numberOfNonComponentStates, uint64_t
const& currentSspChoice,
276 std::map<uint64_t, ValueType> auxiliaryStateToProbabilityMap;
278 for (
auto const& transition : inputTransitionMatrix.
getRow(inputMatrixChoice)) {
280 auto const& sspTransitionTarget = inputToSspStateMap[transition.getColumn()];
283 if (sspTransitionTarget < numberOfNonComponentStates) {
285 sspMatrixBuilder.
addNextValue(currentSspChoice, sspTransitionTarget, transition.getValue());
289 auto insertionRes = auxiliaryStateToProbabilityMap.emplace(sspTransitionTarget, transition.getValue());
290 if (!insertionRes.second) {
293 insertionRes.first->second += transition.getValue();
300 for (
auto const& componentToProbEntry : auxiliaryStateToProbabilityMap) {
301 sspMatrixBuilder.
addNextValue(currentSspChoice, componentToProbEntry.first, componentToProbEntry.second);
305template<
typename ValueType>
307 std::vector<ValueType>
const& mecLraValues, std::vector<uint64_t>
const& inputToSspStateMap,
storm::storage::BitVector const& statesNotInComponent,
308 uint64_t numberOfNonComponentStates, std::vector<std::pair<uint64_t, uint64_t>>* sspComponentExitChoicesToOriginalMap) {
309 auto const& choiceIndices = this->_transitionMatrix.getRowGroupIndices();
311 std::vector<ValueType> rhs;
312 uint64_t numberOfSspStates = numberOfNonComponentStates + this->_longRunComponentDecomposition->size();
315 uint64_t currentSspChoice = 0;
316 for (
auto nonComponentState : statesNotInComponent) {
318 for (uint64_t choice = choiceIndices[nonComponentState]; choice < choiceIndices[nonComponentState + 1]; ++choice, ++currentSspChoice) {
319 rhs.push_back(storm::utility::zero<ValueType>());
320 addSspMatrixChoice(choice, this->_transitionMatrix, inputToSspStateMap, numberOfNonComponentStates, currentSspChoice, sspMatrixBuilder);
324 for (uint64_t componentIndex = 0; componentIndex < this->_longRunComponentDecomposition->size(); ++componentIndex) {
325 auto const& component = (*this->_longRunComponentDecomposition)[componentIndex];
329 for (
auto const& element : component) {
331 for (uint64_t choice = choiceIndices[componentState]; choice < choiceIndices[componentState + 1]; ++choice) {
334 rhs.push_back(storm::utility::zero<ValueType>());
335 addSspMatrixChoice(choice, this->_transitionMatrix, inputToSspStateMap, numberOfNonComponentStates, currentSspChoice, sspMatrixBuilder);
336 if (sspComponentExitChoicesToOriginalMap) {
338 sspComponentExitChoicesToOriginalMap->emplace_back(componentState, choice - choiceIndices[componentState]);
345 rhs.push_back(mecLraValues[componentIndex]);
346 if (sspComponentExitChoicesToOriginalMap) {
348 sspComponentExitChoicesToOriginalMap->emplace_back(std::numeric_limits<uint_fast64_t>::max(), std::numeric_limits<uint_fast64_t>::max());
352 return std::make_pair(sspMatrixBuilder.
build(currentSspChoice, numberOfSspStates, numberOfSspStates), std::move(rhs));
355template<
typename ValueType>
359 std::vector<std::pair<uint64_t, uint64_t>>
const& sspComponentExitChoicesToOriginalMap) {
367 uint64_t exitChoiceOffset = sspMatrix.
getRowGroupIndices()[numberOfNonComponentStates];
368 for (
auto const& mec : *this->_longRunComponentDecomposition) {
370 auto const& sspState = inputToSspStateMap[mec.begin()->first];
371 uint64_t sspChoiceIndex = sspMatrix.
getRowGroupIndices()[sspState] + sspChoices[sspState];
373 auto const& originalStateChoice = sspComponentExitChoicesToOriginalMap[sspChoiceIndex - exitChoiceOffset];
375 if (originalStateChoice.first == std::numeric_limits<uint_fast64_t>::max()) {
382 this->_producedOptimalChoices.get()[originalStateChoice.first] = originalStateChoice.second;
386 for (
auto const& stateActions : mec) {
387 if (stateActions.first != originalStateChoice.first) {
388 this->_producedOptimalChoices.get()[stateActions.first] = std::numeric_limits<uint64_t>::max();
392 this->createBackwardTransitions();
394 std::vector<uint64_t> stack = {originalStateChoice.first};
395 while (!stack.empty()) {
396 uint64_t currentState = stack.back();
398 for (
auto const& backwardsTransition : this->_backwardTransitions->getRowGroup(currentState)) {
399 uint64_t predecessorState = backwardsTransition.getColumn();
400 if (mec.containsState(predecessorState)) {
401 auto& selectedPredChoice = this->_producedOptimalChoices.get()[predecessorState];
402 if (selectedPredChoice == std::numeric_limits<uint64_t>::max()) {
405 for (
auto const& predChoice : mec.getChoicesForState(predecessorState)) {
406 for (
auto const& forwardTransition : this->_transitionMatrix.getRow(predChoice)) {
407 if (forwardTransition.getColumn() == currentState && !
storm::utility::isZero(forwardTransition.getValue())) {
409 selectedPredChoice = predChoice - this->_transitionMatrix.getRowGroupIndices()[predecessorState];
410 stack.push_back(predecessorState);
414 if (selectedPredChoice != std::numeric_limits<uint64_t>::max()) {
426template<
typename ValueType>
428 std::vector<ValueType>
const& componentLraValues) {
429 STORM_LOG_ASSERT(this->_longRunComponentDecomposition !=
nullptr,
"Decomposition not computed, yet.");
438 std::vector<uint64_t> inputToSspStateMap(this->_transitionMatrix.getRowGroupCount(), std::numeric_limits<uint64_t>::max());
439 for (uint64_t currentComponentIndex = 0; currentComponentIndex < this->_longRunComponentDecomposition->size(); ++currentComponentIndex) {
440 for (
auto const& element : (*this->_longRunComponentDecomposition)[currentComponentIndex]) {
442 statesInComponents.
set(state);
443 inputToSspStateMap[state] = currentComponentIndex;
447 uint64_t numberOfNonComponentStates = 0;
449 for (
auto nonComponentState : statesNotInComponent) {
450 inputToSspStateMap[nonComponentState] = numberOfNonComponentStates;
451 ++numberOfNonComponentStates;
456 for (
auto mecState : statesInComponents) {
457 inputToSspStateMap[mecState] += numberOfNonComponentStates;
462 std::vector<std::pair<uint_fast64_t, uint_fast64_t>> sspComponentExitChoicesToOriginalMap;
465 auto sspMatrixVector = buildSspMatrixVector(componentLraValues, inputToSspStateMap, statesNotInComponent, numberOfNonComponentStates,
466 this->isProduceSchedulerSet() ? &sspComponentExitChoicesToOriginalMap :
nullptr);
471 minMaxLinearEquationSolverFactory.getRequirements(env,
true,
true, this->getOptimizationDirection(),
false, this->isProduceSchedulerSet());
475 std::unique_ptr<storm::solver::MinMaxLinearEquationSolver<ValueType>> solver = minMaxLinearEquationSolverFactory.create(env, sspMatrixVector.first);
476 solver->setHasUniqueSolution();
477 solver->setHasNoEndComponents();
478 solver->setTrackScheduler(this->isProduceSchedulerSet());
479 auto lowerUpperBounds = std::minmax_element(componentLraValues.begin(), componentLraValues.end());
480 solver->setLowerBound(*lowerUpperBounds.first);
481 solver->setUpperBound(*lowerUpperBounds.second);
482 solver->setRequirementsChecked();
485 std::vector<ValueType> x(sspMatrixVector.first.getRowGroupCount());
486 solver->solveEquations(env, this->getOptimizationDirection(), x, sspMatrixVector.second);
489 if (this->isProduceSchedulerSet() && solver->hasScheduler()) {
491 constructOptimalChoices(solver->getSchedulerChoices(), sspMatrixVector.first, inputToSspStateMap, statesNotInComponent, numberOfNonComponentStates,
492 sspComponentExitChoicesToOriginalMap);
494 STORM_LOG_ERROR_COND(!this->isProduceSchedulerSet(),
"Requested to produce a scheduler, but no scheduler was generated.");
499 std::vector<ValueType> result = std::move(sspMatrixVector.second);
500 result.resize(this->_transitionMatrix.getRowGroupCount());
501 result.shrink_to_fit();
SolverEnvironment & solver()
storm::RationalNumber const & getAperiodicFactor() const
storm::solver::LraMethod const & getNondetLraMethod() const
bool const & isNondetLraMethodSetFromDefault() const
bool isForceExact() const
LongRunAverageSolverEnvironment & lra()
bool isForceSoundness() const
Expression rational(double value) const
Creates an expression that characterizes the given rational literal.
ExpressionManager const & getManager() const
Retrieves the manager responsible for this variable.
Helper class for model checking queries that depend on the long run behavior of the (nondeterministic...
Helper class for model checking queries that depend on the long run behavior of the (nondeterministic...
storm::storage::Scheduler< ValueType > extractScheduler() const
virtual void createDecomposition() override
virtual std::vector< ValueType > buildAndSolveSsp(Environment const &env, std::vector< ValueType > const &mecLraValues) override
std::pair< storm::storage::SparseMatrix< ValueType >, std::vector< ValueType > > buildSspMatrixVector(std::vector< ValueType > const &mecLraValues, std::vector< uint64_t > const &inputToSspStateMap, storm::storage::BitVector const &statesNotInComponent, uint64_t numberOfNonComponentStates, std::vector< std::pair< uint64_t, uint64_t > > *sspComponentExitChoicesToOriginalMap)
SparseInfiniteHorizonHelper< ValueType, true >::ValueGetter ValueGetter
Function mapping from indices to values.
std::vector< uint64_t > const & getProducedOptimalChoices() const
virtual ValueType computeLraForComponent(Environment const &env, ValueGetter const &stateValuesGetter, ValueGetter const &actionValuesGetter, storm::storage::MaximalEndComponent const &component) override
void constructOptimalChoices(std::vector< uint64_t > const &sspChoices, storm::storage::SparseMatrix< ValueType > const &sspMatrix, std::vector< uint64_t > const &inputToSspStateMap, storm::storage::BitVector const &statesNotInComponent, uint64_t numberOfNonComponentStates, std::vector< std::pair< uint64_t, uint64_t > > const &sspComponentExitChoicesToOriginalMap)
std::pair< bool, ValueType > computeLraForTrivialMec(Environment const &env, ValueGetter const &stateValuesGetter, ValueGetter const &actionValuesGetter, storm::storage::MaximalEndComponent const &mec)
SparseNondeterministicInfiniteHorizonHelper(storm::storage::SparseMatrix< ValueType > const &transitionMatrix)
Initializes the helper for a discrete time model (i.e.
ValueType computeLraForMecLp(Environment const &env, ValueGetter const &stateValuesGetter, ValueGetter const &actionValuesGetter, storm::storage::MaximalEndComponent const &mec)
As computeLraForMec but uses linear programming as a solution method (independent of what is set in e...
ValueType computeLraForMecVi(Environment const &env, ValueGetter const &stateValuesGetter, ValueGetter const &actionValuesGetter, storm::storage::MaximalEndComponent const &mec)
As computeLraForMec but uses value iteration as a solution method (independent of what is set in env)
Helper class that performs iterations of the value iteration method.
bool hasEnabledCriticalRequirement() const
std::string getEnabledRequirementsAsString() const
Returns a string that enumerates the enabled requirements.
A bit vector that is internally represented as a vector of 64-bit values.
void set(uint_fast64_t index, bool value=true)
Sets the given truth value at the given index.
This class represents a maximal end-component of a nondeterministic model.
iterator begin()
Retrieves an iterator that points to the first state and its choices in the MEC.
This class defines which action is chosen in a particular state of a non-deterministic model.
void setChoice(SchedulerChoice< ValueType > const &choice, uint_fast64_t modelState, uint_fast64_t memoryState=0)
Sets the choice defined by the scheduler for the given state.
index_type getNumberOfEntries() const
Retrieves the number of entries in the rows.
A class that can be used to build a sparse matrix by adding value by value.
void addNextValue(index_type row, index_type column, value_type const &value)
Sets the matrix entry at the given row and column to the given value.
void newRowGroup(index_type startingRow)
Starts a new row group in the matrix.
SparseMatrix< value_type > build(index_type overriddenRowCount=0, index_type overriddenColumnCount=0, index_type overriddenRowGroupCount=0)
A class that holds a possibly non-square matrix in the compressed row storage format.
const_rows getRow(index_type row) const
Returns an object representing the given row.
std::vector< index_type > const & getRowGroupIndices() const
Returns the grouping of rows of this matrix.
#define STORM_LOG_INFO(message)
#define STORM_LOG_ASSERT(cond, message)
#define STORM_LOG_ERROR_COND(cond, message)
#define STORM_LOG_THROW(cond, exception, message)
Expression sum(std::vector< storm::expressions::Expression > const &expressions)
bool componentElementChoicesContains(typename storm::storage::StronglyConnectedComponent::value_type const &element, uint64_t choice)
uint64_t getComponentElementState(typename storm::storage::StronglyConnectedComponent::value_type const &element)
Auxiliary functions that deal with the different kinds of components (MECs on potentially nondetermin...
uint64_t const * getComponentElementChoicesEnd(typename storm::storage::StronglyConnectedComponent::value_type const &element)
constexpr uint64_t getComponentElementChoiceCount(typename storm::storage::StronglyConnectedComponent::value_type const &)
@ NondetTsNoIs
deterministic choice at timed states, deterministic choice at instant states (as in Markov Automata w...
@ DetTsNondetIs
deterministic choice at timed states, no instant states (as in DTMCs and CTMCs)
uint64_t const * getComponentElementChoicesBegin(typename storm::storage::StronglyConnectedComponent::value_type const &element)
void addSspMatrixChoice(uint64_t const &inputMatrixChoice, storm::storage::SparseMatrix< ValueType > const &inputTransitionMatrix, std::vector< uint64_t > const &inputToSspStateMap, uint64_t const &numberOfNonComponentStates, uint64_t const ¤tSspChoice, storm::storage::SparseMatrixBuilder< ValueType > &sspMatrixBuilder)
Auxiliary function that adds the entries of the Ssp Matrix for a single choice (i....
void setVectorValues(std::vector< T > &vector, storm::storage::BitVector const &positions, std::vector< T > const &values)
Sets the provided values at the provided positions in the given vector.
void selectVectorValues(std::vector< T > &vector, storm::storage::BitVector const &positions, std::vector< T > const &values)
Selects the elements from a vector at the specified positions and writes them consecutively into anot...
bool isZero(ValueType const &a)