Browse Source

More steps towards integrating LRA in pcaa

tempestpy_adaptions
Tim Quatmann 4 years ago
parent
commit
360c3877b7
  1. 312
      src/storm/modelchecker/multiobjective/pcaa/StandardPcaaWeightVectorChecker.cpp
  2. 10
      src/storm/modelchecker/multiobjective/pcaa/StandardPcaaWeightVectorChecker.h

312
src/storm/modelchecker/multiobjective/pcaa/StandardPcaaWeightVectorChecker.cpp

@ -192,46 +192,60 @@ namespace storm {
}
template <typename ValueType>
std::vector<uint64_t> computeValidInitialScheduler(storm::storage::SparseMatrix<ValueType> const& matrix, storm::storage::BitVector const& rowsWithSumLessOne) {
std::vector<uint64_t> result(matrix.getRowGroupCount());
auto const& groups = matrix.getRowGroupIndices();
auto backwardsTransitions = matrix.transpose(true);
storm::storage::BitVector processedStates(result.size(), false);
for (uint64_t state = 0; state < result.size(); ++state) {
if (rowsWithSumLessOne.getNextSetIndex(groups[state]) < groups[state + 1]) {
result[state] = rowsWithSumLessOne.getNextSetIndex(groups[state]) - groups[state];
processedStates.set(state, true);
}
}
std::vector<uint64_t> stack(processedStates.begin(), processedStates.end());
void computeSchedulerProb1(storm::storage::SparseMatrix<ValueType> const& transitionMatrix, storm::storage::SparseMatrix<ValueType> const& backwardTransitions, storm::storage::BitVector const& consideredStates, storm::storage::BitVector& statesToReach, std::vector<uint64_t>& choices) {
std::vector<uint64_t> stack;
storm::storage::BitVector& processedStates = statesToReach;
stack.insert(stack.end(), processedStates.begin(), processedStates.end());
uint64_t currentState = 0;
while (!stack.empty()) {
uint64_t current = stack.back();
currentState = stack.back();
stack.pop_back();
STORM_LOG_ASSERT(processedStates.get(current), "states on the stack shall be processed.");
for (auto const& entry : backwardsTransitions.getRow(current)) {
uint64_t pred = entry.getColumn();
if (!processedStates.get(pred)) {
// Find a choice that leads to a processed state
uint64_t predChoice = groups[pred];
bool foundSuccessor = false;
for (; predChoice < groups[pred + 1]; ++predChoice) {
for (auto const& predEntry : matrix.getRow(predChoice)) {
if (processedStates.get(predEntry.getColumn())) {
foundSuccessor = true;
for (auto const& predecessorEntry : backwardTransitions.getRow(currentState)) {
auto predecessor = predecessorEntry.getColumn();
if (consideredStates.get(predecessor) & !processedStates.get(predecessor)) {
// Find a choice leading to an already processed state (such a choice has to exist since this is a predecessor of the currentState)
auto const& groupStart = transitionMatrix.getRowGroupIndices()[predecessor];
auto const& groupEnd = transitionMatrix.getRowGroupIndices()[predecessor + 1];
uint64_t row = groupStart;
for (; row < groupEnd; ++row) {
bool hasSuccessorInProcessedStates = false;
for (auto const& successorOfPredecessor : transitionMatrix.getRow(row)) {
if (processedStates.get(successorOfPredecessor.getColumn())) {
hasSuccessorInProcessedStates = true;
break;
}
}
if (foundSuccessor) {
if (hasSuccessorInProcessedStates) {
choices[predecessor] = row - groupStart;
processedStates.set(predecessor, true);
stack.push_back(predecessor);
break;
}
}
STORM_LOG_ASSERT(foundSuccessor && predChoice < groups[pred + 1], "Predecessor of a processed state should have a processed successor");
result[pred] = predChoice - groups[pred];
processedStates.set(pred, true);
stack.push_back(pred);
STORM_LOG_ASSERT(row < groupEnd, "Unable to find choice at a predecessor of a processed state that leads to a processed state.");
}
}
}
STORM_LOG_ASSERT(consideredStates.isSubsetOf(processedStates), "Not all states have been processed.");
}
template <typename ValueType>
std::vector<uint64_t> computeValidInitialScheduler(storm::storage::SparseMatrix<ValueType> const& matrix, storm::storage::BitVector const& rowsWithSumLessOne) {
std::vector<uint64_t> result(matrix.getRowGroupCount());
auto const& groups = matrix.getRowGroupIndices();
auto backwardsTransitions = matrix.transpose(true);
storm::storage::BitVector processedStates(result.size(), false);
for (uint64_t state = 0; state < result.size(); ++state) {
if (rowsWithSumLessOne.getNextSetIndex(groups[state]) < groups[state + 1]) {
result[state] = rowsWithSumLessOne.getNextSetIndex(groups[state]) - groups[state];
processedStates.set(state, true);
}
}
computeSchedulerProb1(matrix, backwardsTransitions, ~processedStates, processedStates, result);
return result;
}
@ -269,39 +283,7 @@ namespace storm {
// Finally, we need to take care of states that will reach a bad state with prob greater 0 (including the bad states themselves).
// due to the precondition, we know that it has to be possible to eventually avoid the bad states for ever.
// Perform a backwards search from the avoid states and store choices with prob. 1
std::vector<uint64_t> stack;
storm::storage::BitVector processedStates(avoidBadStates);
stack.insert(stack.end(), processedStates.begin(), processedStates.end());
uint64_t currentState = 0;
while (!stack.empty()) {
currentState = stack.back();
stack.pop_back();
for (auto predecessorEntryIt = backwardTransitions.begin(currentState), predecessorEntryIte = backwardTransitions.end(currentState); predecessorEntryIt != predecessorEntryIte; ++predecessorEntryIt) {
if (!processedStates.get(predecessorEntryIt->getColumn())) {
// Find a choice leading to an already processed state (such a choice has to exist since this is a predecessor of the currentState)
auto const& groupStart = transitionMatrix.getRowGroupIndices()[predecessorEntryIt->getColumn()];
auto const& groupEnd = transitionMatrix.getRowGroupIndices()[predecessorEntryIt->getColumn() + 1];
for (uint64_t row = groupStart; row < groupEnd; ++row) {
bool hasSuccessorInProcessedStates = false;
for (auto const& successorOfPredecessor : transitionMatrix.getRow(row)) {
if (processedStates.get(successorOfPredecessor.getColumn())) {
hasSuccessorInProcessedStates = true;
break;
}
}
if (hasSuccessorInProcessedStates) {
result[predecessorEntryIt->getColumn()] = row - groupStart;
processedStates.set(predecessorEntryIt->getColumn(), true);
stack.push_back(predecessorEntryIt->getColumn());
break;
}
}
}
}
}
STORM_LOG_ASSERT(processedStates.full(), "Not all states have been processed.");
computeSchedulerProb1(transitionMatrix, backwardTransitions, reachBadWithProbGreater0AStates, avoidBadStates, result);
return result;
}
@ -341,6 +323,7 @@ namespace storm {
// Set up the choice values
storm::utility::vector::selectVectorValues(ecQuotient->auxChoiceValues, ecQuotient->ecqToOriginalChoiceMapping, weightedRewardVector);
std::map<uint64_t, uint64_t> ecqStateToOptimalMecMap;
if (!lraObjectives.empty()) {
// We also need to assign a value for each ecQuotientChoice that corresponds to "staying" in the eliminated EC. (at this point these choices should all have a value of zero).
// Since each of the eliminated ECs has to contain *at least* one LRA EC, we need to find the largest value among the contained LRA ECs
@ -352,18 +335,20 @@ namespace storm {
uint64_t ecqChoice = ecQuotient->ecqStayInEcChoices.getNextSetIndex(ecQuotient->matrix.getRowGroupIndices()[ecqState]);
STORM_LOG_ASSERT(ecqChoice < ecQuotient->matrix.getRowGroupIndices()[ecqState + 1], "Unable to find choice that represents staying inside the (eliminated) ec.");
auto& ecqChoiceValue = ecQuotient->auxChoiceValues[ecqChoice];
if (foundEcqChoices.get(ecqChoice)) {
ecqChoiceValue = std::max(ecqChoiceValue, mecValue);
} else {
auto insertionRes = ecqStateToOptimalMecMap.emplace(ecqState, mecIndex);
if (insertionRes.second) {
// We have seen this ecqState for the first time.
STORM_LOG_ASSERT(storm::utility::isZero(ecqChoiceValue), "Expected a total reward of zero for choices that represent staying in an EC for ever.");
ecqChoiceValue = mecValue;
foundEcqChoices.set(ecqChoice, true);
} else {
if (mecValue > ecqChoiceValue) { // found a larger value
ecqChoiceValue = mecValue;
insertionRes.first->second = mecIndex;
}
}
}
}
// TODO: Continue to integrate LRA here (solver bounds need to be set correctly)
storm::solver::GeneralMinMaxLinearEquationSolverFactory<ValueType> solverFactory;
std::unique_ptr<storm::solver::MinMaxLinearEquationSolver<ValueType>> solver = solverFactory.create(env, ecQuotient->matrix);
solver->setTrackScheduler(true);
@ -390,7 +375,7 @@ namespace storm {
solver->solveEquations(env, ecQuotient->auxStateValues, ecQuotient->auxChoiceValues);
this->weightedResult = std::vector<ValueType>(transitionMatrix.getRowGroupCount());
transformReducedSolutionToOriginalModel(ecQuotient->matrix, ecQuotient->auxStateValues, solver->getSchedulerChoices(), ecQuotient->ecqToOriginalChoiceMapping, ecQuotient->originalToEcqStateMapping, this->weightedResult, this->optimalChoices);
transformEcqSolutionToOriginalModel(ecQuotient->auxStateValues, solver->getSchedulerChoices(), ecqStateToOptimalMecMap, this->weightedResult, this->optimalChoices);
}
template <class SparseModelType>
@ -534,11 +519,19 @@ namespace storm {
ecQuotient->matrix = std::move(ecElimResult.matrix);
ecQuotient->ecqToOriginalChoiceMapping = std::move(ecElimResult.newToOldRowMapping);
ecQuotient->originalToEcqStateMapping = std::move(ecElimResult.oldToNewStateMapping);
ecQuotient->ecqToOriginalStateMapping.resize(ecQuotient->matrix.getRowGroupCount());
for (uint64_t state = 0; state < ecQuotient->originalToEcqStateMapping.size(); ++state) {
uint64_t ecqState = ecQuotient->originalToEcqStateMapping[state];
if (ecqState < ecQuotient->matrix.getRowGroupCount()) {
ecQuotient->ecqToOriginalStateMapping[ecqState].insert(state);
}
}
ecQuotient->ecqStayInEcChoices = std::move(ecElimResult.sinkRows);
ecQuotient->origReward0Choices = std::move(newReward0Choices);
ecQuotient->rowsWithSumLessOne = std::move(rowsWithSumLessOne);
ecQuotient->auxStateValues.resize(ecQuotient->matrix.getRowGroupCount());
ecQuotient->auxChoiceValues.resize(ecQuotient->matrix.getRowCount());
}
}
@ -560,14 +553,26 @@ namespace storm {
template <class SparseModelType>
void StandardPcaaWeightVectorChecker<SparseModelType>::setBoundsToSolver(storm::solver::AbstractEquationSolver<ValueType>& solver, bool requiresLower, bool requiresUpper, std::vector<ValueType> const& weightVector, storm::storage::BitVector const& objectiveFilter, storm::storage::SparseMatrix<ValueType> const& transitions, storm::storage::BitVector const& rowsWithSumLessOne, std::vector<ValueType> const& rewards) const {
// Check whether bounds are already available
boost::optional<ValueType> lowerBound = this->computeWeightedResultBound(true, weightVector, objectiveFilter);
boost::optional<ValueType> lowerBound = this->computeWeightedResultBound(true, weightVector, objectiveFilter & ~lraObjectives);
if (lowerBound) {
if (!lraObjectives.empty()) {
auto min = std::min_element(lraMecDecomposition->auxMecValues.begin(), lraMecDecomposition->auxMecValues.end());
if (min != lraMecDecomposition->auxMecValues.end()) {
lowerBound.get() += *min;
}
}
solver.setLowerBound(lowerBound.get());
}
boost::optional<ValueType> upperBound = this->computeWeightedResultBound(false, weightVector, objectiveFilter);
if (upperBound) {
if (!lraObjectives.empty()) {
auto max = std::max_element(lraMecDecomposition->auxMecValues.begin(), lraMecDecomposition->auxMecValues.end());
if (max != lraMecDecomposition->auxMecValues.end()) {
upperBound.get() += *max;
}
}
solver.setUpperBound(upperBound.get());
}
@ -581,7 +586,7 @@ namespace storm {
// Compute the one step target probs
std::vector<ValueType> oneStepTargetProbs(transitions.getRowCount(), storm::utility::zero<ValueType>());
for (auto const& row : rowsWithSumLessOne) {
for (auto row : rowsWithSumLessOne) {
oneStepTargetProbs[row] = storm::utility::one<ValueType>() - transitions.getRowSum(row);
}
@ -629,101 +634,98 @@ namespace storm {
}
template <class SparseModelType>
void StandardPcaaWeightVectorChecker<SparseModelType>::transformReducedSolutionToOriginalModel(storm::storage::SparseMatrix<ValueType> const& reducedMatrix,
std::vector<ValueType> const& reducedSolution,
std::vector<uint_fast64_t> const& reducedOptimalChoices,
std::vector<uint_fast64_t> const& reducedToOriginalChoiceMapping,
std::vector<uint_fast64_t> const& originalToReducedStateMapping,
std::vector<ValueType>& originalSolution,
std::vector<uint_fast64_t>& originalOptimalChoices) const {
void StandardPcaaWeightVectorChecker<SparseModelType>::transformEcqSolutionToOriginalModel(std::vector<ValueType> const& ecqSolution,
std::vector<uint_fast64_t> const& ecqOptimalChoices,
std::map<uint64_t, uint64_t> const& ecqStateToOptimalMecMap,
std::vector<ValueType>& originalSolution,
std::vector<uint_fast64_t>& originalOptimalChoices) const {
storm::storage::BitVector bottomStates(transitionMatrix.getRowGroupCount(), false);
storm::storage::BitVector statesThatShouldStayInTheirEC(transitionMatrix.getRowGroupCount(), false);
storm::storage::BitVector statesWithUndefSched(transitionMatrix.getRowGroupCount(), false);
auto backwardsTransitions = transitionMatrix.transpose(true);
// Handle all the states for which the choice in the original model is uniquely given by the choice in the reduced model
// Also store some information regarding the remaining states
for (uint_fast64_t state = 0; state < transitionMatrix.getRowGroupCount(); ++state) {
// Check if the state exists in the reduced model, i.e., the mapping retrieves a valid index
uint_fast64_t stateInReducedModel = originalToReducedStateMapping[state];
if (stateInReducedModel < reducedMatrix.getRowGroupCount()) {
originalSolution[state] = reducedSolution[stateInReducedModel];
uint_fast64_t chosenRowInReducedModel = reducedMatrix.getRowGroupIndices()[stateInReducedModel] + reducedOptimalChoices[stateInReducedModel];
uint_fast64_t chosenRowInOriginalModel = reducedToOriginalChoiceMapping[chosenRowInReducedModel];
// Check if the state is a bottom state, i.e., the chosen row stays inside its EC.
bool stateIsBottom = reward0EStates.get(state);
for (auto const& entry : transitionMatrix.getRow(chosenRowInOriginalModel)) {
stateIsBottom &= originalToReducedStateMapping[entry.getColumn()] == stateInReducedModel;
}
if (stateIsBottom) {
bottomStates.set(state);
statesThatShouldStayInTheirEC.set(state);
// Keep track of states for which no choice has been set yet.
storm::storage::BitVector unprocessedStates(transitionMatrix.getRowGroupCount(), true);
// For each eliminated ec, keep track of the states (within the ec) that we want to reach and the states for which a choice needs to be set
// (Declared already at this point to avoid expensive allocations in each loop iteration)
storm::storage::BitVector ecStatesToReach(transitionMatrix.getRowGroupCount(), false);
storm::storage::BitVector ecStatesToProcess(transitionMatrix.getRowGroupCount(), false);
// Run through each state of the ec quotient as well as the associated state(s) of the original model
for (uint64_t ecqState = 0; ecqState < ecqSolution.size(); ++ecqState) {
uint64_t ecqChoice = ecQuotient->matrix.getRowGroupIndices()[ecqState] + ecqOptimalChoices[ecqState];
uint_fast64_t origChoice = ecQuotient->ecqToOriginalChoiceMapping[ecqChoice];
auto const& origStates = ecQuotient->ecqToOriginalStateMapping[ecqState];
STORM_LOG_ASSERT(!origStates.empty(), "Unexpected empty set of original states.");
bool origStatesNeedSchedulerComputation = false;
if (ecQuotient->ecqStayInEcChoices.get(origChoice) && !ecqStateToOptimalMecMap.empty()) {
// The current ecqState represents an elimnated EC and we need to stay in this EC and we need to make sure that optimal MEC decisions are performed within this EC.
STORM_LOG_ASSERT(ecqStateToOptimalMecMap.count(ecqState) > 0, "No Lra Mec associated to given eliminated EC");
auto const& lraMec = lraMecDecomposition->mecs[ecqStateToOptimalMecMap.at(ecqState)];
if (lraMec.size() == origStates.size()) {
// LRA mec and eliminated EC coincide
for (auto const& state : origStates) {
STORM_LOG_ASSERT(lraMec.containsState(state), "Expected state to be contained in the lra mec.");
// Note that the optimal choice for this state has already been set in the infinite horizon phase.
unprocessedStates.set(state, false);
originalSolution[state] = ecqSolution[ecqState];
}
} else {
// Check if the chosen row originaly belonged to the current state (and not to another state of the EC)
if (chosenRowInOriginalModel >= transitionMatrix.getRowGroupIndices()[state] &&
chosenRowInOriginalModel < transitionMatrix.getRowGroupIndices()[state+1]) {
originalOptimalChoices[state] = chosenRowInOriginalModel - transitionMatrix.getRowGroupIndices()[state];
} else {
statesWithUndefSched.set(state);
statesThatShouldStayInTheirEC.set(state);
// LRA mec is proper subset of eliminated ec. There are also other states for which we have to set choices leading to the LRA MEC inside.
STORM_LOG_ASSERT(lraMec.size() < origStates.size(), "Lra Mec should be a proper subset of the eliminated ec.");
origStatesNeedSchedulerComputation = true;
for (auto const& state : origStates) {
if (lraMec.containsState(state)) {
ecStatesToReach.set(state, true);
// Note that the optimal choice for this state has already been set in the infinite horizon phase.
} else {
ecStatesToProcess.set(state, true);
}
unprocessedStates.set(state, false);
originalSolution[state] = ecqSolution[ecqState];
}
}
} else {
// if the state does not exist in the reduced model, it means that the (weighted) result is always zero, independent of the scheduler.
originalSolution[state] = storm::utility::zero<ValueType>();
// However, it might be the case that infinite reward is induced for an objective with weight 0.
// To avoid this, all possible bottom states are made bottom and the remaining states have to reach a bottom state with prob. one
if (reward0EStates.get(state)) {
bottomStates.set(state);
// In this case, we can safely take the origChoice at the corresponding state (say 's').
// For all other origStates associated with ecqState (if there are any), we make sure that the state 's' is reached almost surely.
if (origStates.size() > 1) {
origStatesNeedSchedulerComputation = true;
for (auto const& state : origStates) {
// Check if the orig choice originates from this state
auto groupStart = transitionMatrix.getRowGroupIndices()[state];
auto groupEnd = transitionMatrix.getRowGroupIndices()[state + 1];
if (origChoice >= groupStart && origChoice < groupEnd) {
originalOptimalChoices[state] = origChoice - groupStart;
ecStatesToReach.set(state, true);
} else {
STORM_LOG_ASSERT(origStates.size() > 1, "Multiple original states expected.");
ecStatesToProcess.set(state, true);
}
unprocessedStates.set(state, false);
originalSolution[state] = ecqSolution[ecqState];
}
} else {
statesWithUndefSched.set(state);
// There is just one state so we take the associated choice.
auto state = *origStates.begin();
auto groupStart = transitionMatrix.getRowGroupIndices()[state];
originalOptimalChoices[state] = origChoice - groupStart;
STORM_LOG_ASSERT(originalOptimalChoices[state] > 0 && originalOptimalChoices[state] < transitionMatrix.getRowGroupSize(state), "Invalid choice.");
originalSolution[state] = ecqSolution[ecqState];
unprocessedStates.set(state, false);
}
}
}
// Handle bottom states
for (auto state : bottomStates) {
bool foundRowForState = false;
// Find a row with zero rewards that only leads to bottom states.
// If the state should stay in its EC, we also need to make sure that all successors map to the same state in the reduced model
uint_fast64_t stateInReducedModel = originalToReducedStateMapping[state];
for (uint_fast64_t row = transitionMatrix.getRowGroupIndices()[state]; row < transitionMatrix.getRowGroupIndices()[state+1]; ++row) {
bool rowOnlyLeadsToBottomStates = true;
bool rowStaysInEC = true;
for ( auto const& entry : transitionMatrix.getRow(row)) {
rowOnlyLeadsToBottomStates &= bottomStates.get(entry.getColumn());
rowStaysInEC &= originalToReducedStateMapping[entry.getColumn()] == stateInReducedModel;
}
if (rowOnlyLeadsToBottomStates && (rowStaysInEC || !statesThatShouldStayInTheirEC.get(state)) && actionsWithoutRewardInUnboundedPhase.get(row)) {
foundRowForState = true;
originalOptimalChoices[state] = row - transitionMatrix.getRowGroupIndices()[state];
break;
}
if (origStatesNeedSchedulerComputation) {
computeSchedulerProb1(transitionMatrix, backwardsTransitions, ecStatesToProcess, ecStatesToReach, originalOptimalChoices);
// These states have only been altered if the origStatesNeedSchedulerComputation flag has been set to true. Hence, they only need to be cleared in this if-branch.
ecStatesToProcess.clear();
ecStatesToReach.clear();
}
STORM_LOG_ASSERT(foundRowForState, "Could not find a suitable choice for a bottom state.");
}
// Handle remaining states with still undef. scheduler (either EC states or non-subsystem states)
while(!statesWithUndefSched.empty()) {
for (auto state : statesWithUndefSched) {
// Iteratively Try to find a choice such that at least one successor has a defined scheduler.
uint_fast64_t stateInReducedModel = originalToReducedStateMapping[state];
for (uint_fast64_t row = transitionMatrix.getRowGroupIndices()[state]; row < transitionMatrix.getRowGroupIndices()[state+1]; ++row) {
bool rowStaysInEC = true;
bool rowLeadsToDefinedScheduler = false;
for (auto const& entry : transitionMatrix.getRow(row)) {
rowStaysInEC &= ( stateInReducedModel == originalToReducedStateMapping[entry.getColumn()]);
rowLeadsToDefinedScheduler |= !statesWithUndefSched.get(entry.getColumn());
}
if (rowLeadsToDefinedScheduler && (rowStaysInEC || !statesThatShouldStayInTheirEC.get(state))) {
originalOptimalChoices[state] = row - transitionMatrix.getRowGroupIndices()[state];
statesWithUndefSched.set(state, false);
break;
}
}
}
}
// The states that still not have been processed, there is no associated state of the ec quotient.
// This is because the value for these states will be 0 under all (lra optimal-) schedulers.
// We thus do not change the choices at these states (to remain LRA optimality) and just set the result to 0
storm::utility::vector::setVectorValues(originalSolution, unprocessedStates, storm::utility::zero<ValueType>());
}

10
src/storm/modelchecker/multiobjective/pcaa/StandardPcaaWeightVectorChecker.h

@ -12,6 +12,7 @@
#include "storm/modelchecker/multiobjective/pcaa/PcaaWeightVectorChecker.h"
#include "storm/modelchecker/multiobjective/preprocessing/SparseMultiObjectivePreprocessorResult.h"
#include "storm/utility/vector.h"
#include "storm/storage/BoostTypes.h"
namespace storm {
namespace modelchecker {
@ -106,11 +107,9 @@ namespace storm {
/*!
* Transforms the results of a min-max-solver that considers a reduced model (without end components) to a result for the original (unreduced) model
*/
void transformReducedSolutionToOriginalModel(storm::storage::SparseMatrix<ValueType> const& reducedMatrix,
std::vector<ValueType> const& reducedSolution,
std::vector<uint_fast64_t> const& reducedOptimalChoices,
std::vector<uint_fast64_t> const& reducedToOriginalChoiceMapping,
std::vector<uint_fast64_t> const& originalToReducedStateMapping,
void transformEcqSolutionToOriginalModel(std::vector<ValueType> const& ecqSolution,
std::vector<uint_fast64_t> const& ecqOptimalChoices,
std::map<uint64_t, uint64_t> const& ecqStateToOptimalMecMap,
std::vector<ValueType>& originalSolution,
std::vector<uint_fast64_t>& originalOptimalChoices) const;
@ -155,6 +154,7 @@ namespace storm {
storm::storage::SparseMatrix<ValueType> matrix;
std::vector<uint_fast64_t> ecqToOriginalChoiceMapping;
std::vector<uint_fast64_t> originalToEcqStateMapping;
std::vector<storm::storage::FlatSetStateContainer> ecqToOriginalStateMapping;
storm::storage::BitVector ecqStayInEcChoices;
storm::storage::BitVector origReward0Choices;
storm::storage::BitVector rowsWithSumLessOne;

Loading…
Cancel
Save