Browse Source

first version of long-run-average for parametric DTMCs

Former-commit-id: 85253c4f05
main
dehnert 10 years ago
parent
commit
2a5780d5be
  1. 304
      src/modelchecker/reachability/SparseDtmcEliminationModelChecker.cpp
  2. 9
      src/modelchecker/reachability/SparseDtmcEliminationModelChecker.h

304
src/modelchecker/reachability/SparseDtmcEliminationModelChecker.cpp

@ -105,16 +105,23 @@ namespace storm {
return true; return true;
} }
} else if (formula.isReachabilityRewardFormula()) { } else if (formula.isReachabilityRewardFormula()) {
storm::logic::ReachabilityRewardFormula reachabilityRewardFormula = formula.asReachabilityRewardFormula(); storm::logic::ReachabilityRewardFormula const& reachabilityRewardFormula = formula.asReachabilityRewardFormula();
if (reachabilityRewardFormula.getSubformula().isPropositionalFormula()) { if (reachabilityRewardFormula.getSubformula().isPropositionalFormula()) {
return true; return true;
} }
} else if (formula.isConditionalPathFormula()) { } else if (formula.isConditionalPathFormula()) {
storm::logic::ConditionalPathFormula conditionalPathFormula = formula.asConditionalPathFormula(); storm::logic::ConditionalPathFormula const& conditionalPathFormula = formula.asConditionalPathFormula();
if (conditionalPathFormula.getLeftSubformula().isEventuallyFormula() && conditionalPathFormula.getRightSubformula().isEventuallyFormula()) { if (conditionalPathFormula.getLeftSubformula().isEventuallyFormula() && conditionalPathFormula.getRightSubformula().isEventuallyFormula()) {
return this->canHandle(conditionalPathFormula.getLeftSubformula()) && this->canHandle(conditionalPathFormula.getRightSubformula()); return this->canHandle(conditionalPathFormula.getLeftSubformula()) && this->canHandle(conditionalPathFormula.getRightSubformula());
} }
} else if (formula.isPropositionalFormula()) { } else if (formula.isLongRunAverageOperatorFormula()) {
storm::logic::LongRunAverageOperatorFormula const& longRunAverageOperatorFormula = formula.asLongRunAverageOperatorFormula();
if (longRunAverageOperatorFormula.getSubformula().isPropositionalFormula()) {
return true;
}
}
else if (formula.isPropositionalFormula()) {
return true; return true;
} }
return false; return false;
@ -156,32 +163,20 @@ namespace storm {
if (furtherComputationNeeded) { if (furtherComputationNeeded) {
// Start by decomposing the DTMC into its BSCCs. // Start by decomposing the DTMC into its BSCCs.
storm::storage::StronglyConnectedComponentDecomposition<ValueType> bsccDecomposition(transitionMatrix, storm::storage::BitVector(transitionMatrix.getRowCount(), true), false, true); storm::storage::BitVector allStates(transitionMatrix.getRowCount(), true);
storm::storage::StronglyConnectedComponentDecomposition<ValueType> bsccDecomposition(transitionMatrix, allStates, false, true);
storm::storage::BitVector statesInRelevantBsccs(maybeStates.size()); if (computeResultsForInitialStatesOnly) {
for (auto const& scc : bsccDecomposition) { // Determine the set of states that is reachable from the initial state without jumping over a target state.
// Since all states in an SCC can reach all other states, we only need to check whether an arbitrary storm::storage::BitVector reachableStates = storm::utility::graph::getReachableStates(transitionMatrix, initialStates, allStates, allStates, true);
// state is a maybe state. // Subtract from the maybe states the set of states that is not reachable (on a path from the initial to a target state).
if (maybeStates.get(*scc.begin())) { maybeStates &= reachableStates;
for (auto const& state : scc) {
statesInRelevantBsccs.set(state, true);
}
}
} }
// Determine the set of initial states of the sub-model. std::vector<ValueType> stateValues(maybeStates.size(), storm::utility::zero<ValueType>());
storm::storage::BitVector newInitialStates = initialStates % maybeStates; storm::utility::vector::setVectorValues(stateValues, psiStates, storm::utility::one<ValueType>());
storm::storage::BitVector newPsiStates = psiStates % maybeStates; result = computeLongRunValues(transitionMatrix, this->getModel().getBackwardTransitions(), initialStates, maybeStates, bsccDecomposition, computeResultsForInitialStatesOnly, stateValues);
storm::storage::BitVector newStatesInRelevantBsccs = statesInRelevantBsccs % maybeStates;
// We then build the submatrix that only has the transitions of the maybe states.
storm::storage::SparseMatrix<ValueType> submatrix = transitionMatrix.getSubmatrix(true, maybeStates, maybeStates);
storm::storage::SparseMatrix<ValueType> submatrixTransposed = submatrix.transpose();
std::vector<ValueType> stateValues(maybeStates.getNumberOfSetBits(), storm::utility::zero<ValueType>());
storm::utility::vector::setVectorValues(stateValues, newPsiStates, storm::utility::one<ValueType>());
std::vector<ValueType> subresult = computeLongRunValues(submatrix, submatrixTransposed, newInitialStates, newStatesInRelevantBsccs, computeResultsForInitialStatesOnly, stateValues);
storm::utility::vector::setVectorValues<ValueType>(result, maybeStates, subresult);
} }
// Construct check result based on whether we have computed values for all states or just the initial states. // Construct check result based on whether we have computed values for all states or just the initial states.
@ -195,92 +190,151 @@ namespace storm {
} }
template<typename SparseDtmcModelType> template<typename SparseDtmcModelType>
std::vector<typename SparseDtmcEliminationModelChecker<SparseDtmcModelType>::ValueType> SparseDtmcEliminationModelChecker<SparseDtmcModelType>::computeLongRunValues(storm::storage::SparseMatrix<ValueType> const& transitionMatrix, storm::storage::SparseMatrix<ValueType> const& backwardTransitions, storm::storage::BitVector const& initialStates, storm::storage::BitVector const& statesInBsccs, bool computeResultsForInitialStatesOnly, std::vector<ValueType>& stateValues) { std::vector<typename SparseDtmcEliminationModelChecker<SparseDtmcModelType>::ValueType> SparseDtmcEliminationModelChecker<SparseDtmcModelType>::computeLongRunValues(storm::storage::SparseMatrix<ValueType> const& transitionMatrix, storm::storage::SparseMatrix<ValueType> const& backwardTransitions, storm::storage::BitVector const& initialStates, storm::storage::BitVector const& statesWithProbabilityGreater0, storm::storage::StronglyConnectedComponentDecomposition<ValueType> const& bsccDecomposition, bool computeResultsForInitialStatesOnly, std::vector<ValueType>& stateValues) {
std::chrono::high_resolution_clock::time_point totalTimeStart = std::chrono::high_resolution_clock::now();
// std::chrono::high_resolution_clock::time_point totalTimeStart = std::chrono::high_resolution_clock::now(); std::chrono::high_resolution_clock::time_point conversionStart = std::chrono::high_resolution_clock::now();
// // Then, we convert the reduced matrix to a more flexible format to be able to perform state elimination more easily.
// std::chrono::high_resolution_clock::time_point conversionStart = std::chrono::high_resolution_clock::now(); FlexibleSparseMatrix flexibleMatrix = getFlexibleSparseMatrix(transitionMatrix);
// // Then, we convert the reduced matrix to a more flexible format to be able to perform state elimination more easily. flexibleMatrix.filter(statesWithProbabilityGreater0, statesWithProbabilityGreater0);
// FlexibleSparseMatrix flexibleMatrix = getFlexibleSparseMatrix(transitionMatrix); FlexibleSparseMatrix flexibleBackwardTransitions = getFlexibleSparseMatrix(backwardTransitions);
// FlexibleSparseMatrix flexibleBackwardTransitions = getFlexibleSparseMatrix(backwardTransitions); flexibleBackwardTransitions.filter(statesWithProbabilityGreater0, statesWithProbabilityGreater0);
// auto conversionEnd = std::chrono::high_resolution_clock::now(); auto conversionEnd = std::chrono::high_resolution_clock::now();
// std::chrono::high_resolution_clock::time_point modelCheckingStart = std::chrono::high_resolution_clock::now();
// std::chrono::high_resolution_clock::time_point modelCheckingStart = std::chrono::high_resolution_clock::now(); storm::settings::modules::SparseDtmcEliminationModelCheckerSettings::EliminationOrder order = storm::settings::sparseDtmcEliminationModelCheckerSettings().getEliminationOrder();
// boost::optional<std::vector<uint_fast64_t>> distanceBasedPriorities;
// storm::settings::modules::SparseDtmcEliminationModelCheckerSettings::EliminationOrder order = storm::settings::sparseDtmcEliminationModelCheckerSettings().getEliminationOrder(); if (eliminationOrderNeedsDistances(order)) {
// boost::optional<std::vector<uint_fast64_t>> distanceBasedPriorities; distanceBasedPriorities = getDistanceBasedPriorities(transitionMatrix, backwardTransitions, initialStates, stateValues,
// if (eliminationOrderNeedsDistances(order)) { eliminationOrderNeedsForwardDistances(order), eliminationOrderNeedsReversedDistances(order));
// distanceBasedPriorities = getDistanceBasedPriorities(transitionMatrix, backwardTransitions, initialStates, stateValues, }
// eliminationOrderNeedsForwardDistances(order), eliminationOrderNeedsReversedDistances(order)); uint_fast64_t numberOfStates = transitionMatrix.getRowCount();
// } storm::storage::BitVector statesInBsccs(numberOfStates);
// storm::storage::BitVector relevantBsccs(bsccDecomposition.size());
// // Compute the average time to stay in each state for all states in BSCCs. storm::storage::BitVector bsccRepresentativesAsBitVector(numberOfStates);
// boost::optional<std::vector<ValueType>> averageTimeInStates = std::vector<ValueType>(stateValues.size(), storm::utility::zero<ValueType>()); std::vector<storm::storage::sparse::state_type> bsccRepresentatives;
// for (auto state : statesInBsccs) { uint_fast64_t currentIndex = 0;
// for (auto const& entry : transitionMatrix.getRow(state)) { for (auto const& bscc : bsccDecomposition) {
// bool hasSelfLoop = false; // Since all states in an SCC can reach all other states, we only need to check whether an arbitrary
// if (entry.getColumn() == state) { // state is a maybe state.
// hasSelfLoop = true; if (statesWithProbabilityGreater0.get(*bscc.cbegin())) {
// averageTimeInStates.get()[state] = storm::utility::one<ValueType>() / (storm::utility::one<ValueType>() - entry.getValue()); relevantBsccs.set(currentIndex);
// } bsccRepresentatives.push_back(*bscc.cbegin());
// if (!hasSelfLoop) { bsccRepresentativesAsBitVector.set(*bscc.cbegin(), true);
// averageTimeInStates.get()[state] = storm::utility::one<ValueType>(); for (auto const& state : bscc) {
// } statesInBsccs.set(state, true);
// } }
// } }
// ++currentIndex;
// performOrdinaryStateElimination(flexibleMatrix, flexibleBackwardTransitions, statesInBsccs, initialStates, computeResultsForInitialStatesOnly, stateValues, averageTimeInStates, distanceBasedPriorities); }
// statesInBsccs &= ~bsccRepresentativesAsBitVector;
// // Compute the average time to stay in each state for all states in BSCCs.
// std::vector<ValueType> averageTimeInStates(stateValues.size(), storm::utility::one<ValueType>());
// // First, we eliminate all states in BSCCs (except for the representative states).
// {
// std::unique_ptr<StatePriorityQueue> priorityQueue = createStatePriorityQueue(distanceBasedPriorities, flexibleMatrix, flexibleBackwardTransitions, stateValues, statesInBsccs);
// // Create a bit vector that represents the subsystem of states we still have to eliminate. ValueUpdateCallback valueUpdateCallback = [&stateValues,&averageTimeInStates] (storm::storage::sparse::state_type const& state, ValueType const& loopProbability) {
// storm::storage::BitVector subsystem = storm::storage::BitVector(transitionMatrix.getRowCount(), true); stateValues[state] = storm::utility::simplify(loopProbability * stateValues[state]);
// averageTimeInStates[state] = storm::utility::simplify(loopProbability * averageTimeInStates[state]);
// uint_fast64_t maximalDepth = 0; };
// if (storm::settings::sparseDtmcEliminationModelCheckerSettings().getEliminationMethod() == storm::settings::modules::SparseDtmcEliminationModelCheckerSettings::EliminationMethod::State) { PredecessorUpdateCallback predecessorCallback = [&stateValues,&averageTimeInStates] (storm::storage::sparse::state_type const& predecessor, ValueType const& probability, storm::storage::sparse::state_type const& state) {
// performOrdinaryStateElimination(flexibleMatrix, flexibleBackwardTransitions, subsystem, initialStates, computeResultsForInitialStatesOnly, values, distanceBasedPriorities); stateValues[predecessor] = storm::utility::simplify(stateValues[predecessor] + storm::utility::simplify(probability * stateValues[state]));
// } else if (storm::settings::sparseDtmcEliminationModelCheckerSettings().getEliminationMethod() == storm::settings::modules::SparseDtmcEliminationModelCheckerSettings::EliminationMethod::Hybrid) { averageTimeInStates[predecessor] = storm::utility::simplify(averageTimeInStates[predecessor] + storm::utility::simplify(probability * averageTimeInStates[state]));
// maximalDepth = performHybridStateElimination(transitionMatrix, flexibleMatrix, flexibleBackwardTransitions, subsystem, initialStates, computeResultsForInitialStatesOnly, values, distanceBasedPriorities); };
// } boost::optional<PriorityUpdateCallback> priorityUpdateCallback = PriorityUpdateCallback([&flexibleMatrix,&flexibleBackwardTransitions,&stateValues,&priorityQueue] (storm::storage::sparse::state_type const& state) {
// priorityQueue->update(state, flexibleMatrix, flexibleBackwardTransitions, stateValues);
// STORM_LOG_ASSERT(flexibleMatrix.empty(), "Not all transitions were eliminated."); });
// STORM_LOG_ASSERT(flexibleBackwardTransitions.empty(), "Not all transitions were eliminated."); boost::optional<PredecessorFilterCallback> predecessorFilterCallback = boost::none;
// while (priorityQueue->hasNextState()) {
// std::chrono::high_resolution_clock::time_point modelCheckingEnd = std::chrono::high_resolution_clock::now(); storm::storage::sparse::state_type state = priorityQueue->popNextState();
// std::chrono::high_resolution_clock::time_point totalTimeEnd = std::chrono::high_resolution_clock::now(); eliminateState(state, flexibleMatrix, flexibleBackwardTransitions, valueUpdateCallback, predecessorCallback, priorityUpdateCallback, predecessorFilterCallback, true);
// STORM_LOG_ASSERT(checkConsistent(flexibleMatrix, flexibleBackwardTransitions), "The forward and backward transition matrices became inconsistent.");
// if (storm::settings::generalSettings().isShowStatisticsSet()) { }
// std::chrono::high_resolution_clock::duration conversionTime = conversionEnd - conversionStart; }
// std::chrono::milliseconds conversionTimeInMilliseconds = std::chrono::duration_cast<std::chrono::milliseconds>(conversionTime); // Now, we set the values of all states in BSCCs to that of the representative value (and clear the
// std::chrono::high_resolution_clock::duration modelCheckingTime = modelCheckingEnd - modelCheckingStart; // transitions of the representative states while doing so).
// std::chrono::milliseconds modelCheckingTimeInMilliseconds = std::chrono::duration_cast<std::chrono::milliseconds>(modelCheckingTime); auto representativeIt = bsccRepresentatives.begin();
// std::chrono::high_resolution_clock::duration totalTime = totalTimeEnd - totalTimeStart; for (auto sccIndex : relevantBsccs) {
// std::chrono::milliseconds totalTimeInMilliseconds = std::chrono::duration_cast<std::chrono::milliseconds>(totalTime); // We only need to set the values for all states of the BSCC if we are not computing the values for the
// // initial states only.
// STORM_PRINT_AND_LOG(std::endl); ValueType bsccValue = stateValues[*representativeIt] / averageTimeInStates[*representativeIt];
// STORM_PRINT_AND_LOG("Time breakdown:" << std::endl); auto const& bscc = bsccDecomposition[sccIndex];
// STORM_PRINT_AND_LOG(" * time for conversion: " << conversionTimeInMilliseconds.count() << "ms" << std::endl); if (!computeResultsForInitialStatesOnly) {
// STORM_PRINT_AND_LOG(" * time for checking: " << modelCheckingTimeInMilliseconds.count() << "ms" << std::endl); for (auto const& state : bscc) {
// STORM_PRINT_AND_LOG("------------------------------------------" << std::endl); stateValues[state] = bsccValue;
// STORM_PRINT_AND_LOG(" * total time: " << totalTimeInMilliseconds.count() << "ms" << std::endl); }
// STORM_PRINT_AND_LOG(std::endl); } else {
// STORM_PRINT_AND_LOG("Other:" << std::endl); for (auto const& state : bscc) {
// STORM_PRINT_AND_LOG(" * number of states eliminated: " << transitionMatrix.getRowCount() << std::endl); stateValues[state] = storm::utility::zero<ValueType>();
// if (storm::settings::sparseDtmcEliminationModelCheckerSettings().getEliminationMethod() == storm::settings::modules::SparseDtmcEliminationModelCheckerSettings::EliminationMethod::Hybrid) { }
// STORM_PRINT_AND_LOG(" * maximal depth of SCC decomposition: " << maximalDepth << std::endl); stateValues[*representativeIt] = bsccValue;
// } }
// } typename SparseDtmcEliminationModelChecker<SparseDtmcModelType>::FlexibleSparseMatrix::row_type& representativeForwardRow = flexibleMatrix.getRow(*representativeIt);
// representativeForwardRow.clear();
// // Now, we return the value for the only initial state. representativeForwardRow.shrink_to_fit();
// STORM_LOG_DEBUG("Simplifying and returning result."); typename SparseDtmcEliminationModelChecker<SparseDtmcModelType>::FlexibleSparseMatrix::row_type& representativeBackwardRow = flexibleBackwardTransitions.getRow(*representativeIt);
// for (auto& value : values) { auto it = representativeBackwardRow.begin(), ite = representativeBackwardRow.end();
// value = storm::utility::simplify(value); for (; it != ite; ++it) {
// } if (it->getColumn() == *representativeIt) {
// return values; break;
return std::vector<ValueType>(); }
}
representativeBackwardRow.erase(it);
++representativeIt;
}
// If there are states remaining that are not in BSCCs, we need to eliminate them now.
storm::storage::BitVector remainingStates = statesWithProbabilityGreater0 & ~statesInBsccs;
// Reset the values of all remaining non-BSCC states (they might have been erroneously changed by the
// previous state elimination.
for (auto state : remainingStates) {
if (!bsccRepresentativesAsBitVector.get(state)) {
stateValues[state] = storm::utility::zero<ValueType>();
}
}
performOrdinaryStateElimination(flexibleMatrix, flexibleBackwardTransitions, remainingStates, initialStates, computeResultsForInitialStatesOnly, stateValues, distanceBasedPriorities);
std::chrono::high_resolution_clock::time_point modelCheckingEnd = std::chrono::high_resolution_clock::now();
std::chrono::high_resolution_clock::time_point totalTimeEnd = std::chrono::high_resolution_clock::now();
if (storm::settings::generalSettings().isShowStatisticsSet()) {
std::chrono::high_resolution_clock::duration conversionTime = conversionEnd - conversionStart;
std::chrono::milliseconds conversionTimeInMilliseconds = std::chrono::duration_cast<std::chrono::milliseconds>(conversionTime);
std::chrono::high_resolution_clock::duration modelCheckingTime = modelCheckingEnd - modelCheckingStart;
std::chrono::milliseconds modelCheckingTimeInMilliseconds = std::chrono::duration_cast<std::chrono::milliseconds>(modelCheckingTime);
std::chrono::high_resolution_clock::duration totalTime = totalTimeEnd - totalTimeStart;
std::chrono::milliseconds totalTimeInMilliseconds = std::chrono::duration_cast<std::chrono::milliseconds>(totalTime);
STORM_PRINT_AND_LOG(std::endl);
STORM_PRINT_AND_LOG("Time breakdown:" << std::endl);
STORM_PRINT_AND_LOG(" * time for conversion: " << conversionTimeInMilliseconds.count() << "ms" << std::endl);
STORM_PRINT_AND_LOG(" * time for checking: " << modelCheckingTimeInMilliseconds.count() << "ms" << std::endl);
STORM_PRINT_AND_LOG("------------------------------------------" << std::endl);
STORM_PRINT_AND_LOG(" * total time: " << totalTimeInMilliseconds.count() << "ms" << std::endl);
}
// Now, we return the value for the only initial state.
STORM_LOG_DEBUG("Simplifying and returning result.");
for (auto& value : stateValues) {
value = storm::utility::simplify(value);
}
return stateValues;
} }
template<typename SparseDtmcModelType> template<typename SparseDtmcModelType>
@ -708,8 +762,6 @@ namespace storm {
} }
} }
flexibleMatrix.print();
ValueType numerator = storm::utility::zero<ValueType>(); ValueType numerator = storm::utility::zero<ValueType>();
ValueType denominator = storm::utility::zero<ValueType>(); ValueType denominator = storm::utility::zero<ValueType>();
@ -825,8 +877,9 @@ namespace storm {
while (priorityQueue->hasNextState()) { while (priorityQueue->hasNextState()) {
storm::storage::sparse::state_type state = priorityQueue->popNextState(); storm::storage::sparse::state_type state = priorityQueue->popNextState();
eliminateState(state, transitionMatrix, backwardTransitions, valueUpdateCallback, predecessorCallback, priorityUpdateCallback, predecessorFilterCallback, computeResultsForInitialStatesOnly && !initialStates.get(state)); bool removeForwardTransitions = computeResultsForInitialStatesOnly && !initialStates.get(state);
if (computeResultsForInitialStatesOnly && !initialStates.get(state)) { eliminateState(state, transitionMatrix, backwardTransitions, valueUpdateCallback, predecessorCallback, priorityUpdateCallback, predecessorFilterCallback, removeForwardTransitions);
if (removeForwardTransitions) {
values[state] = storm::utility::zero<ValueType>(); values[state] = storm::utility::zero<ValueType>();
} }
STORM_LOG_ASSERT(checkConsistent(transitionMatrix, backwardTransitions), "The forward and backward transition matrices became inconsistent."); STORM_LOG_ASSERT(checkConsistent(transitionMatrix, backwardTransitions), "The forward and backward transition matrices became inconsistent.");
@ -1350,6 +1403,25 @@ namespace storm {
return true; return true;
} }
template<typename SparseDtmcModelType>
void SparseDtmcEliminationModelChecker<SparseDtmcModelType>::FlexibleSparseMatrix::filter(storm::storage::BitVector const& rowFilter, storm::storage::BitVector const& columnFilter) {
for (uint_fast64_t rowIndex = 0; rowIndex < this->data.size(); ++rowIndex) {
auto& row = this->data[rowIndex];
if (!rowFilter.get(rowIndex)) {
row.clear();
row.shrink_to_fit();
continue;
}
row_type newRow;
for (auto const& element : row) {
if (columnFilter.get(element.getColumn())) {
newRow.push_back(element);
}
}
row = std::move(newRow);
}
}
template<typename SparseDtmcModelType> template<typename SparseDtmcModelType>
typename SparseDtmcEliminationModelChecker<SparseDtmcModelType>::FlexibleSparseMatrix SparseDtmcEliminationModelChecker<SparseDtmcModelType>::getFlexibleSparseMatrix(storm::storage::SparseMatrix<ValueType> const& matrix, bool setAllValuesToOne) { typename SparseDtmcEliminationModelChecker<SparseDtmcModelType>::FlexibleSparseMatrix SparseDtmcEliminationModelChecker<SparseDtmcModelType>::getFlexibleSparseMatrix(storm::storage::SparseMatrix<ValueType> const& matrix, bool setAllValuesToOne) {
FlexibleSparseMatrix flexibleMatrix(matrix.getRowCount()); FlexibleSparseMatrix flexibleMatrix(matrix.getRowCount());

9
src/modelchecker/reachability/SparseDtmcEliminationModelChecker.h

@ -7,6 +7,11 @@
#include "src/utility/constants.h" #include "src/utility/constants.h"
namespace storm { namespace storm {
namespace storage {
template<typename ValueType>
class StronglyConnectedComponentDecomposition;
}
namespace modelchecker { namespace modelchecker {
template<typename SparseDtmcModelType> template<typename SparseDtmcModelType>
@ -54,6 +59,8 @@ namespace storm {
bool empty() const; bool empty() const;
void filter(storm::storage::BitVector const& rowFilter, storm::storage::BitVector const& columnFilter);
/*! /*!
* Checks whether the given state has a self-loop with an arbitrary probability in the probability matrix. * Checks whether the given state has a self-loop with an arbitrary probability in the probability matrix.
* *
@ -110,7 +117,7 @@ namespace storm {
PenaltyFunctionType penaltyFunction; PenaltyFunctionType penaltyFunction;
}; };
static std::vector<ValueType> computeLongRunValues(storm::storage::SparseMatrix<ValueType> const& transitionMatrix, storm::storage::SparseMatrix<ValueType> const& backwardTransitions, storm::storage::BitVector const& initialStates, storm::storage::BitVector const& statesInBsccs, bool computeResultsForInitialStatesOnly, std::vector<ValueType>& stateValues); static std::vector<ValueType> computeLongRunValues(storm::storage::SparseMatrix<ValueType> const& transitionMatrix, storm::storage::SparseMatrix<ValueType> const& backwardTransitions, storm::storage::BitVector const& initialStates, storm::storage::BitVector const& statesWithProbabilityGreater0, storm::storage::StronglyConnectedComponentDecomposition<ValueType> const& bsccDecomposition, bool computeResultsForInitialStatesOnly, std::vector<ValueType>& stateValues);
static std::vector<ValueType> computeReachabilityValues(storm::storage::SparseMatrix<ValueType> const& transitionMatrix, std::vector<ValueType>& values, storm::storage::SparseMatrix<ValueType> const& backwardTransitions, storm::storage::BitVector const& initialStates, bool computeResultsForInitialStatesOnly, storm::storage::BitVector const& phiStates, storm::storage::BitVector const& psiStates, std::vector<ValueType> const& oneStepProbabilitiesToTarget); static std::vector<ValueType> computeReachabilityValues(storm::storage::SparseMatrix<ValueType> const& transitionMatrix, std::vector<ValueType>& values, storm::storage::SparseMatrix<ValueType> const& backwardTransitions, storm::storage::BitVector const& initialStates, bool computeResultsForInitialStatesOnly, storm::storage::BitVector const& phiStates, storm::storage::BitVector const& psiStates, std::vector<ValueType> const& oneStepProbabilitiesToTarget);

|||||||
100:0
Loading…
Cancel
Save