Browse Source

Merge remote-tracking branch 'origin/multi-objective' into multi-objective

tempestpy_adaptions
TimQu 8 years ago
parent
commit
1797a63757
  1. 5
      resources/3rdparty/CMakeLists.txt
  2. 4
      src/storm/cli/cli.cpp
  3. 90
      src/storm/modelchecker/multiobjective/pcaa/SparseMaPcaaWeightVectorChecker.cpp
  4. 4
      src/storm/modelchecker/multiobjective/pcaa/SparseMaPcaaWeightVectorChecker.h
  5. 44
      src/storm/modelchecker/multiobjective/pcaa/SparsePcaaQuantitativeQuery.cpp
  6. 81
      src/storm/modelchecker/multiobjective/pcaa/SparsePcaaWeightVectorChecker.cpp

5
resources/3rdparty/CMakeLists.txt

@ -43,8 +43,8 @@ list(APPEND STORM_DEP_TARGETS gmm)
##
#############################################################
add_imported_library_interface(Eigen33 "${PROJECT_SOURCE_DIR}/resources/3rdparty/eigen-3.3-beta1")
list(APPEND STORM_DEP_TARGETS Eigen33)
# add_imported_library_interface(Eigen33 "${PROJECT_SOURCE_DIR}/resources/3rdparty/eigen-3.3-beta1")
# list(APPEND STORM_DEP_TARGETS Eigen33)
#############################################################
@ -589,3 +589,4 @@ if(ENABLE_CUDA)
endif()
add_custom_target(copy_resources_headers DEPENDS ${CMAKE_BINARY_DIR}/include/resources/3rdparty/sparsepp/sparsepp.h ${CMAKE_BINARY_DIR}/include/resources/3rdparty/sparsepp/sparsepp.h)

4
src/storm/cli/cli.cpp

