16 changed files with 430 additions and 45 deletions
-
9src/builder/ExplicitPrismModelBuilder.cpp
-
9src/counterexamples/MILPMinimalLabelSetGenerator.h
-
7src/counterexamples/SMTMinimalCommandSetGenerator.h
-
30src/modelchecker/AbstractModelChecker.cpp
-
160src/modelchecker/csl/SparseCtmcCslModelChecker.cpp
-
49src/modelchecker/csl/SparseCtmcCslModelChecker.h
-
25src/modelchecker/prctl/SparseDtmcPrctlModelChecker.cpp
-
6src/modelchecker/prctl/SparseDtmcPrctlModelChecker.h
-
12src/modelchecker/prctl/SparseMdpPrctlModelChecker.cpp
-
6src/modelchecker/prctl/SparseMdpPrctlModelChecker.h
-
8src/modelchecker/results/ExplicitQuantitativeCheckResult.cpp
-
7src/storage/expressions/BinaryNumericalFunctionExpression.cpp
-
2src/storage/expressions/UnaryNumericalFunctionExpression.cpp
-
6src/utility/cli.h
-
131src/utility/numerical.h
@ -0,0 +1,160 @@ |
|||||
|
#include "src/modelchecker/csl/SparseCtmcCslModelChecker.h"
|
||||
|
|
||||
|
#include <vector>
|
||||
|
|
||||
|
#include "src/utility/macros.h"
|
||||
|
#include "src/utility/vector.h"
|
||||
|
#include "src/utility/graph.h"
|
||||
|
#include "src/utility/solver.h"
|
||||
|
|
||||
|
#include "src/modelchecker/prctl/SparseDtmcPrctlModelChecker.h"
|
||||
|
#include "src/modelchecker/results/ExplicitQualitativeCheckResult.h"
|
||||
|
#include "src/modelchecker/results/ExplicitQuantitativeCheckResult.h"
|
||||
|
|
||||
|
#include "src/exceptions/InvalidStateException.h"
|
||||
|
#include "src/exceptions/InvalidPropertyException.h"
|
||||
|
#include "src/exceptions/NotImplementedException.h"
|
||||
|
|
||||
|
namespace storm { |
||||
|
namespace modelchecker { |
||||
|
template<class ValueType> |
||||
|
SparseCtmcCslModelChecker<ValueType>::SparseCtmcCslModelChecker(storm::models::sparse::Ctmc<ValueType> const& model) : SparsePropositionalModelChecker<ValueType>(model), linearEquationSolver(storm::utility::solver::getLinearEquationSolver<ValueType>()) { |
||||
|
// Intentionally left empty.
|
||||
|
} |
||||
|
|
||||
|
template<class ValueType> |
||||
|
SparseCtmcCslModelChecker<ValueType>::SparseCtmcCslModelChecker(storm::models::sparse::Ctmc<ValueType> const& model, std::unique_ptr<storm::solver::LinearEquationSolver<ValueType>>&& linearEquationSolver) : SparsePropositionalModelChecker<ValueType>(model), linearEquationSolver(std::move(linearEquationSolver)) { |
||||
|
// Intentionally left empty.
|
||||
|
} |
||||
|
|
||||
|
template<class ValueType> |
||||
|
bool SparseCtmcCslModelChecker<ValueType>::canHandle(storm::logic::Formula const& formula) const { |
||||
|
// FIXME: refine.
|
||||
|
return formula.isCslStateFormula() || formula.isCslPathFormula() || formula.isRewardPathFormula(); |
||||
|
} |
||||
|
|
||||
|
template<class ValueType> |
||||
|
std::unique_ptr<CheckResult> SparseCtmcCslModelChecker<ValueType>::computeBoundedUntilProbabilities(storm::logic::BoundedUntilFormula const& pathFormula, bool qualitative, boost::optional<storm::logic::OptimalityType> const& optimalityType) { |
||||
|
STORM_LOG_THROW(pathFormula.isIntervalBounded(), storm::exceptions::InvalidPropertyException, "Cannot treat non-interval bounded until."); |
||||
|
|
||||
|
std::unique_ptr<CheckResult> leftResultPointer = this->check(pathFormula.getLeftSubformula()); |
||||
|
std::unique_ptr<CheckResult> rightResultPointer = this->check(pathFormula.getRightSubformula()); |
||||
|
ExplicitQualitativeCheckResult const& leftResult = leftResultPointer->asExplicitQualitativeCheckResult();; |
||||
|
ExplicitQualitativeCheckResult const& rightResult = rightResultPointer->asExplicitQualitativeCheckResult(); |
||||
|
std::pair<double, double> const& intervalBounds = pathFormula.getIntervalBounds(); |
||||
|
std::unique_ptr<CheckResult> result = std::unique_ptr<CheckResult>(new ExplicitQuantitativeCheckResult<ValueType>(this->computeBoundedUntilProbabilitiesHelper(leftResult.getTruthValuesVector(), rightResult.getTruthValuesVector(), qualitative, intervalBounds.first, intervalBounds.second))); |
||||
|
|
||||
|
return result; |
||||
|
} |
||||
|
|
||||
|
template<class ValueType> |
||||
|
std::unique_ptr<CheckResult> SparseCtmcCslModelChecker<ValueType>::computeNextProbabilities(storm::logic::NextFormula const& pathFormula, bool qualitative, boost::optional<storm::logic::OptimalityType> const& optimalityType) { |
||||
|
std::unique_ptr<CheckResult> subResultPointer = this->check(pathFormula.getSubformula()); |
||||
|
ExplicitQualitativeCheckResult const& subResult = subResultPointer->asExplicitQualitativeCheckResult(); |
||||
|
std::vector<ValueType> result = SparseDtmcPrctlModelChecker<ValueType>::computeNextProbabilitiesHelper(this->getModel().getTransitionMatrix(), subResult.getTruthValuesVector(), *this->linearEquationSolver); |
||||
|
return std::unique_ptr<CheckResult>(new ExplicitQuantitativeCheckResult<ValueType>(std::move(result))); |
||||
|
} |
||||
|
|
||||
|
template<class ValueType> |
||||
|
std::unique_ptr<CheckResult> SparseCtmcCslModelChecker<ValueType>::computeUntilProbabilities(storm::logic::UntilFormula const& pathFormula, bool qualitative, boost::optional<storm::logic::OptimalityType> const& optimalityType) { |
||||
|
std::unique_ptr<CheckResult> leftResultPointer = this->check(pathFormula.getLeftSubformula()); |
||||
|
std::unique_ptr<CheckResult> rightResultPointer = this->check(pathFormula.getRightSubformula()); |
||||
|
ExplicitQualitativeCheckResult const& leftResult = leftResultPointer->asExplicitQualitativeCheckResult(); |
||||
|
ExplicitQualitativeCheckResult const& rightResult = rightResultPointer->asExplicitQualitativeCheckResult(); |
||||
|
return std::unique_ptr<CheckResult>(new ExplicitQuantitativeCheckResult<ValueType>(this->computeUntilProbabilitiesHelper(this->getModel().getTransitionMatrix(), this->getModel().getBackwardTransitions(), leftResult.getTruthValuesVector(), rightResult.getTruthValuesVector(), qualitative, *this->linearEquationSolver))); |
||||
|
} |
||||
|
|
||||
|
template<class ValueType> |
||||
|
storm::models::sparse::Ctmc<ValueType> const& SparseCtmcCslModelChecker<ValueType>::getModel() const { |
||||
|
return this->template getModelAs<storm::models::sparse::Ctmc<ValueType>>(); |
||||
|
} |
||||
|
|
||||
|
template<class ValueType> |
||||
|
std::vector<ValueType> SparseCtmcCslModelChecker<ValueType>::computeBoundedUntilProbabilitiesHelper(storm::storage::BitVector const& phiStates, storm::storage::BitVector const& psiStates, std::vector<ValueType> const& exitRates, bool qualitative, double lowerBound, double upperBound) const { |
||||
|
// If the time bounds are [0, inf], we rather call untimed reachability.
|
||||
|
storm::utility::ConstantsComparator<ValueType> comparator; |
||||
|
if (comparator.isZero(lowerBound) && comparator.isInfinity(upperBound)) { |
||||
|
return std::unique_ptr<CheckResult>(new ExplicitQuantitativeCheckResult<ValueType>(this->computeUntilProbabilitiesHelper(this->getModel().getTransitionMatrix(), this->getModel().getBackwardTransitions(), phiStates, psiStates, qualitative, *this->linearEquationSolver))); |
||||
|
} |
||||
|
|
||||
|
// From this point on, we know that we have to solve a more complicated problem [t, t'] with either t != 0
|
||||
|
// or t' != inf.
|
||||
|
|
||||
|
// Create the result vector.
|
||||
|
std::vector<ValueType> result(this->getModel().getNumberOfStates(), storm::utility::zero<ValueType>()); |
||||
|
|
||||
|
// If we identify the states that have probability 0 of reaching the target states, we can exclude them from the
|
||||
|
// further computations.
|
||||
|
storm::storage::BitVector statesWithProbabilityGreater0 = storm::utility::graph::performProbGreater0(this->getModel(), this->getModel().getBackwardTransitions(), phiStates, psiStates); |
||||
|
STORM_LOG_INFO("Found " << statesWithProbabilityGreater0.getNumberOfSetBits() << " states with probability greater 0."); |
||||
|
storm::storage::BitVector statesWithProbabilityGreater0NonPsi = statesWithProbabilityGreater0 & ~psiStates; |
||||
|
STORM_LOG_INFO("Found " << statesWithProbabilityGreater0NonPsi.getNumberOfSetBits() << " 'maybe' states."); |
||||
|
|
||||
|
if (!statesWithProbabilityGreater0NonPsi.empty()) { |
||||
|
if (comparator.isZero(upperBound)) { |
||||
|
// In this case, the interval is of the form [0, 0].
|
||||
|
storm::utility::vector::setVectorValues<ValueType>(result, psiStates, storm::utility::one<ValueType>()); |
||||
|
} else { |
||||
|
// Find the maximal rate of all 'maybe' states to take it as the uniformization rate.
|
||||
|
ValueType uniformizationRate = 0; |
||||
|
for (auto const& state : statesWithProbabilityGreater0NonPsi) { |
||||
|
uniformizationRate = std::max(uniformizationRate, exitRates[state]); |
||||
|
} |
||||
|
STORM_LOG_THROW(uniformizationRate > 0, storm::exceptions::InvalidStateException, "The uniformization rate must be positive."); |
||||
|
|
||||
|
if (comparator.isZero(lowerBound)) { |
||||
|
// In this case, the interval is of the form [0, t].
|
||||
|
// Note that this excludes [0, inf] since this is untimed reachability and we considered this case earlier.
|
||||
|
|
||||
|
// Compute the uniformized matrix.
|
||||
|
storm::storage::SparseMatrix<ValueType> uniformizedMatrix = this->computeUniformizedMatrix(this->getModel().getTransitionMatrix(), statesWithProbabilityGreater0, psiStates, uniformizationRate, exitRates); |
||||
|
|
||||
|
} else if (comparator.isInfinity(upperBound)) { |
||||
|
// In this case, the interval is of the form [t, inf] with t != 0.
|
||||
|
|
||||
|
} else { |
||||
|
// In this case, the interval is of the form [t, t'] with t != 0 and t' != inf.
|
||||
|
|
||||
|
} |
||||
|
} |
||||
|
} |
||||
|
|
||||
|
return result; |
||||
|
} |
||||
|
|
||||
|
template<class ValueType> |
||||
|
storm::storage::SparseMatrix<ValueType> SparseCtmcCslModelChecker<ValueType>::computeUniformizedMatrix(storm::storage::SparseMatrix<ValueType> const& transitionMatrix, storm::storage::BitVector const& maybeStates, storm::storage::BitVector const& absorbingStates, ValueType uniformizationRate, std::vector<ValueType> const& exitRates) { |
||||
|
// Create the submatrix that only contains the states with a positive probability (including the
|
||||
|
// psi states) and reserve space for elements on the diagonal.
|
||||
|
storm::storage::SparseMatrix<ValueType> uniformizedMatrix = transitionMatrix.getSubmatrix(false, maybeStates, maybeStates, true); |
||||
|
|
||||
|
// Make the appropriate states absorbing.
|
||||
|
uniformizedMatrix.makeRowsAbsorbing(absorbingStates % maybeStates); |
||||
|
|
||||
|
// Now we need to perform the actual uniformization. That is, all entries need to be divided by
|
||||
|
// the uniformization rate, and the diagonal needs to be set to the negative exit rate of the
|
||||
|
// state plus the self-loop rate and then increased by one.
|
||||
|
uint_fast64_t currentRow = 0; |
||||
|
for (auto const& state : maybeStates) { |
||||
|
for (auto& element : uniformizedMatrix.getRow(currentRow)) { |
||||
|
if (element.getColumn() == currentRow) { |
||||
|
if (absorbingStates.get(state)) { |
||||
|
// Nothing to do here, since the state has already been made absorbing.
|
||||
|
} else { |
||||
|
element.setValue(-(exitRates[state] + element.getValue()) / uniformizationRate + storm::utility::one<ValueType>()); |
||||
|
} |
||||
|
} else { |
||||
|
element.setValue(element.getValue() / uniformizationRate); |
||||
|
} |
||||
|
} |
||||
|
++currentRow; |
||||
|
} |
||||
|
} |
||||
|
|
||||
|
template<class ValueType> |
||||
|
std::vector<ValueType> SparseCtmcCslModelChecker<ValueType>::computeUntilProbabilitiesHelper(storm::storage::SparseMatrix<ValueType> const& transitionMatrix, storm::storage::SparseMatrix<ValueType> const& backwardTransitions, storm::storage::BitVector const& phiStates, storm::storage::BitVector const& psiStates, bool qualitative, storm::solver::LinearEquationSolver<ValueType> const& linearEquationSolver) { |
||||
|
return SparseDtmcPrctlModelChecker<ValueType>::computeUntilProbabilitiesHelper(transitionMatrix, backwardTransitions, phiStates, psiStates, qualitative, linearEquationSolver); |
||||
|
} |
||||
|
|
||||
|
} // namespace modelchecker
|
||||
|
} // namespace storm
|
@ -0,0 +1,49 @@ |
|||||
|
#ifndef STORM_MODELCHECKER_SPARSECTMCCSLMODELCHECKER_H_ |
||||
|
#define STORM_MODELCHECKER_SPARSECTMCCSLMODELCHECKER_H_ |
||||
|
|
||||
|
#include "src/modelchecker/propositional/SparsePropositionalModelChecker.h" |
||||
|
#include "src/models/sparse/Ctmc.h" |
||||
|
#include "src/solver/LinearEquationSolver.h" |
||||
|
|
||||
|
namespace storm { |
||||
|
namespace modelchecker { |
||||
|
|
||||
|
template<class ValueType> |
||||
|
class SparseCtmcCslModelChecker : public SparsePropositionalModelChecker<ValueType> { |
||||
|
public: |
||||
|
explicit SparseCtmcCslModelChecker(storm::models::sparse::Ctmc<ValueType> const& model); |
||||
|
explicit SparseCtmcCslModelChecker(storm::models::sparse::Ctmc<ValueType> const& model, std::unique_ptr<storm::solver::LinearEquationSolver<ValueType>>&& linearEquationSolver); |
||||
|
|
||||
|
// The implemented methods of the AbstractModelChecker interface. |
||||
|
virtual bool canHandle(storm::logic::Formula const& formula) const override; |
||||
|
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; |
||||
|
|
||||
|
protected: |
||||
|
storm::models::sparse::Ctmc<ValueType> const& getModel() const override; |
||||
|
|
||||
|
private: |
||||
|
// The methods that perform the actual checking. |
||||
|
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; |
||||
|
static std::vector<ValueType> computeUntilProbabilitiesHelper(storm::storage::SparseMatrix<ValueType> const& transitionMatrix, storm::storage::SparseMatrix<ValueType> const& backwardTransitions, storm::storage::BitVector const& phiStates, storm::storage::BitVector const& psiStates, bool qualitative, storm::solver::LinearEquationSolver<ValueType> const& linearEquationSolver); |
||||
|
|
||||
|
/*! |
||||
|
* 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 absorbingStates The states that need to be made absorbing. |
||||
|
* @param uniformizationRate The rate to be used for uniformization. |
||||
|
* @param exitRates The exit rates of all states. |
||||
|
*/ |
||||
|
storm::storage::SparseMatrix<ValueType> computeUniformizedMatrix(storm::storage::SparseMatrix<ValueType> const& transitionMatrix, storm::storage::BitVector const& maybeStates, storm::storage::BitVector const& absorbingStates, ValueType uniformizationRate, std::vector<ValueType> const& exitRates); |
||||
|
|
||||
|
// An object that is used for solving linear equations and performing matrix-vector multiplication. |
||||
|
std::unique_ptr<storm::solver::LinearEquationSolver<ValueType>> linearEquationSolver; |
||||
|
}; |
||||
|
|
||||
|
} // namespace modelchecker |
||||
|
} // namespace storm |
||||
|
|
||||
|
#endif /* STORM_MODELCHECKER_SPARSECTMCCSLMODELCHECKER_H_ */ |
@ -0,0 +1,131 @@ |
|||||
|
#include <vector> |
||||
|
#include <tuple> |
||||
|
#include <cmath> |
||||
|
|
||||
|
#include "src/utility/macros.h" |
||||
|
#include "src/utility/constants.h" |
||||
|
#include "src/exceptions/InvalidArgumentException.h" |
||||
|
|
||||
|
namespace storm { |
||||
|
namespace utility { |
||||
|
namespace numerical { |
||||
|
template<typename ValueType> |
||||
|
std::tuple<uint_fast64_t, uint_fast64_t, ValueType, std::vector<ValueType>> getFoxGlynnCutoff(ValueType lambda, ValueType underflow, ValueType overflow, ValueType accuracy) { |
||||
|
storm::utility::ConstantsComparator<ValueType> comparator; |
||||
|
STORM_LOG_THROW(!comparator.isZero(lambda), storm::exceptions::InvalidArgumentException, "Error in Fox-Glynn algorithm: lambda must not be zero."); |
||||
|
|
||||
|
// This code is a modified version of the one in PRISM. According to their implementation, for lambda |
||||
|
// smaller than 400, we compute the result using the naive method. |
||||
|
if (lambda < 25) { |
||||
|
ValueType eToPowerMinusLambda = std::exp(-lambda); |
||||
|
ValueType targetValue = (1 - accuracy) / eToPowerMinusLambda; |
||||
|
std::vector<ValueType> weights; |
||||
|
|
||||
|
ValueType exactlyKEvents = 1; |
||||
|
ValueType atMostKEvents = exactlyKEvents; |
||||
|
weights.push_back(exactlyKEvents * eToPowerMinusLambda); |
||||
|
|
||||
|
uint_fast64_t k = 1; |
||||
|
do { |
||||
|
exactlyKEvents *= lambda / k; |
||||
|
atMostKEvents += exactlyKEvents; |
||||
|
weights.push_back(exactlyKEvents * eToPowerMinusLambda); |
||||
|
++k; |
||||
|
} while (atMostKEvents < targetValue); |
||||
|
|
||||
|
return std::make_tuple(0, k - 1, 1.0, weights); |
||||
|
} else { |
||||
|
STORM_LOG_THROW(accuracy >= 1e-10, storm::exceptions::InvalidArgumentException, "Error in Fox-Glynn algorithm: the accuracy must not be below 1e-10."); |
||||
|
|
||||
|
// Factor from Fox&Glynn's paper. The paper does not explain where it comes from. |
||||
|
ValueType factor = 1e+10; |
||||
|
|
||||
|
// Now start the Finder algorithm to find the truncation points. |
||||
|
ValueType m = std::floor(lambda); |
||||
|
uint_fast64_t leftTruncationPoint = 0, rightTruncationPoint = 0; |
||||
|
{ |
||||
|
// Factors used by the corollaries explained in Fox & Glynns paper. |
||||
|
// Square root of pi. |
||||
|
ValueType sqrtpi = 1.77245385090551602729; |
||||
|
|
||||
|
// Square root of 2. |
||||
|
ValueType sqrt2 = 1.41421356237309504880; |
||||
|
|
||||
|
// Set up a_\lambda, b_\lambda, and the square root of lambda. |
||||
|
ValueType aLambda = 0, bLambda = 0, sqrtLambda = 0; |
||||
|
if (m < 400) { |
||||
|
sqrtLambda = std::sqrt(400.0); |
||||
|
aLambda = (1.0 + 1.0 / 400.0) * exp(0.0625) * sqrt2; |
||||
|
bLambda = (1.0 + 1.0 / 400.0) * exp(0.125 / 400.0); |
||||
|
} else { |
||||
|
sqrtLambda = std::sqrt(lambda); |
||||
|
aLambda = (1.0 + 1.0 / lambda) * exp(0.0625) * sqrt2; |
||||
|
bLambda = (1.0 + 1.0 / lambda) * exp(0.125 / lambda); |
||||
|
} |
||||
|
|
||||
|
// Use Corollary 1 from the paper to find the right truncation point. |
||||
|
uint_fast64_t k = 3; |
||||
|
|
||||
|
ValueType dkl = 1.0 / (1 - std::exp(-(2.0 / 9.0) * (k * sqrt2 * sqrtLambda + 1.5))); |
||||
|
|
||||
|
// According to David Jansen the upper bound can be ignored to achieve more accurate results. |
||||
|
// Right hand side of the equation in Corollary 1. |
||||
|
while ((accuracy / 2.0) < (aLambda * dkl * std::exp(-(k*k / 2.0)) / (k * sqrt2 * sqrtpi))) { |
||||
|
++k; |
||||
|
|
||||
|
// d(k,Lambda) from the paper. |
||||
|
dkl = 1.0 / (1 - std::exp(-(2.0 / 9.0)*(k * sqrt2 * sqrtLambda + 1.5))); |
||||
|
} |
||||
|
|
||||
|
// Left hand side of the equation in Corollary 1. |
||||
|
rightTruncationPoint = static_cast<uint_fast64_t>(std::ceil((m + k * sqrt2 * sqrtLambda + 1.5))); |
||||
|
|
||||
|
// Use Corollary 2 to find left truncation point. |
||||
|
k = 3; |
||||
|
|
||||
|
// Right hand side of the equation in Corollary 2. |
||||
|
while ((accuracy / 2.0) < ((bLambda * exp(-(k*k / 2.0))) / (k * sqrt2 * sqrtpi))) { |
||||
|
++k; |
||||
|
} |
||||
|
|
||||
|
// Left hand side of the equation in Corollary 2. |
||||
|
leftTruncationPoint = std::max(static_cast<uint_fast64_t>(std::trunc(m - k * sqrtLambda - 1.5)), uint_fast64_t(0)); |
||||
|
|
||||
|
// Check for underflow. |
||||
|
STORM_LOG_THROW(static_cast<uint_fast64_t>(std::trunc((m - k * sqrtLambda - 1.5))) >= 0, storm::exceptions::OutOfRangeException, "Error in Fox-Glynn algorithm: Underflow of left truncation point."); |
||||
|
} |
||||
|
|
||||
|
std::vector<ValueType> weights(rightTruncationPoint - leftTruncationPoint + 1); |
||||
|
weights[m - leftTruncationPoint] = overflow / (factor * (rightTruncationPoint - leftTruncationPoint)); |
||||
|
for (uint_fast64_t j = m; j > leftTruncationPoint; --j) { |
||||
|
weights[j - 1 - leftTruncationPoint] = (j / lambda) * weights[j - leftTruncationPoint]; |
||||
|
} |
||||
|
|
||||
|
for (uint_fast64_t j = m; j < rightTruncationPoint; ++j) { |
||||
|
weights[j + 1 - leftTruncationPoint] = (lambda / (j + 1)) * weights[j - leftTruncationPoint]; |
||||
|
} |
||||
|
|
||||
|
|
||||
|
// Compute the total weight and start with smallest to avoid roundoff errors. |
||||
|
ValueType totalWeight = storm::utility::zero<ValueType>(); |
||||
|
|
||||
|
uint_fast64_t s = leftTruncationPoint; |
||||
|
uint_fast64_t t = rightTruncationPoint; |
||||
|
while (s < t) { |
||||
|
if (weights[s - leftTruncationPoint] <= weights[t - leftTruncationPoint]) { |
||||
|
totalWeight += weights[s - leftTruncationPoint]; |
||||
|
++s; |
||||
|
} else { |
||||
|
totalWeight += weights[t - leftTruncationPoint]; |
||||
|
--t; |
||||
|
} |
||||
|
} |
||||
|
|
||||
|
totalWeight += weights[s - leftTruncationPoint]; |
||||
|
|
||||
|
return std::make_tuple(leftTruncationPoint, rightTruncationPoint, totalWeight, weights); |
||||
|
} |
||||
|
} |
||||
|
} |
||||
|
} |
||||
|
} |
Reference in new issue
xxxxxxxxxx