Browse Source

Finalized hybrid CTMC model checker.

Former-commit-id: c217e11b06
tempestpy_adaptions
dehnert 10 years ago
parent
commit
be66ef2751
  1. 8
      src/builder/DdPrismModelBuilder.cpp
  2. 29
      src/builder/ExplicitPrismModelBuilder.cpp
  3. 229
      src/modelchecker/csl/HybridCtmcCslModelChecker.cpp
  4. 39
      src/modelchecker/csl/HybridCtmcCslModelChecker.h
  5. 24
      src/modelchecker/csl/SparseCtmcCslModelChecker.cpp
  6. 6
      src/settings/Argument.h
  7. 15
      src/storage/prism/Program.cpp
  8. 12
      src/utility/cli.h
  9. 11
      test/functional/modelchecker/GmmxxCtmcCslModelCheckerTest.cpp
  10. 100
      test/functional/modelchecker/GmmxxHybridCtmcCslModelCheckerTest.cpp
  11. 30
      test/functional/modelchecker/SparseCtmcCslModelCheckerTest.cpp

8
src/builder/DdPrismModelBuilder.cpp

@ -562,10 +562,12 @@ namespace storm {
}
storm::dd::Add<Type> transitionRewardDd = synchronization * states * rewards;
if (generationInfo.program.getModelType() == storm::prism::Program::ModelType::MDP) {
transitionRewardDd = transitions.notZero().toAdd() * transitionRewardDd;
} else {
if (generationInfo.program.getModelType() == storm::prism::Program::ModelType::DTMC) {
// For DTMCs we need to keep the weighting for the scaling that follows.
transitionRewardDd = transitions * transitionRewardDd;
} else {
// For all other model types, we do not scale the rewards.
transitionRewardDd = transitions.notZero().toAdd() * transitionRewardDd;
}
// Perform some sanity checks.

29
src/builder/ExplicitPrismModelBuilder.cpp

@ -127,21 +127,6 @@ namespace storm {
#endif
}
storm::prism::RewardModel rewardModel = storm::prism::RewardModel();
// Select the appropriate reward model.
if (options.buildRewards) {
// If a specific reward model was selected or one with the empty name exists, select it.
if (options.rewardModelName) {
rewardModel = preparedProgram.getRewardModel(options.rewardModelName.get());
} else if (preparedProgram.hasRewardModel("")) {
rewardModel = preparedProgram.getRewardModel("");
} else if (preparedProgram.hasRewardModel()) {
// Otherwise, we select the first one.
rewardModel = preparedProgram.getRewardModel(0);
}
}
// If the set of labels we are supposed to built is restricted, we need to remove the other labels from the program.
if (options.labelsToBuild) {
preparedProgram.filterLabels(options.labelsToBuild.get());
@ -162,6 +147,20 @@ namespace storm {
// Now that the program is fixed, we we need to substitute all constants with their concrete value.
preparedProgram = preparedProgram.substituteConstants();
// Select the appropriate reward model (after the constants have been substituted).
storm::prism::RewardModel rewardModel = storm::prism::RewardModel();
if (options.buildRewards) {
// If a specific reward model was selected or one with the empty name exists, select it.
if (options.rewardModelName) {
rewardModel = preparedProgram.getRewardModel(options.rewardModelName.get());
} else if (preparedProgram.hasRewardModel("")) {
rewardModel = preparedProgram.getRewardModel("");
} else if (preparedProgram.hasRewardModel()) {
// Otherwise, we select the first one.
rewardModel = preparedProgram.getRewardModel(0);
}
}
ModelComponents modelComponents = buildModelComponents(preparedProgram, rewardModel, options);

229
src/modelchecker/csl/HybridCtmcCslModelChecker.cpp

@ -169,11 +169,8 @@ namespace storm {
// staying in phi states.
std::unique_ptr<CheckResult> unboundedResult = HybridDtmcPrctlModelChecker<DdType, ValueType>::computeUntilProbabilitiesHelper(this->getModel(), this->computeProbabilityMatrix(this->getModel(), this->getModel().getTransitionMatrix(), this->getModel().getExitRateVector()), phiStates, psiStates, qualitative, *this->linearEquationSolverFactory);
// Determine the set of states that achieved a strictly positive probability.
std::unique_ptr<CheckResult> statesWithValuesGreaterZero = unboundedResult->asQuantitativeCheckResult().compareAgainstBound(storm::logic::ComparisonType::Greater, storm::utility::zero<ValueType>());
// And use it to compute the set of relevant states.
storm::dd::Bdd<DdType> relevantStates = storm::utility::graph::performProbGreater0(this->getModel(), this->getModel().getTransitionMatrix().notZero(), phiStates, statesWithValuesGreaterZero->asSymbolicQualitativeCheckResult<DdType>().getTruthValuesVector());
// Compute the set of relevant states.
storm::dd::Bdd<DdType> relevantStates = statesWithProbabilityGreater0 && phiStates;
// Filter the unbounded result such that it only contains values for the relevant states.
unboundedResult->filter(SymbolicQualitativeCheckResult<DdType>(this->getModel().getReachableStates(), relevantStates));
@ -204,28 +201,21 @@ namespace storm {
} else {
// In this case, the interval is of the form [t, t'] with t != 0 and t' != inf.
// Prepare some variables that are used by the two following blocks.
storm::dd::Bdd<DdType> relevantStates;
ValueType uniformizationRate = 0;
storm::dd::Add<DdType> uniformizedMatrix;
std::vector<ValueType> newSubresult;
storm::dd::Odd<DdType> odd;
if (lowerBound == upperBound) {
relevantStates = statesWithProbabilityGreater0;
} else {
if (lowerBound != upperBound) {
// In this case, the interval is of the form [t, t'] with t != 0, t' != inf and t != t'.
// Find the maximal rate of all 'maybe' states to take it as the uniformization rate.
ValueType uniformizationRate = 1.02 * (statesWithProbabilityGreater0NonPsi.toAdd() * exitRates).getMax();
STORM_LOG_THROW(uniformizationRate > 0, storm::exceptions::InvalidStateException, "The uniformization rate must be positive.");
// Compute the (first) uniformized matrix.
uniformizedMatrix = this->computeUniformizedMatrix(this->getModel(), this->getModel().getTransitionMatrix(), exitRates, statesWithProbabilityGreater0NonPsi, uniformizationRate);
storm::dd::Add<DdType> uniformizedMatrix = this->computeUniformizedMatrix(this->getModel(), this->getModel().getTransitionMatrix(), exitRates, statesWithProbabilityGreater0NonPsi, uniformizationRate);
// Create the one-step vector.
storm::dd::Add<DdType> b = (statesWithProbabilityGreater0NonPsi.toAdd() * this->getModel().getTransitionMatrix() * psiStates.swapVariables(this->getModel().getRowColumnMetaVariablePairs()).toAdd()).sumAbstract(this->getModel().getColumnVariables()) / this->getModel().getManager().getConstant(uniformizationRate);
// Build an ODD for the relevant states and translate the symbolic parts to their explicit representation.
odd = storm::dd::Odd<DdType>(statesWithProbabilityGreater0NonPsi);
storm::dd::Odd<DdType> odd = storm::dd::Odd<DdType>(statesWithProbabilityGreater0NonPsi);
storm::storage::SparseMatrix<ValueType> explicitUniformizedMatrix = uniformizedMatrix.toMatrix(odd, odd);
std::vector<ValueType> explicitB = b.template toVector<ValueType>(odd);
@ -241,8 +231,8 @@ namespace storm {
// Determine the set of states that achieved a strictly positive probability.
std::unique_ptr<CheckResult> statesWithValuesGreaterZero = hybridResult.compareAgainstBound(storm::logic::ComparisonType::Greater, storm::utility::zero<ValueType>());
// And use it to compute the set of relevant states.
storm::dd::Bdd<DdType> relevantStates = storm::utility::graph::performProbGreater0(this->getModel(), this->getModel().getTransitionMatrix().notZero(), phiStates, statesWithValuesGreaterZero->asSymbolicQualitativeCheckResult<DdType>().getTruthValuesVector());
// Compute the set of relevant states.
storm::dd::Bdd<DdType> relevantStates = statesWithProbabilityGreater0 && phiStates;
// Filter the unbounded result such that it only contains values for the relevant states.
hybridResult.filter(SymbolicQualitativeCheckResult<DdType>(this->getModel().getReachableStates(), relevantStates));
@ -252,28 +242,48 @@ namespace storm {
std::vector<ValueType> result;
std::unique_ptr<CheckResult> explicitResult = hybridResult.toExplicitQuantitativeCheckResult();
newSubresult = std::move(explicitResult->asExplicitQuantitativeCheckResult<ValueType>().getValueVector());
}
// Then compute the transient probabilities of being in such a state after t time units. For this,
// we must re-uniformize the CTMC, so we need to compute the second uniformized matrix.
uniformizationRate = 1.02 * (relevantStates.toAdd() * exitRates).getMax();
STORM_LOG_THROW(uniformizationRate > 0, storm::exceptions::InvalidStateException, "The uniformization rate must be positive.");
// If the lower and upper bounds coincide, we have only determined the relevant states at this
// point, but we still need to construct the starting vector.
if (lowerBound == upperBound) {
odd = storm::dd::Odd<DdType>(relevantStates);
newSubresult = psiStates.toAdd().template toVector<ValueType>(odd);
std::vector<ValueType> newSubresult = std::move(explicitResult->asExplicitQuantitativeCheckResult<ValueType>().getValueVector());
// Then compute the transient probabilities of being in such a state after t time units. For this,
// we must re-uniformize the CTMC, so we need to compute the second uniformized matrix.
uniformizationRate = 1.02 * (relevantStates.toAdd() * exitRates).getMax();
STORM_LOG_THROW(uniformizationRate > 0, storm::exceptions::InvalidStateException, "The uniformization rate must be positive.");
// If the lower and upper bounds coincide, we have only determined the relevant states at this
// point, but we still need to construct the starting vector.
if (lowerBound == upperBound) {
odd = storm::dd::Odd<DdType>(relevantStates);
newSubresult = psiStates.toAdd().template toVector<ValueType>(odd);
}
// Finally, we compute the second set of transient probabilities.
uniformizedMatrix = this->computeUniformizedMatrix(this->getModel(), this->getModel().getTransitionMatrix(), exitRates, relevantStates, uniformizationRate);
explicitUniformizedMatrix = uniformizedMatrix.toMatrix(odd, odd);
newSubresult = SparseCtmcCslModelChecker<ValueType>::computeTransientProbabilities(explicitUniformizedMatrix, nullptr, lowerBound, uniformizationRate, newSubresult, *this->linearEquationSolverFactory);
return std::unique_ptr<CheckResult>(new HybridQuantitativeCheckResult<DdType>(this->getModel().getReachableStates(), !relevantStates && this->getModel().getReachableStates(), this->getModel().getManager().getAddZero(), relevantStates, odd, newSubresult));
} else {
// In this case, the interval is of the form [t, t] with t != 0, t != inf.
// Build an ODD for the relevant states.
storm::dd::Odd<DdType> odd = storm::dd::Odd<DdType>(statesWithProbabilityGreater0);
std::vector<ValueType> newSubresult = psiStates.toAdd().template toVector<ValueType>(odd);
// Then compute the transient probabilities of being in such a state after t time units. For this,
// we must re-uniformize the CTMC, so we need to compute the second uniformized matrix.
ValueType uniformizationRate = 1.02 * (statesWithProbabilityGreater0.toAdd() * exitRates).getMax();
STORM_LOG_THROW(uniformizationRate > 0, storm::exceptions::InvalidStateException, "The uniformization rate must be positive.");
// Finally, we compute the second set of transient probabilities.
storm::dd::Add<DdType> uniformizedMatrix = this->computeUniformizedMatrix(this->getModel(), this->getModel().getTransitionMatrix(), exitRates, statesWithProbabilityGreater0, uniformizationRate);
storm::storage::SparseMatrix<ValueType> explicitUniformizedMatrix = uniformizedMatrix.toMatrix(odd, odd);
newSubresult = SparseCtmcCslModelChecker<ValueType>::computeTransientProbabilities(explicitUniformizedMatrix, nullptr, lowerBound, uniformizationRate, newSubresult, *this->linearEquationSolverFactory);
return std::unique_ptr<CheckResult>(new HybridQuantitativeCheckResult<DdType>(this->getModel().getReachableStates(), !statesWithProbabilityGreater0 && this->getModel().getReachableStates(), this->getModel().getManager().getAddZero(), statesWithProbabilityGreater0, odd, newSubresult));
}
// Finally, we compute the second set of transient probabilities.
uniformizedMatrix = this->computeUniformizedMatrix(this->getModel(), this->getModel().getTransitionMatrix(), exitRates, relevantStates, uniformizationRate);
storm::storage::SparseMatrix<ValueType> explicitUniformizedMatrix = uniformizedMatrix.toMatrix(odd, odd);
newSubresult = SparseCtmcCslModelChecker<ValueType>::computeTransientProbabilities(explicitUniformizedMatrix, nullptr, lowerBound, uniformizationRate, newSubresult, *this->linearEquationSolverFactory);
return std::unique_ptr<CheckResult>(new HybridQuantitativeCheckResult<DdType>(this->getModel().getReachableStates(), !relevantStates && this->getModel().getReachableStates(), this->getModel().getManager().getAddZero(), relevantStates, odd, newSubresult));
}
}
} else {
@ -281,76 +291,75 @@ namespace storm {
}
}
// template<storm::dd::DdType DdType, class ValueType>
// std::unique_ptr<CheckResult> HybridCtmcCslModelChecker<DdType, ValueType>::computeInstantaneousRewards(storm::logic::InstantaneousRewardFormula const& rewardPathFormula, bool qualitative, boost::optional<storm::logic::OptimalityType> const& optimalityType) {
// return std::unique_ptr<CheckResult>(new ExplicitQuantitativeCheckResult<ValueType>(this->computeInstantaneousRewardsHelper(rewardPathFormula.getContinuousTimeBound())));
// }
//
// template<storm::dd::DdType DdType, class ValueType>
// std::vector<ValueType> HybridCtmcCslModelChecker<DdType, ValueType>::computeInstantaneousRewardsHelper(double timeBound) const {
// // Only compute the result if the model has a state-based reward this->getModel().
// STORM_LOG_THROW(this->getModel().hasStateRewards(), storm::exceptions::InvalidPropertyException, "Missing reward model for formula. Skipping formula.");
//
// // Initialize result to state rewards of the this->getModel().
// std::vector<ValueType> result(this->getModel().getStateRewardVector());
//
// // If the time-bound is not zero, we need to perform a transient analysis.
// if (timeBound > 0) {
// ValueType uniformizationRate = 0;
// for (auto const& rate : this->getModel().getExitRateVector()) {
// uniformizationRate = std::max(uniformizationRate, rate);
// }
// uniformizationRate *= 1.02;
// STORM_LOG_THROW(uniformizationRate > 0, storm::exceptions::InvalidStateException, "The uniformization rate must be positive.");
//
// storm::storage::SparseMatrix<ValueType> uniformizedMatrix = this->computeUniformizedMatrix(this->getModel().getTransitionMatrix(), storm::storage::BitVector(this->getModel().getNumberOfStates(), true), uniformizationRate, this->getModel().getExitRateVector());
// result = this->computeTransientProbabilities(uniformizedMatrix, nullptr, timeBound, uniformizationRate, result, *this->linearEquationSolverFactory);
// }
//
// return result;
// }
//
// template<storm::dd::DdType DdType, class ValueType>
// std::unique_ptr<CheckResult> HybridCtmcCslModelChecker<DdType, ValueType>::computeCumulativeRewards(storm::logic::CumulativeRewardFormula const& rewardPathFormula, bool qualitative, boost::optional<storm::logic::OptimalityType> const& optimalityType) {
// return std::unique_ptr<CheckResult>(new ExplicitQuantitativeCheckResult<ValueType>(this->computeCumulativeRewardsHelper(rewardPathFormula.getContinuousTimeBound())));
// }
//
// template<storm::dd::DdType DdType, class ValueType>
// std::vector<ValueType> HybridCtmcCslModelChecker<DdType, ValueType>::computeCumulativeRewardsHelper(double timeBound) const {
// // Only compute the result if the model has a state-based reward this->getModel().
// STORM_LOG_THROW(this->getModel().hasStateRewards() || this->getModel().hasTransitionRewards(), storm::exceptions::InvalidPropertyException, "Missing reward model for formula. Skipping formula.");
//
// // If the time bound is zero, the result is the constant zero vector.
// if (timeBound == 0) {
// return std::vector<ValueType>(this->getModel().getNumberOfStates(), storm::utility::zero<ValueType>());
// }
//
// // Otherwise, we need to perform some computations.
//
// // Start with the uniformization.
// ValueType uniformizationRate = 0;
// for (auto const& rate : this->getModel().getExitRateVector()) {
// uniformizationRate = std::max(uniformizationRate, rate);
// }
// uniformizationRate *= 1.02;
// STORM_LOG_THROW(uniformizationRate > 0, storm::exceptions::InvalidStateException, "The uniformization rate must be positive.");
//
// storm::storage::SparseMatrix<ValueType> uniformizedMatrix = this->computeUniformizedMatrix(this->getModel().getTransitionMatrix(), storm::storage::BitVector(this->getModel().getNumberOfStates(), true), uniformizationRate, this->getModel().getExitRateVector());
//
// // Compute the total state reward vector.
// std::vector<ValueType> totalRewardVector;
// if (this->getModel().hasTransitionRewards()) {
// totalRewardVector = this->getModel().getTransitionMatrix().getPointwiseProductRowSumVector(this->getModel().getTransitionRewardMatrix());
// if (this->getModel().hasStateRewards()) {
// storm::utility::vector::addVectors(totalRewardVector, this->getModel().getStateRewardVector(), totalRewardVector);
// }
// } else {
// totalRewardVector = std::vector<ValueType>(this->getModel().getStateRewardVector());
// }
//
// // Finally, compute the transient probabilities.
// return this->computeTransientProbabilities<true>(uniformizedMatrix, nullptr, timeBound, uniformizationRate, totalRewardVector, *this->linearEquationSolverFactory);
// }
template<storm::dd::DdType DdType, class ValueType>
std::unique_ptr<CheckResult> HybridCtmcCslModelChecker<DdType, ValueType>::computeInstantaneousRewards(storm::logic::InstantaneousRewardFormula const& rewardPathFormula, bool qualitative, boost::optional<storm::logic::OptimalityType> const& optimalityType) {
return this->computeInstantaneousRewardsHelper(rewardPathFormula.getContinuousTimeBound());
}
template<storm::dd::DdType DdType, class ValueType>
std::unique_ptr<CheckResult> HybridCtmcCslModelChecker<DdType, ValueType>::computeInstantaneousRewardsHelper(double timeBound) const {
// Only compute the result if the model has a state-based reward this->getModel().
STORM_LOG_THROW(this->getModel().hasStateRewards(), storm::exceptions::InvalidPropertyException, "Missing reward model for formula. Skipping formula.");
// Create ODD for the translation.
storm::dd::Odd<DdType> odd(this->getModel().getReachableStates());
// Initialize result to state rewards of the this->getModel().
std::vector<ValueType> result = this->getModel().getStateRewardVector().template toVector<ValueType>(odd);
// If the time-bound is not zero, we need to perform a transient analysis.
if (timeBound > 0) {
ValueType uniformizationRate = 1.02 * this->getModel().getExitRateVector().getMax();
STORM_LOG_THROW(uniformizationRate > 0, storm::exceptions::InvalidStateException, "The uniformization rate must be positive.");
storm::dd::Add<DdType> uniformizedMatrix = this->computeUniformizedMatrix(this->getModel(), this->getModel().getTransitionMatrix(), this->getModel().getExitRateVector(), this->getModel().getReachableStates(), uniformizationRate);
storm::storage::SparseMatrix<ValueType> explicitUniformizedMatrix = uniformizedMatrix.toMatrix(odd, odd);
result = SparseCtmcCslModelChecker<ValueType>::computeTransientProbabilities(explicitUniformizedMatrix, nullptr, timeBound, uniformizationRate, result, *this->linearEquationSolverFactory);
}
return std::unique_ptr<CheckResult>(new HybridQuantitativeCheckResult<DdType>(this->getModel().getReachableStates(), this->getModel().getManager().getBddZero(), this->getModel().getManager().getAddZero(), this->getModel().getReachableStates(), odd, result));
}
template<storm::dd::DdType DdType, class ValueType>
std::unique_ptr<CheckResult> HybridCtmcCslModelChecker<DdType, ValueType>::computeCumulativeRewards(storm::logic::CumulativeRewardFormula const& rewardPathFormula, bool qualitative, boost::optional<storm::logic::OptimalityType> const& optimalityType) {
return this->computeCumulativeRewardsHelper(rewardPathFormula.getContinuousTimeBound());
}
template<storm::dd::DdType DdType, class ValueType>
std::unique_ptr<CheckResult> HybridCtmcCslModelChecker<DdType, ValueType>::computeCumulativeRewardsHelper(double timeBound) const {
// Only compute the result if the model has a state-based reward this->getModel().
STORM_LOG_THROW(this->getModel().hasStateRewards() || this->getModel().hasTransitionRewards(), storm::exceptions::InvalidPropertyException, "Missing reward model for formula. Skipping formula.");
// If the time bound is zero, the result is the constant zero vector.
if (timeBound == 0) {
return std::unique_ptr<CheckResult>(new SymbolicQuantitativeCheckResult<DdType>(this->getModel().getReachableStates(), this->getModel().getManager().getAddZero()));
}
// Otherwise, we need to perform some computations.
// Start with the uniformization.
ValueType uniformizationRate = 1.02 * this->getModel().getExitRateVector().getMax();
STORM_LOG_THROW(uniformizationRate > 0, storm::exceptions::InvalidStateException, "The uniformization rate must be positive.");
// Create ODD for the translation.
storm::dd::Odd<DdType> odd(this->getModel().getReachableStates());
// Compute the uniformized matrix.
storm::dd::Add<DdType> uniformizedMatrix = this->computeUniformizedMatrix(this->getModel(), this->getModel().getTransitionMatrix(), this->getModel().getExitRateVector(), this->getModel().getReachableStates(), uniformizationRate);
storm::storage::SparseMatrix<ValueType> explicitUniformizedMatrix = uniformizedMatrix.toMatrix(odd, odd);
// Then compute the state reward vector to use in the computation.
storm::dd::Add<DdType> totalRewardVector = this->getModel().hasStateRewards() ? this->getModel().getStateRewardVector() : this->getModel().getManager().getAddZero();
if (this->getModel().hasTransitionRewards()) {
totalRewardVector += (this->getModel().getTransitionMatrix() * this->getModel().getTransitionRewardMatrix()).sumAbstract(this->getModel().getColumnVariables());
}
std::vector<ValueType> explicitTotalRewardVector = totalRewardVector.template toVector<ValueType>(odd);
// Finally, compute the transient probabilities.
std::vector<ValueType> result = SparseCtmcCslModelChecker<ValueType>::template computeTransientProbabilities<true>(explicitUniformizedMatrix, nullptr, timeBound, uniformizationRate, explicitTotalRewardVector, *this->linearEquationSolverFactory);
return std::unique_ptr<CheckResult>(new HybridQuantitativeCheckResult<DdType>(this->getModel().getReachableStates(), this->getModel().getManager().getBddZero(), this->getModel().getManager().getAddZero(), this->getModel().getReachableStates(), odd, result));
}
// Explicitly instantiate the model checker.
template class HybridCtmcCslModelChecker<storm::dd::DdType::CUDD, double>;

39
src/modelchecker/csl/HybridCtmcCslModelChecker.h

@ -20,8 +20,8 @@ namespace storm {
virtual std::unique_ptr<CheckResult> computeBoundedUntilProbabilities(storm::logic::BoundedUntilFormula const& pathFormula, bool qualitative = false, boost::optional<storm::logic::OptimalityType> const& optimalityType = boost::optional<storm::logic::OptimalityType>()) override;
virtual std::unique_ptr<CheckResult> computeNextProbabilities(storm::logic::NextFormula const& pathFormula, bool qualitative = false, boost::optional<storm::logic::OptimalityType> const& optimalityType = boost::optional<storm::logic::OptimalityType>()) override;
virtual std::unique_ptr<CheckResult> computeUntilProbabilities(storm::logic::UntilFormula const& pathFormula, bool qualitative = false, boost::optional<storm::logic::OptimalityType> const& optimalityType = boost::optional<storm::logic::OptimalityType>()) override;
// virtual std::unique_ptr<CheckResult> computeInstantaneousRewards(storm::logic::InstantaneousRewardFormula const& rewardPathFormula, bool qualitative = false, boost::optional<storm::logic::OptimalityType> const& optimalityType = boost::optional<storm::logic::OptimalityType>()) override;
// virtual std::unique_ptr<CheckResult> computeCumulativeRewards(storm::logic::CumulativeRewardFormula const& rewardPathFormula, bool qualitative = false, boost::optional<storm::logic::OptimalityType> const& optimalityType = boost::optional<storm::logic::OptimalityType>()) override;
virtual std::unique_ptr<CheckResult> computeInstantaneousRewards(storm::logic::InstantaneousRewardFormula const& rewardPathFormula, bool qualitative = false, boost::optional<storm::logic::OptimalityType> const& optimalityType = boost::optional<storm::logic::OptimalityType>()) override;
virtual std::unique_ptr<CheckResult> computeCumulativeRewards(storm::logic::CumulativeRewardFormula const& rewardPathFormula, bool qualitative = false, boost::optional<storm::logic::OptimalityType> const& optimalityType = boost::optional<storm::logic::OptimalityType>()) override;
virtual std::unique_ptr<CheckResult> computeReachabilityRewards(storm::logic::ReachabilityRewardFormula const& rewardPathFormula, bool qualitative = false, boost::optional<storm::logic::OptimalityType> const& optimalityType = boost::optional<storm::logic::OptimalityType>()) override;
protected:
@ -52,39 +52,8 @@ namespace storm {
// The methods that perform the actual checking.
std::unique_ptr<CheckResult> computeBoundedUntilProbabilitiesHelper(storm::dd::Bdd<DdType> const& phiStates, storm::dd::Bdd<DdType> const& psiStates, storm::dd::Add<DdType> const& exitRates, bool qualitative, double lowerBound, double upperBound) const;
// std::vector<ValueType> computeBoundedUntilProbabilitiesHelper(storm::storage::BitVector const& phiStates, storm::storage::BitVector const& psiStates, std::vector<ValueType> const& exitRates, bool qualitative, double lowerBound, double upperBound) const;
// std::vector<ValueType> computeInstantaneousRewardsHelper(double timeBound) const;
// std::vector<ValueType> computeCumulativeRewardsHelper(double timeBound) const;
// std::vector<ValueType> computeReachabilityRewardsHelper(storm::storage::BitVector const& targetStates, bool qualitative) const;
//
// /*!
// * Computes the matrix representing the transitions of the uniformized CTMC.
// *
// * @param transitionMatrix The matrix to uniformize.
// * @param maybeStates The states that need to be considered.
// * @param uniformizationRate The rate to be used for uniformization.
// * @param exitRates The exit rates of all states.
// * @return The uniformized matrix.
// */
// static storm::storage::SparseMatrix<ValueType> computeUniformizedMatrix(storm::storage::SparseMatrix<ValueType> const& transitionMatrix, storm::storage::BitVector const& maybeStates, ValueType uniformizationRate, std::vector<ValueType> const& exitRates);
//
// /*!
// * Computes the transient probabilities for lambda time steps.
// *
// * @param uniformizedMatrix The uniformized transition matrix.
// * @param addVector A vector that is added in each step as a possible compensation for removing absorbing states
// * with a non-zero initial value. If this is not supposed to be used, it can be set to nullptr.
// * @param timeBound The time bound to use.
// * @param uniformizationRate The used uniformization rate.
// * @param values A vector mapping each state to an initial probability.
// * @param linearEquationSolverFactory The factory to use when instantiating new linear equation solvers.
// * @param useMixedPoissonProbabilities If set to true, instead of taking the poisson probabilities, mixed
// * poisson probabilities are used.
// * @return The vector of transient probabilities.
// */
// template<bool useMixedPoissonProbabilities = false>
// std::vector<ValueType> computeTransientProbabilities(storm::storage::SparseMatrix<ValueType> const& uniformizedMatrix, std::vector<ValueType> const* addVector, ValueType timeBound, ValueType uniformizationRate, std::vector<ValueType> values, storm::utility::solver::LinearEquationSolverFactory<ValueType> const& linearEquationSolverFactory) const;
std::unique_ptr<CheckResult> computeInstantaneousRewardsHelper(double timeBound) const;
std::unique_ptr<CheckResult> computeCumulativeRewardsHelper(double timeBound) const;
// An object that is used for solving linear equations and performing matrix-vector multiplication.
std::unique_ptr<storm::utility::solver::LinearEquationSolverFactory<ValueType>> linearEquationSolverFactory;

24
src/modelchecker/csl/SparseCtmcCslModelChecker.cpp

@ -161,14 +161,11 @@ namespace storm {
} else {
// In this case, the interval is of the form [t, t'] with t != 0 and t' != inf.
// Prepare some variables that are used by the two following blocks.
ValueType uniformizationRate = 0;
storm::storage::SparseMatrix<ValueType> uniformizedMatrix;
std::vector<ValueType> newSubresult;
if (lowerBound != upperBound) {
// In this case, the interval is of the form [t, t'] with t != 0, t' != inf and t != t'.
// Find the maximal rate of all 'maybe' states to take it as the uniformization rate.
uniformizationRate = 0;
ValueType uniformizationRate = storm::utility::zero<ValueType>();
for (auto const& state : statesWithProbabilityGreater0NonPsi) {
uniformizationRate = std::max(uniformizationRate, exitRates[state]);
}
@ -176,7 +173,7 @@ namespace storm {
STORM_LOG_THROW(uniformizationRate > 0, storm::exceptions::InvalidStateException, "The uniformization rate must be positive.");
// Compute the (first) uniformized matrix.
uniformizedMatrix = this->computeUniformizedMatrix(this->getModel().getTransitionMatrix(), statesWithProbabilityGreater0NonPsi, uniformizationRate, exitRates);
storm::storage::SparseMatrix<ValueType> uniformizedMatrix = this->computeUniformizedMatrix(this->getModel().getTransitionMatrix(), statesWithProbabilityGreater0NonPsi, uniformizationRate, exitRates);
// Compute the vector that is to be added as a compensation for removing the absorbing states.
std::vector<ValueType> b = this->getModel().getTransitionMatrix().getConstrainedRowSumVector(statesWithProbabilityGreater0NonPsi, psiStates);
@ -189,13 +186,13 @@ namespace storm {
std::vector<ValueType> subresult = computeTransientProbabilities(uniformizedMatrix, &b, upperBound - lowerBound, uniformizationRate, values, *this->linearEquationSolverFactory);
storm::storage::BitVector relevantStates = statesWithProbabilityGreater0 & phiStates;
newSubresult = std::vector<ValueType>(relevantStates.getNumberOfSetBits());
std::vector<ValueType> newSubresult = std::vector<ValueType>(relevantStates.getNumberOfSetBits());
storm::utility::vector::selectVectorValues(newSubresult, statesWithProbabilityGreater0NonPsi % relevantStates, subresult);
storm::utility::vector::setVectorValues(newSubresult, psiStates % relevantStates, storm::utility::one<ValueType>());
// Then compute the transient probabilities of being in such a state after t time units. For this,
// we must re-uniformize the CTMC, so we need to compute the second uniformized matrix.
uniformizationRate = 0;
uniformizationRate = storm::utility::zero<ValueType>();
for (auto const& state : relevantStates) {
uniformizationRate = std::max(uniformizationRate, exitRates[state]);
}
@ -211,12 +208,14 @@ namespace storm {
storm::utility::vector::setVectorValues(result, ~relevantStates, storm::utility::zero<ValueType>());
storm::utility::vector::setVectorValues(result, relevantStates, newSubresult);
} else {
newSubresult = std::vector<ValueType>(statesWithProbabilityGreater0.getNumberOfSetBits());
// In this case, the interval is of the form [t, t] with t != 0, t != inf.
std::vector<ValueType> newSubresult = std::vector<ValueType>(statesWithProbabilityGreater0.getNumberOfSetBits());
storm::utility::vector::setVectorValues(newSubresult, psiStates % statesWithProbabilityGreater0, storm::utility::one<ValueType>());
// Then compute the transient probabilities of being in such a state after t time units. For this,
// we must re-uniformize the CTMC, so we need to compute the second uniformized matrix.
uniformizationRate = 0;
ValueType uniformizationRate = storm::utility::zero<ValueType>();
for (auto const& state : statesWithProbabilityGreater0) {
uniformizationRate = std::max(uniformizationRate, exitRates[state]);
}
@ -224,14 +223,13 @@ namespace storm {
STORM_LOG_THROW(uniformizationRate > 0, storm::exceptions::InvalidStateException, "The uniformization rate must be positive.");
// Finally, we compute the second set of transient probabilities.
uniformizedMatrix = this->computeUniformizedMatrix(this->getModel().getTransitionMatrix(), statesWithProbabilityGreater0, uniformizationRate, exitRates);
storm::storage::SparseMatrix<ValueType> uniformizedMatrix = this->computeUniformizedMatrix(this->getModel().getTransitionMatrix(), statesWithProbabilityGreater0, uniformizationRate, exitRates);
newSubresult = computeTransientProbabilities(uniformizedMatrix, nullptr, lowerBound, uniformizationRate, newSubresult, *this->linearEquationSolverFactory);
// Fill in the correct values.
result = std::vector<ValueType>(this->getModel().getNumberOfStates(), storm::utility::zero<ValueType>());
storm::utility::vector::setVectorValues(result, ~statesWithProbabilityGreater0, storm::utility::zero<ValueType>());
storm::utility::vector::setVectorValues(result, statesWithProbabilityGreater0, newSubresult);
}
}
}

6
src/settings/Argument.h

@ -73,12 +73,12 @@ namespace storm {
return this->setFromTypeValue(newValue);
}
bool setFromTypeValue(T const& newValue) {
bool setFromTypeValue(T const& newValue, bool hasBeenSet = true) {
if (!this->validate(newValue)) {
return false;
}
this->argumentValue = newValue;
this->hasBeenSet = true;
this->hasBeenSet = hasBeenSet;
return true;
}
@ -120,7 +120,7 @@ namespace storm {
void setFromDefaultValue() override {
STORM_LOG_THROW(this->hasDefaultValue, storm::exceptions::IllegalFunctionCallException, "Unable to set value from default value, because the argument has none.");
bool result = this->setFromTypeValue(this->defaultValue);
bool result = this->setFromTypeValue(this->defaultValue, false);
STORM_LOG_THROW(result, storm::exceptions::IllegalArgumentValueException, "Unable to assign default value to argument, because it was rejected.");
}

15
src/storage/prism/Program.cpp

@ -595,6 +595,7 @@ namespace storm {
// Check the commands of the modules.
bool hasProbabilisticCommand = false;
bool hasMarkovianCommand = false;
bool hasLabeledMarkovianCommand = false;
for (auto const& module : this->getModules()) {
std::set<storm::expressions::Variable> legalVariables = globalVariables;
for (auto const& variable : module.getBooleanVariables()) {
@ -621,11 +622,7 @@ namespace storm {
// If the command is Markovian and labeled, we throw an error or raise a warning, depending on
// whether or not the PRISM compatibility mode was enabled.
if (command.isMarkovian() && command.isLabeled()) {
if (storm::settings::generalSettings().isPrismCompatibilityEnabled()) {
STORM_LOG_WARN_COND(false, "The model uses synchronizing Markovian commands. This may lead to unexpected verification results, because of unclear semantics.");
} else {
STORM_LOG_THROW(false, storm::exceptions::WrongFormatException, "The model uses synchronizing Markovian commands. This may lead to unexpected verification results, because of unclear semantics.");
}
hasLabeledMarkovianCommand = true;
}
// Check all updates.
@ -665,6 +662,14 @@ namespace storm {
}
}
if (hasLabeledMarkovianCommand) {
if (storm::settings::generalSettings().isPrismCompatibilityEnabled()) {
STORM_LOG_WARN_COND(false, "The model uses synchronizing Markovian commands. This may lead to unexpected verification results, because of unclear semantics.");
} else {
STORM_LOG_THROW(false, storm::exceptions::WrongFormatException, "The model uses synchronizing Markovian commands. This may lead to unexpected verification results, because of unclear semantics.");
}
}
if (this->getModelType() == Program::ModelType::DTMC || this->getModelType() == Program::ModelType::MDP) {
STORM_LOG_THROW(!hasMarkovianCommand, storm::exceptions::WrongFormatException, "Discrete-time model must not have Markovian commands.");
} else if (this->getModelType() == Program::ModelType::CTMC) {

12
src/utility/cli.h

@ -318,8 +318,12 @@ namespace storm {
std::string constants = settings.getConstantDefinitionString();
bool buildRewards = false;
if (formula) {
boost::optional<std::string> rewardModelName;
if (formula || settings.isSymbolicRewardModelNameSet()) {
buildRewards = formula.get()->isRewardOperatorFormula() || formula.get()->isRewardPathFormula();
if (settings.isSymbolicRewardModelNameSet()) {
rewardModelName = settings.getSymbolicRewardModelName();
}
}
// Customize and perform model-building.
@ -329,6 +333,8 @@ namespace storm {
options = typename storm::builder::ExplicitPrismModelBuilder<ValueType>::Options(*formula.get());
}
options.addConstantDefinitionsFromString(program, settings.getConstantDefinitionString());
options.buildRewards = buildRewards;
options.rewardModelName = rewardModelName;
// Generate command labels if we are going to build a counterexample later.
if (storm::settings::counterexampleGeneratorSettings().isMinimalCommandSetGenerationSet()) {
@ -342,7 +348,9 @@ namespace storm {
options = typename storm::builder::DdPrismModelBuilder<storm::dd::DdType::CUDD>::Options(*formula.get());
}
options.addConstantDefinitionsFromString(program, settings.getConstantDefinitionString());
options.buildRewards = buildRewards;
options.rewardModelName = rewardModelName;
result = storm::builder::DdPrismModelBuilder<storm::dd::DdType::CUDD>::translateProgram(program, options);
}

11
test/functional/modelchecker/GmmxxCtmcCslModelCheckerTest.cpp

@ -69,12 +69,19 @@ TEST(GmmxxCtmcCslModelCheckerTest, Cluster) {
storm::modelchecker::ExplicitQuantitativeCheckResult<double> quantitativeCheckResult5 = checkResult->asExplicitQuantitativeCheckResult<double>();
EXPECT_NEAR(0, quantitativeCheckResult5[initialState], storm::settings::generalSettings().getPrecision());
formula = formulaParser.parseFromString("R=? [C<=100]");
formula = formulaParser.parseFromString("P=? [ \"minimum\" U[1,inf] !\"minimum\"]");
checkResult = modelchecker.check(*formula);
ASSERT_TRUE(checkResult->isExplicitQuantitativeCheckResult());
storm::modelchecker::ExplicitQuantitativeCheckResult<double> quantitativeCheckResult6 = checkResult->asExplicitQuantitativeCheckResult<double>();
EXPECT_NEAR(0.8602815057967503, quantitativeCheckResult6[initialState], storm::settings::generalSettings().getPrecision());
EXPECT_NEAR(0.9999999033633374, quantitativeCheckResult6[initialState], storm::settings::generalSettings().getPrecision());
formula = formulaParser.parseFromString("R=? [C<=100]");
checkResult = modelchecker.check(*formula);
ASSERT_TRUE(checkResult->isExplicitQuantitativeCheckResult());
storm::modelchecker::ExplicitQuantitativeCheckResult<double> quantitativeCheckResult7 = checkResult->asExplicitQuantitativeCheckResult<double>();
EXPECT_NEAR(0.8602815057967503, quantitativeCheckResult7[initialState], storm::settings::generalSettings().getPrecision());
}
TEST(GmmxxCtmcCslModelCheckerTest, Embedded) {

100
test/functional/modelchecker/GmmxxHybridCtmcCslModelCheckerTest.cpp

@ -50,36 +50,54 @@ TEST(GmmxxHybridCtmcCslModelCheckerTest, Cluster) {
ASSERT_TRUE(checkResult->isHybridQuantitativeCheckResult());
checkResult->filter(storm::modelchecker::SymbolicQualitativeCheckResult<storm::dd::DdType::CUDD>(ctmc->getReachableStates(), ctmc->getInitialStates()));
storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::CUDD> quantitativeCheckResult2 = checkResult->asHybridQuantitativeCheckResult<storm::dd::DdType::CUDD>();
storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::CUDD> quantitativeCheckResult2 = checkResult->asHybridQuantitativeCheckResult<storm::dd::DdType::CUDD>();
EXPECT_NEAR(2.3397873548343415E-6, quantitativeCheckResult2.getMin(), storm::settings::generalSettings().getPrecision());
EXPECT_NEAR(2.3397873548343415E-6, quantitativeCheckResult2.getMax(), storm::settings::generalSettings().getPrecision());
formula = formulaParser.parseFromString("P=? [ F[100,2000] !\"minimum\"]");
checkResult = modelchecker.check(*formula);
ASSERT_TRUE(checkResult->isHybridQuantitativeCheckResult());
checkResult->filter(storm::modelchecker::SymbolicQualitativeCheckResult<storm::dd::DdType::CUDD>(ctmc->getReachableStates(), ctmc->getInitialStates()));
storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::CUDD> quantitativeCheckResult3 = checkResult->asHybridQuantitativeCheckResult<storm::dd::DdType::CUDD>();
EXPECT_NEAR(0.001112543139248814, quantitativeCheckResult3.getMin(), storm::settings::generalSettings().getPrecision());
EXPECT_NEAR(0.001112543139248814, quantitativeCheckResult3.getMax(), storm::settings::generalSettings().getPrecision());
formula = formulaParser.parseFromString("P=? [ \"minimum\" U<=10 \"premium\"]");
checkResult = modelchecker.check(*formula);
ASSERT_TRUE(checkResult->isHybridQuantitativeCheckResult());
checkResult->filter(storm::modelchecker::SymbolicQualitativeCheckResult<storm::dd::DdType::CUDD>(ctmc->getReachableStates(), ctmc->getInitialStates()));
storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::CUDD> quantitativeCheckResult3 = checkResult->asHybridQuantitativeCheckResult<storm::dd::DdType::CUDD>();
EXPECT_NEAR(1, quantitativeCheckResult3.getMin(), storm::settings::generalSettings().getPrecision());
EXPECT_NEAR(1, quantitativeCheckResult3.getMax(), storm::settings::generalSettings().getPrecision());
storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::CUDD> quantitativeCheckResult4 = checkResult->asHybridQuantitativeCheckResult<storm::dd::DdType::CUDD>();
EXPECT_NEAR(1, quantitativeCheckResult4.getMin(), storm::settings::generalSettings().getPrecision());
EXPECT_NEAR(1, quantitativeCheckResult4.getMax(), storm::settings::generalSettings().getPrecision());
formula = formulaParser.parseFromString("P=? [ !\"minimum\" U[1,inf] \"minimum\"]");
checkResult = modelchecker.check(*formula);
ASSERT_TRUE(checkResult->isHybridQuantitativeCheckResult());
checkResult->filter(storm::modelchecker::SymbolicQualitativeCheckResult<storm::dd::DdType::CUDD>(ctmc->getReachableStates(), ctmc->getInitialStates()));
storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::CUDD> quantitativeCheckResult4 = checkResult->asHybridQuantitativeCheckResult<storm::dd::DdType::CUDD>();
EXPECT_NEAR(0, quantitativeCheckResult4.getMin(), storm::settings::generalSettings().getPrecision());
EXPECT_NEAR(0, quantitativeCheckResult4.getMax(), storm::settings::generalSettings().getPrecision());
// formula = formulaParser.parseFromString("R=? [C<=100]");
// checkResult = modelchecker.check(*formula);
//
// ASSERT_TRUE(checkResult->isHybridQuantitativeCheckResult());
// checkResult->filter(storm::modelchecker::SymbolicQualitativeCheckResult<storm::dd::DdType::CUDD>(ctmc->getReachableStates(), ctmc->getInitialStates()));
// storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::CUDD> quantitativeCheckResult5 = checkResult->asHybridQuantitativeCheckResult<storm::dd::DdType::CUDD>();
// EXPECT_NEAR(0.8602815057967503, quantitativeCheckResult5.getMin(), storm::settings::generalSettings().getPrecision());
// EXPECT_NEAR(0.8602815057967503, quantitativeCheckResult5.getMax(), storm::settings::generalSettings().getPrecision());
storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::CUDD> quantitativeCheckResult5 = checkResult->asHybridQuantitativeCheckResult<storm::dd::DdType::CUDD>();
EXPECT_NEAR(0, quantitativeCheckResult5.getMin(), storm::settings::generalSettings().getPrecision());
EXPECT_NEAR(0, quantitativeCheckResult5.getMax(), storm::settings::generalSettings().getPrecision());
formula = formulaParser.parseFromString("P=? [ \"minimum\" U[1,inf] !\"minimum\"]");
checkResult = modelchecker.check(*formula);
ASSERT_TRUE(checkResult->isHybridQuantitativeCheckResult());
checkResult->filter(storm::modelchecker::SymbolicQualitativeCheckResult<storm::dd::DdType::CUDD>(ctmc->getReachableStates(), ctmc->getInitialStates()));
storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::CUDD> quantitativeCheckResult6 = checkResult->asHybridQuantitativeCheckResult<storm::dd::DdType::CUDD>();
EXPECT_NEAR(0.9999999033633374, quantitativeCheckResult6.getMin(), storm::settings::generalSettings().getPrecision());
EXPECT_NEAR(0.9999999033633374, quantitativeCheckResult6.getMax(), storm::settings::generalSettings().getPrecision());
formula = formulaParser.parseFromString("R=? [C<=100]");
checkResult = modelchecker.check(*formula);
ASSERT_TRUE(checkResult->isHybridQuantitativeCheckResult());
checkResult->filter(storm::modelchecker::SymbolicQualitativeCheckResult<storm::dd::DdType::CUDD>(ctmc->getReachableStates(), ctmc->getInitialStates()));
storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::CUDD> quantitativeCheckResult7 = checkResult->asHybridQuantitativeCheckResult<storm::dd::DdType::CUDD>();
EXPECT_NEAR(0.8602815057967503, quantitativeCheckResult7.getMin(), storm::settings::generalSettings().getPrecision());
EXPECT_NEAR(0.8602815057967503, quantitativeCheckResult7.getMax(), storm::settings::generalSettings().getPrecision());
}
TEST(GmmxxHybridCtmcCslModelCheckerTest, Embedded) {
@ -139,14 +157,14 @@ TEST(GmmxxHybridCtmcCslModelCheckerTest, Embedded) {
EXPECT_NEAR(4.429620626755424E-5, quantitativeCheckResult4.getMin(), storm::settings::generalSettings().getPrecision());
EXPECT_NEAR(4.429620626755424E-5, quantitativeCheckResult4.getMax(), storm::settings::generalSettings().getPrecision());
// formula = formulaParser.parseFromString("R=? [C<=10000]");
// checkResult = modelchecker.check(*formula);
//
// ASSERT_TRUE(checkResult->isHybridQuantitativeCheckResult());
// checkResult->filter(storm::modelchecker::SymbolicQualitativeCheckResult<storm::dd::DdType::CUDD>(ctmc->getReachableStates(), ctmc->getInitialStates()));
// storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::CUDD> quantitativeCheckResult5 = checkResult->asHybridQuantitativeCheckResult<storm::dd::DdType::CUDD>();
// EXPECT_NEAR(2.7720429852369946, quantitativeCheckResult5.getMin(), storm::settings::generalSettings().getPrecision());
// EXPECT_NEAR(2.7720429852369946, quantitativeCheckResult5.getMax(), storm::settings::generalSettings().getPrecision());
formula = formulaParser.parseFromString("R=? [C<=10000]");
checkResult = modelchecker.check(*formula);
ASSERT_TRUE(checkResult->isHybridQuantitativeCheckResult());
checkResult->filter(storm::modelchecker::SymbolicQualitativeCheckResult<storm::dd::DdType::CUDD>(ctmc->getReachableStates(), ctmc->getInitialStates()));
storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::CUDD> quantitativeCheckResult5 = checkResult->asHybridQuantitativeCheckResult<storm::dd::DdType::CUDD>();
EXPECT_NEAR(2.7720429852369946, quantitativeCheckResult5.getMin(), storm::settings::generalSettings().getPrecision());
EXPECT_NEAR(2.7720429852369946, quantitativeCheckResult5.getMax(), storm::settings::generalSettings().getPrecision());
}
TEST(GmmxxHybridCtmcCslModelCheckerTest, Polling) {
@ -232,23 +250,23 @@ TEST(GmmxxHybridCtmcCslModelCheckerTest, Tandem) {
EXPECT_NEAR(1, quantitativeCheckResult3.getMin(), storm::settings::generalSettings().getPrecision());
EXPECT_NEAR(1, quantitativeCheckResult3.getMax(), storm::settings::generalSettings().getPrecision());
// formula = formulaParser.parseFromString("R=? [I=10]");
// checkResult = modelchecker.check(*formula);
//
// ASSERT_TRUE(checkResult->isHybridQuantitativeCheckResult());
// checkResult->filter(storm::modelchecker::SymbolicQualitativeCheckResult<storm::dd::DdType::CUDD>(ctmc->getReachableStates(), ctmc->getInitialStates()));
// storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::CUDD> quantitativeCheckResult4 = checkResult->asHybridQuantitativeCheckResult<storm::dd::DdType::CUDD>();
// EXPECT_NEAR(5.679243850315877, quantitativeCheckResult4.getMin(), storm::settings::generalSettings().getPrecision());
// EXPECT_NEAR(5.679243850315877, quantitativeCheckResult4.getMax(), storm::settings::generalSettings().getPrecision());
//
// formula = formulaParser.parseFromString("R=? [C<=10]");
// checkResult = modelchecker.check(*formula);
//
// ASSERT_TRUE(checkResult->isHybridQuantitativeCheckResult());
// checkResult->filter(storm::modelchecker::SymbolicQualitativeCheckResult<storm::dd::DdType::CUDD>(ctmc->getReachableStates(), ctmc->getInitialStates()));
// storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::CUDD> quantitativeCheckResult5 = checkResult->asHybridQuantitativeCheckResult<storm::dd::DdType::CUDD>();
// EXPECT_NEAR(55.44792186036232, quantitativeCheckResult5.getMin(), storm::settings::generalSettings().getPrecision());
// EXPECT_NEAR(55.44792186036232, quantitativeCheckResult5.getMax(), storm::settings::generalSettings().getPrecision());
formula = formulaParser.parseFromString("R=? [I=10]");
checkResult = modelchecker.check(*formula);
ASSERT_TRUE(checkResult->isHybridQuantitativeCheckResult());
checkResult->filter(storm::modelchecker::SymbolicQualitativeCheckResult<storm::dd::DdType::CUDD>(ctmc->getReachableStates(), ctmc->getInitialStates()));
storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::CUDD> quantitativeCheckResult4 = checkResult->asHybridQuantitativeCheckResult<storm::dd::DdType::CUDD>();
EXPECT_NEAR(5.679243850315877, quantitativeCheckResult4.getMin(), storm::settings::generalSettings().getPrecision());
EXPECT_NEAR(5.679243850315877, quantitativeCheckResult4.getMax(), storm::settings::generalSettings().getPrecision());
formula = formulaParser.parseFromString("R=? [C<=10]");
checkResult = modelchecker.check(*formula);
ASSERT_TRUE(checkResult->isHybridQuantitativeCheckResult());
checkResult->filter(storm::modelchecker::SymbolicQualitativeCheckResult<storm::dd::DdType::CUDD>(ctmc->getReachableStates(), ctmc->getInitialStates()));
storm::modelchecker::HybridQuantitativeCheckResult<storm::dd::DdType::CUDD> quantitativeCheckResult5 = checkResult->asHybridQuantitativeCheckResult<storm::dd::DdType::CUDD>();
EXPECT_NEAR(55.44792186036232, quantitativeCheckResult5.getMin(), storm::settings::generalSettings().getPrecision());
EXPECT_NEAR(55.44792186036232, quantitativeCheckResult5.getMax(), storm::settings::generalSettings().getPrecision());
formula = formulaParser.parseFromString("R=? [F \"first_queue_full\"&\"second_queue_full\"]");
checkResult = modelchecker.check(*formula);

30
test/functional/modelchecker/SparseCtmcCslModelCheckerTest.cpp

@ -47,27 +47,41 @@ TEST(SparseCtmcCslModelCheckerTest, Cluster) {
ASSERT_TRUE(checkResult->isExplicitQuantitativeCheckResult());
storm::modelchecker::ExplicitQuantitativeCheckResult<double> quantitativeCheckResult2 = checkResult->asExplicitQuantitativeCheckResult<double>();
EXPECT_NEAR(2.3397873548343415E-6, quantitativeCheckResult2[initialState], storm::settings::generalSettings().getPrecision());
formula = formulaParser.parseFromString("P=? [ \"minimum\" U<=10 \"premium\"]");
formula = formulaParser.parseFromString("P=? [ F[100,2000] !\"minimum\"]");
checkResult = modelchecker.check(*formula);
ASSERT_TRUE(checkResult->isExplicitQuantitativeCheckResult());
storm::modelchecker::ExplicitQuantitativeCheckResult<double> quantitativeCheckResult3 = checkResult->asExplicitQuantitativeCheckResult<double>();
EXPECT_NEAR(1, quantitativeCheckResult3[initialState], storm::settings::generalSettings().getPrecision());
formula = formulaParser.parseFromString("P=? [ !\"minimum\" U[1,inf] \"minimum\"]");
EXPECT_NEAR(0.001105335651650576, quantitativeCheckResult3[initialState], storm::settings::generalSettings().getPrecision());
formula = formulaParser.parseFromString("P=? [ \"minimum\" U<=10 \"premium\"]");
checkResult = modelchecker.check(*formula);
ASSERT_TRUE(checkResult->isExplicitQuantitativeCheckResult());
storm::modelchecker::ExplicitQuantitativeCheckResult<double> quantitativeCheckResult4 = checkResult->asExplicitQuantitativeCheckResult<double>();
EXPECT_NEAR(0, quantitativeCheckResult4[initialState], storm::settings::generalSettings().getPrecision());
EXPECT_NEAR(1, quantitativeCheckResult4[initialState], storm::settings::generalSettings().getPrecision());
formula = formulaParser.parseFromString("R=? [C<=100]");
formula = formulaParser.parseFromString("P=? [ !\"minimum\" U[1,inf] \"minimum\"]");
checkResult = modelchecker.check(*formula);
ASSERT_TRUE(checkResult->isExplicitQuantitativeCheckResult());
storm::modelchecker::ExplicitQuantitativeCheckResult<double> quantitativeCheckResult5 = checkResult->asExplicitQuantitativeCheckResult<double>();
EXPECT_NEAR(0.8602815057967503, quantitativeCheckResult5[initialState], storm::settings::generalSettings().getPrecision());
EXPECT_NEAR(0, quantitativeCheckResult5[initialState], storm::settings::generalSettings().getPrecision());
formula = formulaParser.parseFromString("P=? [ \"minimum\" U[1,inf] !\"minimum\"]");
checkResult = modelchecker.check(*formula);
ASSERT_TRUE(checkResult->isExplicitQuantitativeCheckResult());
storm::modelchecker::ExplicitQuantitativeCheckResult<double> quantitativeCheckResult6 = checkResult->asExplicitQuantitativeCheckResult<double>();
EXPECT_NEAR(0.9999999033633374, quantitativeCheckResult6[initialState], storm::settings::generalSettings().getPrecision());
formula = formulaParser.parseFromString("R=? [C<=100]");
checkResult = modelchecker.check(*formula);
ASSERT_TRUE(checkResult->isExplicitQuantitativeCheckResult());
storm::modelchecker::ExplicitQuantitativeCheckResult<double> quantitativeCheckResult7 = checkResult->asExplicitQuantitativeCheckResult<double>();
EXPECT_NEAR(0.8602815057967503, quantitativeCheckResult7[initialState], storm::settings::generalSettings().getPrecision());
}
TEST(SparseCtmcCslModelCheckerTest, Embedded) {
Loading…
Cancel
Save