@ -135,14 +135,14 @@ namespace storm {
getrusage(RUSAGE_SELF, &ru);
std::cout << "Performance statistics:" << std::endl;
std::cout << " * peak memory usage: " << ru.ru_maxrss/1024/1024 << " mb" << std::endl;
std::cout << " * peak memory usage: " << ru.ru_maxrss/1024 << " mb" << std::endl;
std::cout << " * CPU time: " << ru.ru_utime.tv_sec << "." << std::setw(3) << std::setfill('0') << ru.ru_utime.tv_usec/1000 << " seconds" << std::endl;
if (wallclockMilliseconds != 0) {
std::cout << " * wallclock time: " << (wallclockMilliseconds/1000) << "." << std::setw(3) << std::setfill('0') << (wallclockMilliseconds % 1000) << " seconds" << std::endl;
}
std::cout << "STATISTICS_OVERALL_HEADERS;" << "memory;CPU time;wallclock time;" << std::endl;
std::cout << "STATISTICS_OVERALL_DATA;"
<< ru.ru_maxrss/1024/1024 << ";"
<< ru.ru_maxrss/1024 << ";"
<< ru.ru_utime.tv_sec << "." << std::setw(3) << std::setfill('0') << ru.ru_utime.tv_usec/1000 << ";"
<< (wallclockMilliseconds/1000) << "." << std::setw(3) << std::setfill('0') << (wallclockMilliseconds % 1000) << ";" << std::endl;
#else

90
src/storm/modelchecker/multiobjective/pcaa/SparseMaPcaaWeightVectorChecker.cpp

@ -75,12 +75,12 @@ namespace storm {
updateDataToCurrentEpoch(MS, PS, *minMax, consideredObjectives, currentEpoch, weightVector, lowerTimeBoundIt, lowerTimeBounds, upperTimeBoundIt, upperTimeBounds);
// Compute the values that can be obtained at probabilistic states in the current time epoch
performPSStep(PS, MS, *minMax, *linEq, optimalChoicesAtCurrentEpoch, consideredObjectives);
performPSStep(PS, MS, *minMax, *linEq, optimalChoicesAtCurrentEpoch, consideredObjectives, weightVector);
// Compute values that can be obtained at Markovian states after letting one (digitized) time unit pass.
// Only perform such a step if there is time left.
if(currentEpoch>0) {
performMSStep(MS, PS, consideredObjectives);
performMSStep(MS, PS, consideredObjectives, weightVector);
--currentEpoch;
} else {
break;
@ -355,55 +355,73 @@ namespace storm {
}
template <class SparseMaModelType>
void SparseMaPcaaWeightVectorChecker<SparseMaModelType>::performPSStep(SubModel& PS, SubModel const& MS, MinMaxSolverData& minMax, LinEqSolverData& linEq, std::vector<uint_fast64_t>& optimalChoicesAtCurrentEpoch, storm::storage::BitVector const& consideredObjectives) const {
void SparseMaPcaaWeightVectorChecker<SparseMaModelType>::performPSStep(SubModel& PS, SubModel const& MS, MinMaxSolverData& minMax, LinEqSolverData& linEq, std::vector<uint_fast64_t>& optimalChoicesAtCurrentEpoch, storm::storage::BitVector const& consideredObjectives, std::vector<ValueType> const& weightVector) const {
// compute a choice vector for the probabilistic states that is optimal w.r.t. the weighted reward vector
minMax.solver->solveEquations(PS.weightedSolutionVector, minMax.b);
auto newScheduler = minMax.solver->getScheduler();
// check whether the linEqSolver needs to be updated, i.e., whether the scheduler has changed
if(linEq.solver == nullptr || newScheduler->getChoices() != optimalChoicesAtCurrentEpoch) {
if(consideredObjectives.getNumberOfSetBits() == 1 && !storm::utility::isZero(weightVector[*consideredObjectives.begin()])) {
// In this case there is no need to perform the computation on the individual objectives
optimalChoicesAtCurrentEpoch = newScheduler->getChoices();
linEq.solver = nullptr;
storm::storage::SparseMatrix<ValueType> linEqMatrix = PS.toPS.selectRowsFromRowGroups(optimalChoicesAtCurrentEpoch, true);
linEqMatrix.convertToEquationSystem();
linEq.solver = linEq.factory.create(std::move(linEqMatrix));
linEq.solver->setCachingEnabled(true);
}
// Get the results for the individual objectives.
// Note that we do not consider an estimate for each objective (as done in the unbounded phase) since the results from the previous epoch are already pretty close
for(auto objIndex : consideredObjectives) {
auto const& objectiveRewardVectorPS = PS.objectiveRewardVectors[objIndex];
auto const& objectiveSolutionVectorMS = MS.objectiveSolutionVectors[objIndex];
// compute rhs of equation system, i.e., PS.toMS * x + Rewards
// To safe some time, only do this for the obtained optimal choices
auto itGroupIndex = PS.toPS.getRowGroupIndices().begin();
auto itChoiceOffset = optimalChoicesAtCurrentEpoch.begin();
for(auto& bValue : linEq.b) {
uint_fast64_t row = (*itGroupIndex) + (*itChoiceOffset);
bValue = objectiveRewardVectorPS[row];
for(auto const& entry : PS.toMS.getRow(row)){
bValue += entry.getValue() * objectiveSolutionVectorMS[entry.getColumn()];
auto objIndex = *consideredObjectives.begin();
PS.objectiveSolutionVectors[objIndex] = PS.weightedSolutionVector;
if(!storm::utility::isOne(weightVector[objIndex])) {
storm::utility::vector::scaleVectorInPlace(PS.objectiveSolutionVectors[objIndex], storm::utility::one<ValueType>()/weightVector[objIndex]);
}
} else {
// check whether the linEqSolver needs to be updated, i.e., whether the scheduler has changed
if(linEq.solver == nullptr || newScheduler->getChoices() != optimalChoicesAtCurrentEpoch) {
optimalChoicesAtCurrentEpoch = newScheduler->getChoices();
linEq.solver = nullptr;
storm::storage::SparseMatrix<ValueType> linEqMatrix = PS.toPS.selectRowsFromRowGroups(optimalChoicesAtCurrentEpoch, true);
linEqMatrix.convertToEquationSystem();
linEq.solver = linEq.factory.create(std::move(linEqMatrix));
linEq.solver->setCachingEnabled(true);
}
// Get the results for the individual objectives.
// Note that we do not consider an estimate for each objective (as done in the unbounded phase) since the results from the previous epoch are already pretty close
for(auto objIndex : consideredObjectives) {
auto const& objectiveRewardVectorPS = PS.objectiveRewardVectors[objIndex];
auto const& objectiveSolutionVectorMS = MS.objectiveSolutionVectors[objIndex];
// compute rhs of equation system, i.e., PS.toMS * x + Rewards
// To safe some time, only do this for the obtained optimal choices
auto itGroupIndex = PS.toPS.getRowGroupIndices().begin();
auto itChoiceOffset = optimalChoicesAtCurrentEpoch.begin();
for(auto& bValue : linEq.b) {
uint_fast64_t row = (*itGroupIndex) + (*itChoiceOffset);
bValue = objectiveRewardVectorPS[row];
for(auto const& entry : PS.toMS.getRow(row)){
bValue += entry.getValue() * objectiveSolutionVectorMS[entry.getColumn()];
}
++itGroupIndex;
++itChoiceOffset;
}
++itGroupIndex;
++itChoiceOffset;
linEq.solver->solveEquations(PS.objectiveSolutionVectors[objIndex], linEq.b);
}
linEq.solver->solveEquations(PS.objectiveSolutionVectors[objIndex], linEq.b);
}
}
template <class SparseMaModelType>
void SparseMaPcaaWeightVectorChecker<SparseMaModelType>::performMSStep(SubModel& MS, SubModel const& PS, storm::storage::BitVector const& consideredObjectives) const {
void SparseMaPcaaWeightVectorChecker<SparseMaModelType>::performMSStep(SubModel& MS, SubModel const& PS, storm::storage::BitVector const& consideredObjectives, std::vector<ValueType> const& weightVector) const {
MS.toMS.multiplyWithVector(MS.weightedSolutionVector, MS.auxChoiceValues);
storm::utility::vector::addVectors(MS.weightedRewardVector, MS.auxChoiceValues, MS.weightedSolutionVector);
MS.toPS.multiplyWithVector(PS.weightedSolutionVector, MS.auxChoiceValues);
storm::utility::vector::addVectors(MS.weightedSolutionVector, MS.auxChoiceValues, MS.weightedSolutionVector);
for(auto objIndex : consideredObjectives) {
MS.toMS.multiplyWithVector(MS.objectiveSolutionVectors[objIndex], MS.auxChoiceValues);
storm::utility::vector::addVectors(MS.objectiveRewardVectors[objIndex], MS.auxChoiceValues, MS.objectiveSolutionVectors[objIndex]);
MS.toPS.multiplyWithVector(PS.objectiveSolutionVectors[objIndex], MS.auxChoiceValues);
storm::utility::vector::addVectors(MS.objectiveSolutionVectors[objIndex], MS.auxChoiceValues, MS.objectiveSolutionVectors[objIndex]);
if(consideredObjectives.getNumberOfSetBits() == 1 && !storm::utility::isZero(weightVector[*consideredObjectives.begin()])) {
// In this case there is no need to perform the computation on the individual objectives
auto objIndex = *consideredObjectives.begin();
MS.objectiveSolutionVectors[objIndex] = MS.weightedSolutionVector;
if(!storm::utility::isOne(weightVector[objIndex])) {
storm::utility::vector::scaleVectorInPlace(MS.objectiveSolutionVectors[objIndex], storm::utility::one<ValueType>()/weightVector[objIndex]);
}
} else {
for(auto objIndex : consideredObjectives) {
MS.toMS.multiplyWithVector(MS.objectiveSolutionVectors[objIndex], MS.auxChoiceValues);
storm::utility::vector::addVectors(MS.objectiveRewardVectors[objIndex], MS.auxChoiceValues, MS.objectiveSolutionVectors[objIndex]);
MS.toPS.multiplyWithVector(PS.objectiveSolutionVectors[objIndex], MS.auxChoiceValues);
storm::utility::vector::addVectors(MS.objectiveSolutionVectors[objIndex], MS.auxChoiceValues, MS.objectiveSolutionVectors[objIndex]);
}
}
}

4
src/storm/modelchecker/multiobjective/pcaa/SparseMaPcaaWeightVectorChecker.h

@ -138,7 +138,7 @@ namespace storm {
*
* The resulting values represent the rewards at probabilistic states that are obtained at the current time epoch.
*/
void performPSStep(SubModel& PS, SubModel const& MS, MinMaxSolverData& minMax, LinEqSolverData& linEq, std::vector<uint_fast64_t>& optimalChoicesAtCurrentEpoch, storm::storage::BitVector const& consideredObjectives) const;
void performPSStep(SubModel& PS, SubModel const& MS, MinMaxSolverData& minMax, LinEqSolverData& linEq, std::vector<uint_fast64_t>& optimalChoicesAtCurrentEpoch, storm::storage::BitVector const& consideredObjectives, std::vector<ValueType> const& weightVector) const;
/*
* Performs a step for the Markovian states, that is
@ -146,7 +146,7 @@ namespace storm {
*
* The resulting values represent the rewards at Markovian states that are obtained after one (digitized) time unit has passed.
*/
void performMSStep(SubModel& MS, SubModel const& PS, storm::storage::BitVector const& consideredObjectives) const;
void performMSStep(SubModel& MS, SubModel const& PS, storm::storage::BitVector const& consideredObjectives, std::vector<ValueType> const& weightVector) const;
};

44
src/storm/modelchecker/multiobjective/pcaa/SparsePcaaQuantitativeQuery.cpp

@ -51,8 +51,6 @@ namespace storm {
template <class SparseModelType, typename GeometryValueType>
std::unique_ptr<CheckResult> SparsePcaaQuantitativeQuery<SparseModelType, GeometryValueType>::check() {
// First find one solution that achieves the given thresholds ...
if(this->checkAchievability()) {
// ... then improve it
@ -73,21 +71,26 @@ namespace storm {
template <class SparseModelType, typename GeometryValueType>
bool SparsePcaaQuantitativeQuery<SparseModelType, GeometryValueType>::checkAchievability() {
// We don't care for the optimizing objective at this point
this->diracWeightVectorsToBeChecked.set(indexOfOptimizingObjective, false);
if(this->objectives.size()>1) {
// We don't care for the optimizing objective at this point
this->diracWeightVectorsToBeChecked.set(indexOfOptimizingObjective, false);
while(!this->maxStepsPerformed()){
WeightVector separatingVector = this->findSeparatingVector(thresholds);
this->updateWeightedPrecisionInAchievabilityPhase(separatingVector);
this->performRefinementStep(std::move(separatingVector));
//Pick the threshold for the optimizing objective low enough so valid solutions are not excluded
thresholds[indexOfOptimizingObjective] = std::min(thresholds[indexOfOptimizingObjective], this->refinementSteps.back().lowerBoundPoint[indexOfOptimizingObjective]);
if(!checkIfThresholdsAreSatisfied(this->overApproximation)){
return false;
}
if(checkIfThresholdsAreSatisfied(this->underApproximation)){
return true;
while(!this->maxStepsPerformed()){
WeightVector separatingVector = this->findSeparatingVector(thresholds);
this->updateWeightedPrecisionInAchievabilityPhase(separatingVector);
this->performRefinementStep(std::move(separatingVector));
//Pick the threshold for the optimizing objective low enough so valid solutions are not excluded
thresholds[indexOfOptimizingObjective] = std::min(thresholds[indexOfOptimizingObjective], this->refinementSteps.back().lowerBoundPoint[indexOfOptimizingObjective]);
if(!checkIfThresholdsAreSatisfied(this->overApproximation)){
return false;
}
if(checkIfThresholdsAreSatisfied(this->underApproximation)){
return true;
}
}
} else {
// If there is only one objective than its the optimizing one. Thus the query has to be achievable.
return true;
}
STORM_LOG_ERROR("Could not check whether thresholds are achievable: Exceeded maximum number of refinement steps");
return false;
@ -124,6 +127,12 @@ namespace storm {
// thresholds) and (with very low probability) a scheduler that satisfies all (possibly strict) thresholds.
GeometryValueType result = storm::utility::zero<GeometryValueType>();
while(!this->maxStepsPerformed()) {
if(this->refinementSteps.empty()) {
// We did not make any refinement steps during the checkAchievability phase (e.g., because there is only one objective).
this->weightVectorChecker->setWeightedPrecision(storm::utility::convertNumber<typename SparseModelType::ValueType>(storm::settings::getModule<storm::settings::modules::MultiObjectiveSettings>().getPrecision()));
WeightVector separatingVector = directionOfOptimizingObjective;
this->performRefinementStep(std::move(separatingVector));
}
std::pair<Point, bool> optimizationRes = this->underApproximation->intersection(thresholdsAsPolytope)->optimize(directionOfOptimizingObjective);
STORM_LOG_THROW(optimizationRes.second, storm::exceptions::UnexpectedException, "The underapproximation is either unbounded or empty.");
result = optimizationRes.first[indexOfOptimizingObjective];
@ -155,11 +164,12 @@ namespace storm {
void SparsePcaaQuantitativeQuery<SparseModelType, GeometryValueType>::updateWeightedPrecisionInImprovingPhase(WeightVector const& weights) {
STORM_LOG_THROW(!storm::utility::isZero(weights[this->indexOfOptimizingObjective]), exceptions::UnexpectedException, "The chosen weight-vector gives zero weight for the objective that is to be optimized.");
// If weighs[indexOfOptimizingObjective] is low, the computation of the weightVectorChecker needs to be more precise.
// Our heuristic ensures that if p is the new vertex of the under-approximation, then max{ eps | p' = p + (0..0 eps 0..0) is in the over-approximation } <= multiobjective_precision/2
// Our heuristic ensures that if p is the new vertex of the under-approximation, then max{ eps | p' = p + (0..0 eps 0..0) is in the over-approximation } <= multiobjective_precision/0.9
GeometryValueType weightedPrecision = weights[this->indexOfOptimizingObjective] * storm::utility::convertNumber<GeometryValueType>(storm::settings::getModule<storm::settings::modules::MultiObjectiveSettings>().getPrecision());
// Normalize by division with the Euclidean Norm of the weight-vector
weightedPrecision /= storm::utility::sqrt(storm::utility::vector::dotProduct(weights, weights));
weightedPrecision /= GeometryValueType(2);
weightedPrecision *= storm::utility::convertNumber<GeometryValueType>(0.9);
this->weightVectorChecker->setWeightedPrecision(storm::utility::convertNumber<typename SparseModelType::ValueType>(weightedPrecision));
}

81
src/storm/modelchecker/multiobjective/pcaa/SparsePcaaWeightVectorChecker.cpp

@ -159,49 +159,62 @@ namespace storm {
template <class SparseModelType>
void SparsePcaaWeightVectorChecker<SparseModelType>::unboundedIndividualPhase(std::vector<ValueType> const& weightVector) {
storm::storage::SparseMatrix<ValueType> deterministicMatrix = model.getTransitionMatrix().selectRowsFromRowGroups(this->scheduler.getChoices(), true);
storm::storage::SparseMatrix<ValueType> deterministicBackwardTransitions = deterministicMatrix.transpose();
std::vector<ValueType> deterministicStateRewards(deterministicMatrix.getRowCount());
storm::solver::GeneralLinearEquationSolverFactory<ValueType> linearEquationSolverFactory;
// We compute an estimate for the results of the individual objectives which is obtained from the weighted result and the result of the objectives computed so far.
// Note that weightedResult = Sum_{i=1}^{n} w_i * objectiveResult_i.
std::vector<ValueType> weightedSumOfUncheckedObjectives = weightedResult;
ValueType sumOfWeightsOfUncheckedObjectives = storm::utility::vector::sum_if(weightVector, objectivesWithNoUpperTimeBound);
for(uint_fast64_t const& objIndex : storm::utility::vector::getSortedIndices(weightVector)) {
if(objectivesWithNoUpperTimeBound.get(objIndex)){
offsetsToLowerBound[objIndex] = storm::utility::zero<ValueType>();
offsetsToUpperBound[objIndex] = storm::utility::zero<ValueType>();
storm::utility::vector::selectVectorValues(deterministicStateRewards, this->scheduler.getChoices(), model.getTransitionMatrix().getRowGroupIndices(), discreteActionRewards[objIndex]);
storm::storage::BitVector statesWithRewards = ~storm::utility::vector::filterZero(deterministicStateRewards);
// As target states, we pick the states from which no reward is reachable.
storm::storage::BitVector targetStates = ~storm::utility::graph::performProbGreater0(deterministicBackwardTransitions, storm::storage::BitVector(deterministicMatrix.getRowCount(), true), statesWithRewards);
// Compute the estimate for this objective
if(!storm::utility::isZero(weightVector[objIndex])) {
objectiveResults[objIndex] = weightedSumOfUncheckedObjectives;
storm::utility::vector::scaleVectorInPlace(objectiveResults[objIndex], storm::utility::one<ValueType>() / sumOfWeightsOfUncheckedObjectives);
if(objectivesWithNoUpperTimeBound.getNumberOfSetBits()==1) {
auto objIndex = *objectivesWithNoUpperTimeBound.begin();
objectiveResults[objIndex] = weightedResult;
if(!storm::utility::isZero(weightVector[objIndex])) {
storm::utility::vector::scaleVectorInPlace(objectiveResults[objIndex], storm::utility::one<ValueType>()/weightVector[objIndex]);
}
for (uint_fast64_t objIndex2 = 0; objIndex2 < objectives.size(); ++objIndex2) {
if(objIndex != objIndex2) {
objectiveResults[objIndex2] = std::vector<ValueType>(model.getNumberOfStates(), storm::utility::zero<ValueType>());
}
}
} else {
storm::storage::SparseMatrix<ValueType> deterministicMatrix = model.getTransitionMatrix().selectRowsFromRowGroups(this->scheduler.getChoices(), true);
storm::storage::SparseMatrix<ValueType> deterministicBackwardTransitions = deterministicMatrix.transpose();
std::vector<ValueType> deterministicStateRewards(deterministicMatrix.getRowCount());
storm::solver::GeneralLinearEquationSolverFactory<ValueType> linearEquationSolverFactory;
// We compute an estimate for the results of the individual objectives which is obtained from the weighted result and the result of the objectives computed so far.
// Note that weightedResult = Sum_{i=1}^{n} w_i * objectiveResult_i.
std::vector<ValueType> weightedSumOfUncheckedObjectives = weightedResult;
ValueType sumOfWeightsOfUncheckedObjectives = storm::utility::vector::sum_if(weightVector, objectivesWithNoUpperTimeBound);
for(uint_fast64_t const& objIndex : storm::utility::vector::getSortedIndices(weightVector)) {
if(objectivesWithNoUpperTimeBound.get(objIndex)){
offsetsToLowerBound[objIndex] = storm::utility::zero<ValueType>();
offsetsToUpperBound[objIndex] = storm::utility::zero<ValueType>();
storm::utility::vector::selectVectorValues(deterministicStateRewards, this->scheduler.getChoices(), model.getTransitionMatrix().getRowGroupIndices(), discreteActionRewards[objIndex]);
storm::storage::BitVector statesWithRewards = ~storm::utility::vector::filterZero(deterministicStateRewards);
// As target states, we pick the states from which no reward is reachable.
storm::storage::BitVector targetStates = ~storm::utility::graph::performProbGreater0(deterministicBackwardTransitions, storm::storage::BitVector(deterministicMatrix.getRowCount(), true), statesWithRewards);
// Compute the estimate for this objective
if(!storm::utility::isZero(weightVector[objIndex])) {
objectiveResults[objIndex] = weightedSumOfUncheckedObjectives;
storm::utility::vector::scaleVectorInPlace(objectiveResults[objIndex], storm::utility::one<ValueType>() / sumOfWeightsOfUncheckedObjectives);
}
// Make sure that the objectiveResult is initialized in some way
objectiveResults[objIndex].resize(model.getNumberOfStates(), storm::utility::zero<ValueType>());
// Make sure that the objectiveResult is initialized in some way
objectiveResults[objIndex].resize(model.getNumberOfStates(), storm::utility::zero<ValueType>());
// Invoke the linear equation solver
objectiveResults[objIndex] = storm::modelchecker::helper::SparseDtmcPrctlHelper<ValueType>::computeReachabilityRewards(deterministicMatrix,
// Invoke the linear equation solver
objectiveResults[objIndex] = storm::modelchecker::helper::SparseDtmcPrctlHelper<ValueType>::computeReachabilityRewards(deterministicMatrix,
deterministicBackwardTransitions,
deterministicStateRewards,
targetStates,
false, //no qualitative checking,
linearEquationSolverFactory,
objectiveResults[objIndex]);
// Update the estimate for the next objectives.
if(!storm::utility::isZero(weightVector[objIndex])) {
storm::utility::vector::addScaledVector(weightedSumOfUncheckedObjectives, objectiveResults[objIndex], -weightVector[objIndex]);
sumOfWeightsOfUncheckedObjectives -= weightVector[objIndex];
// Update the estimate for the next objectives.
if(!storm::utility::isZero(weightVector[objIndex])) {
storm::utility::vector::addScaledVector(weightedSumOfUncheckedObjectives, objectiveResults[objIndex], -weightVector[objIndex]);
sumOfWeightsOfUncheckedObjectives -= weightVector[objIndex];
}
} else {
objectiveResults[objIndex] = std::vector<ValueType>(model.getNumberOfStates(), storm::utility::zero<ValueType>());
}
} else {
objectiveResults[objIndex] = std::vector<ValueType>(model.getNumberOfStates(), storm::utility::zero<ValueType>());
}
}
}

Loading…
Cancel
Save