You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

682 lines
50 KiB

  1. #include "src/modelchecker/csl/helper/SparseCtmcCslHelper.h"
  2. #include "src/modelchecker/prctl/helper/SparseDtmcPrctlHelper.h"
  3. #include "src/modelchecker/reachability/SparseDtmcEliminationModelChecker.h"
  4. #include "src/models/sparse/StandardRewardModel.h"
  5. #include "src/settings/SettingsManager.h"
  6. #include "src/settings/modules/GeneralSettings.h"
  7. #include "src/solver/LinearEquationSolver.h"
  8. #include "src/storage/StronglyConnectedComponentDecomposition.h"
  9. #include "src/adapters/CarlAdapter.h"
  10. #include "src/utility/macros.h"
  11. #include "src/utility/vector.h"
  12. #include "src/utility/graph.h"
  13. #include "src/utility/numerical.h"
  14. #include "src/exceptions/InvalidStateException.h"
  15. #include "src/exceptions/InvalidPropertyException.h"
  16. namespace storm {
  17. namespace modelchecker {
  18. namespace helper {
  19. template <typename ValueType>
  20. std::vector<ValueType> SparseCtmcCslHelper<ValueType>::computeBoundedUntilProbabilities(storm::storage::SparseMatrix<ValueType> const& rateMatrix, storm::storage::SparseMatrix<ValueType> const& backwardTransitions, storm::storage::BitVector const& phiStates, storm::storage::BitVector const& psiStates, std::vector<ValueType> const& exitRates, bool qualitative, double lowerBound, double upperBound, storm::solver::LinearEquationSolverFactory<ValueType> const& linearEquationSolverFactory) {
  21. uint_fast64_t numberOfStates = rateMatrix.getRowCount();
  22. // If the time bounds are [0, inf], we rather call untimed reachability.
  23. if (storm::utility::isZero(lowerBound) && upperBound == storm::utility::infinity<ValueType>()) {
  24. return computeUntilProbabilities(rateMatrix, backwardTransitions, exitRates, phiStates, psiStates, qualitative, linearEquationSolverFactory);
  25. }
  26. // From this point on, we know that we have to solve a more complicated problem [t, t'] with either t != 0
  27. // or t' != inf.
  28. // Create the result vector.
  29. std::vector<ValueType> result;
  30. // If we identify the states that have probability 0 of reaching the target states, we can exclude them from the
  31. // further computations.
  32. storm::storage::BitVector statesWithProbabilityGreater0 = storm::utility::graph::performProbGreater0(backwardTransitions, phiStates, psiStates);
  33. STORM_LOG_INFO("Found " << statesWithProbabilityGreater0.getNumberOfSetBits() << " states with probability greater 0.");
  34. storm::storage::BitVector statesWithProbabilityGreater0NonPsi = statesWithProbabilityGreater0 & ~psiStates;
  35. STORM_LOG_INFO("Found " << statesWithProbabilityGreater0NonPsi.getNumberOfSetBits() << " 'maybe' states.");
  36. if (!statesWithProbabilityGreater0NonPsi.empty()) {
  37. if (storm::utility::isZero(upperBound)) {
  38. // In this case, the interval is of the form [0, 0].
  39. result = std::vector<ValueType>(numberOfStates, storm::utility::zero<ValueType>());
  40. storm::utility::vector::setVectorValues<ValueType>(result, psiStates, storm::utility::one<ValueType>());
  41. } else {
  42. if (storm::utility::isZero(lowerBound)) {
  43. // In this case, the interval is of the form [0, t].
  44. // Note that this excludes [0, inf] since this is untimed reachability and we considered this case earlier.
  45. // Find the maximal rate of all 'maybe' states to take it as the uniformization rate.
  46. ValueType uniformizationRate = 0;
  47. for (auto const& state : statesWithProbabilityGreater0NonPsi) {
  48. uniformizationRate = std::max(uniformizationRate, exitRates[state]);
  49. }
  50. uniformizationRate *= 1.02;
  51. STORM_LOG_THROW(uniformizationRate > 0, storm::exceptions::InvalidStateException, "The uniformization rate must be positive.");
  52. // Compute the uniformized matrix.
  53. storm::storage::SparseMatrix<ValueType> uniformizedMatrix = computeUniformizedMatrix(rateMatrix, statesWithProbabilityGreater0NonPsi, uniformizationRate, exitRates);
  54. // Compute the vector that is to be added as a compensation for removing the absorbing states.
  55. std::vector<ValueType> b = rateMatrix.getConstrainedRowSumVector(statesWithProbabilityGreater0NonPsi, psiStates);
  56. for (auto& element : b) {
  57. element /= uniformizationRate;
  58. }
  59. // Finally compute the transient probabilities.
  60. std::vector<ValueType> values(statesWithProbabilityGreater0NonPsi.getNumberOfSetBits(), storm::utility::zero<ValueType>());
  61. std::vector<ValueType> subresult = computeTransientProbabilities(uniformizedMatrix, &b, upperBound, uniformizationRate, values, linearEquationSolverFactory);
  62. result = std::vector<ValueType>(numberOfStates, storm::utility::zero<ValueType>());
  63. storm::utility::vector::setVectorValues(result, statesWithProbabilityGreater0NonPsi, subresult);
  64. storm::utility::vector::setVectorValues(result, psiStates, storm::utility::one<ValueType>());
  65. } else if (upperBound == storm::utility::infinity<ValueType>()) {
  66. // In this case, the interval is of the form [t, inf] with t != 0.
  67. // Start by computing the (unbounded) reachability probabilities of reaching psi states while
  68. // staying in phi states.
  69. result = computeUntilProbabilities(rateMatrix, backwardTransitions, exitRates, phiStates, psiStates, qualitative, linearEquationSolverFactory);
  70. // Determine the set of states that must be considered further.
  71. storm::storage::BitVector relevantStates = statesWithProbabilityGreater0 & phiStates;
  72. std::vector<ValueType> subResult(relevantStates.getNumberOfSetBits());
  73. storm::utility::vector::selectVectorValues(subResult, relevantStates, result);
  74. ValueType uniformizationRate = 0;
  75. for (auto const& state : relevantStates) {
  76. uniformizationRate = std::max(uniformizationRate, exitRates[state]);
  77. }
  78. uniformizationRate *= 1.02;
  79. STORM_LOG_THROW(uniformizationRate > 0, storm::exceptions::InvalidStateException, "The uniformization rate must be positive.");
  80. // Compute the uniformized matrix.
  81. storm::storage::SparseMatrix<ValueType> uniformizedMatrix = computeUniformizedMatrix(rateMatrix, relevantStates, uniformizationRate, exitRates);
  82. // Compute the transient probabilities.
  83. subResult = computeTransientProbabilities(uniformizedMatrix, nullptr, lowerBound, uniformizationRate, subResult, linearEquationSolverFactory);
  84. // Fill in the correct values.
  85. storm::utility::vector::setVectorValues(result, ~relevantStates, storm::utility::zero<ValueType>());
  86. storm::utility::vector::setVectorValues(result, relevantStates, subResult);
  87. } else {
  88. // In this case, the interval is of the form [t, t'] with t != 0 and t' != inf.
  89. if (lowerBound != upperBound) {
  90. // In this case, the interval is of the form [t, t'] with t != 0, t' != inf and t != t'.
  91. // Find the maximal rate of all 'maybe' states to take it as the uniformization rate.
  92. ValueType uniformizationRate = storm::utility::zero<ValueType>();
  93. for (auto const& state : statesWithProbabilityGreater0NonPsi) {
  94. uniformizationRate = std::max(uniformizationRate, exitRates[state]);
  95. }
  96. uniformizationRate *= 1.02;
  97. STORM_LOG_THROW(uniformizationRate > 0, storm::exceptions::InvalidStateException, "The uniformization rate must be positive.");
  98. // Compute the (first) uniformized matrix.
  99. storm::storage::SparseMatrix<ValueType> uniformizedMatrix = computeUniformizedMatrix(rateMatrix, statesWithProbabilityGreater0NonPsi, uniformizationRate, exitRates);
  100. // Compute the vector that is to be added as a compensation for removing the absorbing states.
  101. std::vector<ValueType> b = rateMatrix.getConstrainedRowSumVector(statesWithProbabilityGreater0NonPsi, psiStates);
  102. for (auto& element : b) {
  103. element /= uniformizationRate;
  104. }
  105. // Start by computing the transient probabilities of reaching a psi state in time t' - t.
  106. std::vector<ValueType> values(statesWithProbabilityGreater0NonPsi.getNumberOfSetBits(), storm::utility::zero<ValueType>());
  107. std::vector<ValueType> subresult = computeTransientProbabilities(uniformizedMatrix, &b, upperBound - lowerBound, uniformizationRate, values, linearEquationSolverFactory);
  108. storm::storage::BitVector relevantStates = statesWithProbabilityGreater0 & phiStates;
  109. std::vector<ValueType> newSubresult = std::vector<ValueType>(relevantStates.getNumberOfSetBits());
  110. storm::utility::vector::setVectorValues(newSubresult, statesWithProbabilityGreater0NonPsi % relevantStates, subresult);
  111. storm::utility::vector::setVectorValues(newSubresult, psiStates % relevantStates, storm::utility::one<ValueType>());
  112. // Then compute the transient probabilities of being in such a state after t time units. For this,
  113. // we must re-uniformize the CTMC, so we need to compute the second uniformized matrix.
  114. uniformizationRate = storm::utility::zero<ValueType>();
  115. for (auto const& state : relevantStates) {
  116. uniformizationRate = std::max(uniformizationRate, exitRates[state]);
  117. }
  118. uniformizationRate *= 1.02;
  119. STORM_LOG_THROW(uniformizationRate > 0, storm::exceptions::InvalidStateException, "The uniformization rate must be positive.");
  120. // Finally, we compute the second set of transient probabilities.
  121. uniformizedMatrix = computeUniformizedMatrix(rateMatrix, relevantStates, uniformizationRate, exitRates);
  122. newSubresult = computeTransientProbabilities(uniformizedMatrix, nullptr, lowerBound, uniformizationRate, newSubresult, linearEquationSolverFactory);
  123. // Fill in the correct values.
  124. result = std::vector<ValueType>(numberOfStates, storm::utility::zero<ValueType>());
  125. storm::utility::vector::setVectorValues(result, ~relevantStates, storm::utility::zero<ValueType>());
  126. storm::utility::vector::setVectorValues(result, relevantStates, newSubresult);
  127. } else {
  128. // In this case, the interval is of the form [t, t] with t != 0, t != inf.
  129. std::vector<ValueType> newSubresult = std::vector<ValueType>(statesWithProbabilityGreater0.getNumberOfSetBits());
  130. storm::utility::vector::setVectorValues(newSubresult, psiStates % statesWithProbabilityGreater0, storm::utility::one<ValueType>());
  131. // Then compute the transient probabilities of being in such a state after t time units. For this,
  132. // we must re-uniformize the CTMC, so we need to compute the second uniformized matrix.
  133. ValueType uniformizationRate = storm::utility::zero<ValueType>();
  134. for (auto const& state : statesWithProbabilityGreater0) {
  135. uniformizationRate = std::max(uniformizationRate, exitRates[state]);
  136. }
  137. uniformizationRate *= 1.02;
  138. STORM_LOG_THROW(uniformizationRate > 0, storm::exceptions::InvalidStateException, "The uniformization rate must be positive.");
  139. // Finally, we compute the second set of transient probabilities.
  140. storm::storage::SparseMatrix<ValueType> uniformizedMatrix = computeUniformizedMatrix(rateMatrix, statesWithProbabilityGreater0, uniformizationRate, exitRates);
  141. newSubresult = computeTransientProbabilities(uniformizedMatrix, nullptr, lowerBound, uniformizationRate, newSubresult, linearEquationSolverFactory);
  142. // Fill in the correct values.
  143. result = std::vector<ValueType>(numberOfStates, storm::utility::zero<ValueType>());
  144. storm::utility::vector::setVectorValues(result, ~statesWithProbabilityGreater0, storm::utility::zero<ValueType>());
  145. storm::utility::vector::setVectorValues(result, statesWithProbabilityGreater0, newSubresult);
  146. }
  147. }
  148. }
  149. }
  150. return result;
  151. }
  152. template <typename ValueType>
  153. std::vector<ValueType> SparseCtmcCslHelper<ValueType>::computeUntilProbabilities(storm::storage::SparseMatrix<ValueType> const& rateMatrix, storm::storage::SparseMatrix<ValueType> const& backwardTransitions, std::vector<ValueType> const& exitRateVector, storm::storage::BitVector const& phiStates, storm::storage::BitVector const& psiStates, bool qualitative, storm::solver::LinearEquationSolverFactory<ValueType> const& linearEquationSolverFactory) {
  154. return SparseDtmcPrctlHelper<ValueType>::computeUntilProbabilities(computeProbabilityMatrix(rateMatrix, exitRateVector), backwardTransitions, phiStates, psiStates, qualitative, linearEquationSolverFactory);
  155. }
  156. template <typename ValueType>
  157. std::vector<ValueType> SparseCtmcCslHelper<ValueType>::computeNextProbabilities(storm::storage::SparseMatrix<ValueType> const& rateMatrix, std::vector<ValueType> const& exitRateVector, storm::storage::BitVector const& nextStates, storm::solver::LinearEquationSolverFactory<ValueType> const& linearEquationSolverFactory) {
  158. return SparseDtmcPrctlHelper<ValueType>::computeNextProbabilities(computeProbabilityMatrix(rateMatrix, exitRateVector), nextStates, linearEquationSolverFactory);
  159. }
  160. template <typename ValueType>
  161. storm::storage::SparseMatrix<ValueType> SparseCtmcCslHelper<ValueType>::computeUniformizedMatrix(storm::storage::SparseMatrix<ValueType> const& rateMatrix, storm::storage::BitVector const& maybeStates, ValueType uniformizationRate, std::vector<ValueType> const& exitRates) {
  162. STORM_LOG_DEBUG("Computing uniformized matrix using uniformization rate " << uniformizationRate << ".");
  163. STORM_LOG_DEBUG("Keeping " << maybeStates.getNumberOfSetBits() << " rows.");
  164. // Create the submatrix that only contains the states with a positive probability (including the
  165. // psi states) and reserve space for elements on the diagonal.
  166. storm::storage::SparseMatrix<ValueType> uniformizedMatrix = rateMatrix.getSubmatrix(false, maybeStates, maybeStates, true);
  167. // Now we need to perform the actual uniformization. That is, all entries need to be divided by
  168. // the uniformization rate, and the diagonal needs to be set to the negative exit rate of the
  169. // state plus the self-loop rate and then increased by one.
  170. uint_fast64_t currentRow = 0;
  171. for (auto const& state : maybeStates) {
  172. for (auto& element : uniformizedMatrix.getRow(currentRow)) {
  173. if (element.getColumn() == currentRow) {
  174. element.setValue((element.getValue() - exitRates[state]) / uniformizationRate + storm::utility::one<ValueType>());
  175. } else {
  176. element.setValue(element.getValue() / uniformizationRate);
  177. }
  178. }
  179. ++currentRow;
  180. }
  181. return uniformizedMatrix;
  182. }
  183. template <typename ValueType>
  184. template<bool computeCumulativeReward>
  185. std::vector<ValueType> SparseCtmcCslHelper<ValueType>::computeTransientProbabilities(storm::storage::SparseMatrix<ValueType> const& uniformizedMatrix, std::vector<ValueType> const* addVector, ValueType timeBound, ValueType uniformizationRate, std::vector<ValueType> values, storm::solver::LinearEquationSolverFactory<ValueType> const& linearEquationSolverFactory) {
  186. ValueType lambda = timeBound * uniformizationRate;
  187. // If no time can pass, the current values are the result.
  188. if (storm::utility::isZero(lambda)) {
  189. return values;
  190. }
  191. // Use Fox-Glynn to get the truncation points and the weights.
  192. std::tuple<uint_fast64_t, uint_fast64_t, ValueType, std::vector<ValueType>> foxGlynnResult = storm::utility::numerical::getFoxGlynnCutoff(lambda, 1e-300, 1e+300, storm::settings::getModule<storm::settings::modules::GeneralSettings>().getPrecision() / 8.0);
  193. STORM_LOG_DEBUG("Fox-Glynn cutoff points: left=" << std::get<0>(foxGlynnResult) << ", right=" << std::get<1>(foxGlynnResult));
  194. // Scale the weights so they add up to one.
  195. for (auto& element : std::get<3>(foxGlynnResult)) {
  196. element /= std::get<2>(foxGlynnResult);
  197. }
  198. // If the cumulative reward is to be computed, we need to adjust the weights.
  199. if (computeCumulativeReward) {
  200. ValueType sum = storm::utility::zero<ValueType>();
  201. for (auto& element : std::get<3>(foxGlynnResult)) {
  202. sum += element;
  203. element = (1 - sum) / uniformizationRate;
  204. }
  205. }
  206. STORM_LOG_DEBUG("Starting iterations with " << uniformizedMatrix.getRowCount() << " x " << uniformizedMatrix.getColumnCount() << " matrix.");
  207. // Initialize result.
  208. std::vector<ValueType> result;
  209. uint_fast64_t startingIteration = std::get<0>(foxGlynnResult);
  210. if (startingIteration == 0) {
  211. result = values;
  212. storm::utility::vector::scaleVectorInPlace(result, std::get<3>(foxGlynnResult)[0]);
  213. ++startingIteration;
  214. } else {
  215. if (computeCumulativeReward) {
  216. result = std::vector<ValueType>(values.size());
  217. std::function<ValueType (ValueType const&)> scaleWithUniformizationRate = [&uniformizationRate] (ValueType const& a) -> ValueType { return a / uniformizationRate; };
  218. storm::utility::vector::applyPointwise(values, result, scaleWithUniformizationRate);
  219. } else {
  220. result = std::vector<ValueType>(values.size());
  221. }
  222. }
  223. std::vector<ValueType> multiplicationResult(result.size());
  224. std::unique_ptr<storm::solver::LinearEquationSolver<ValueType>> solver = linearEquationSolverFactory.create(std::move(uniformizedMatrix));
  225. if (!computeCumulativeReward && std::get<0>(foxGlynnResult) > 1) {
  226. // Perform the matrix-vector multiplications (without adding).
  227. solver->performMatrixVectorMultiplication(values, addVector, std::get<0>(foxGlynnResult) - 1, &multiplicationResult);
  228. } else if (computeCumulativeReward) {
  229. std::function<ValueType(ValueType const&, ValueType const&)> addAndScale = [&uniformizationRate] (ValueType const& a, ValueType const& b) { return a + b / uniformizationRate; };
  230. // For the iterations below the left truncation point, we need to add and scale the result with the uniformization rate.
  231. for (uint_fast64_t index = 1; index < startingIteration; ++index) {
  232. solver->performMatrixVectorMultiplication(values, nullptr, 1, &multiplicationResult);
  233. storm::utility::vector::applyPointwise(result, values, result, addAndScale);
  234. }
  235. }
  236. // For the indices that fall in between the truncation points, we need to perform the matrix-vector
  237. // multiplication, scale and add the result.
  238. ValueType weight = 0;
  239. std::function<ValueType(ValueType const&, ValueType const&)> addAndScale = [&weight] (ValueType const& a, ValueType const& b) { return a + weight * b; };
  240. for (uint_fast64_t index = startingIteration; index <= std::get<1>(foxGlynnResult); ++index) {
  241. solver->performMatrixVectorMultiplication(values, addVector, 1, &multiplicationResult);
  242. weight = std::get<3>(foxGlynnResult)[index - std::get<0>(foxGlynnResult)];
  243. storm::utility::vector::applyPointwise(result, values, result, addAndScale);
  244. }
  245. return result;
  246. }
  247. template <typename ValueType>
  248. template <typename RewardModelType>
  249. std::vector<ValueType> SparseCtmcCslHelper<ValueType>::computeInstantaneousRewards(storm::storage::SparseMatrix<ValueType> const& rateMatrix, std::vector<ValueType> const& exitRateVector, RewardModelType const& rewardModel, double timeBound, storm::solver::LinearEquationSolverFactory<ValueType> const& linearEquationSolverFactory) {
  250. // Only compute the result if the model has a state-based reward this->getModel().
  251. STORM_LOG_THROW(!rewardModel.empty(), storm::exceptions::InvalidPropertyException, "Missing reward model for formula. Skipping formula.");
  252. uint_fast64_t numberOfStates = rateMatrix.getRowCount();
  253. // Initialize result to state rewards of the this->getModel().
  254. std::vector<ValueType> result(rewardModel.getStateRewardVector());
  255. // If the time-bound is not zero, we need to perform a transient analysis.
  256. if (timeBound > 0) {
  257. ValueType uniformizationRate = 0;
  258. for (auto const& rate : exitRateVector) {
  259. uniformizationRate = std::max(uniformizationRate, rate);
  260. }
  261. uniformizationRate *= 1.02;
  262. STORM_LOG_THROW(uniformizationRate > 0, storm::exceptions::InvalidStateException, "The uniformization rate must be positive.");
  263. storm::storage::SparseMatrix<ValueType> uniformizedMatrix = computeUniformizedMatrix(rateMatrix, storm::storage::BitVector(numberOfStates, true), uniformizationRate, exitRateVector);
  264. result = computeTransientProbabilities(uniformizedMatrix, nullptr, timeBound, uniformizationRate, result, linearEquationSolverFactory);
  265. }
  266. return result;
  267. }
  268. template <typename ValueType>
  269. template <typename RewardModelType>
  270. std::vector<ValueType> SparseCtmcCslHelper<ValueType>::computeCumulativeRewards(storm::storage::SparseMatrix<ValueType> const& rateMatrix, std::vector<ValueType> const& exitRateVector, RewardModelType const& rewardModel, double timeBound, storm::solver::LinearEquationSolverFactory<ValueType> const& linearEquationSolverFactory) {
  271. // Only compute the result if the model has a state-based reward this->getModel().
  272. STORM_LOG_THROW(!rewardModel.empty(), storm::exceptions::InvalidPropertyException, "Missing reward model for formula. Skipping formula.");
  273. uint_fast64_t numberOfStates = rateMatrix.getRowCount();
  274. // If the time bound is zero, the result is the constant zero vector.
  275. if (timeBound == 0) {
  276. return std::vector<ValueType>(numberOfStates, storm::utility::zero<ValueType>());
  277. }
  278. // Otherwise, we need to perform some computations.
  279. // Start with the uniformization.
  280. ValueType uniformizationRate = 0;
  281. for (auto const& rate : exitRateVector) {
  282. uniformizationRate = std::max(uniformizationRate, rate);
  283. }
  284. uniformizationRate *= 1.02;
  285. STORM_LOG_THROW(uniformizationRate > 0, storm::exceptions::InvalidStateException, "The uniformization rate must be positive.");
  286. storm::storage::SparseMatrix<ValueType> uniformizedMatrix = computeUniformizedMatrix(rateMatrix, storm::storage::BitVector(numberOfStates, true), uniformizationRate, exitRateVector);
  287. // Compute the total state reward vector.
  288. std::vector<ValueType> totalRewardVector = rewardModel.getTotalRewardVector(rateMatrix, exitRateVector);
  289. // Finally, compute the transient probabilities.
  290. return computeTransientProbabilities<true>(uniformizedMatrix, nullptr, timeBound, uniformizationRate, totalRewardVector, linearEquationSolverFactory);
  291. }
  292. template <typename ValueType>
  293. template <typename RewardModelType>
  294. std::vector<ValueType> SparseCtmcCslHelper<ValueType>::computeReachabilityRewards(storm::storage::SparseMatrix<ValueType> const& rateMatrix, storm::storage::SparseMatrix<ValueType> const& backwardTransitions, std::vector<ValueType> const& exitRateVector, RewardModelType const& rewardModel, storm::storage::BitVector const& targetStates, bool qualitative, storm::solver::LinearEquationSolverFactory<ValueType> const& linearEquationSolverFactory) {
  295. STORM_LOG_THROW(!rewardModel.empty(), storm::exceptions::InvalidPropertyException, "Missing reward model for formula. Skipping formula.");
  296. storm::storage::SparseMatrix<ValueType> probabilityMatrix = computeProbabilityMatrix(rateMatrix, exitRateVector);
  297. std::vector<ValueType> totalRewardVector;
  298. if (rewardModel.hasStateRewards()) {
  299. totalRewardVector = rewardModel.getStateRewardVector();
  300. typename std::vector<ValueType>::const_iterator it2 = exitRateVector.begin();
  301. for (typename std::vector<ValueType>::iterator it1 = totalRewardVector.begin(), ite1 = totalRewardVector.end(); it1 != ite1; ++it1, ++it2) {
  302. *it1 /= *it2;
  303. }
  304. if (rewardModel.hasStateActionRewards()) {
  305. storm::utility::vector::addVectors(totalRewardVector, rewardModel.getStateActionRewardVector(), totalRewardVector);
  306. }
  307. if (rewardModel.hasTransitionRewards()) {
  308. storm::utility::vector::addVectors(totalRewardVector, probabilityMatrix.getPointwiseProductRowSumVector(rewardModel.getTransitionRewardMatrix()), totalRewardVector);
  309. }
  310. } else if (rewardModel.hasTransitionRewards()) {
  311. totalRewardVector = probabilityMatrix.getPointwiseProductRowSumVector(rewardModel.getTransitionRewardMatrix());
  312. if (rewardModel.hasStateActionRewards()) {
  313. storm::utility::vector::addVectors(totalRewardVector, rewardModel.getStateActionRewardVector(), totalRewardVector);
  314. }
  315. } else {
  316. totalRewardVector = rewardModel.getStateActionRewardVector();
  317. }
  318. return storm::modelchecker::helper::SparseDtmcPrctlHelper<ValueType>::computeReachabilityRewards(probabilityMatrix, backwardTransitions, totalRewardVector, targetStates, qualitative, linearEquationSolverFactory);
  319. }
  320. template <typename ValueType>
  321. storm::storage::SparseMatrix<ValueType> SparseCtmcCslHelper<ValueType>::computeProbabilityMatrix(storm::storage::SparseMatrix<ValueType> const& rateMatrix, std::vector<ValueType> const& exitRates) {
  322. // Turn the rates into probabilities by scaling each row with the exit rate of the state.
  323. storm::storage::SparseMatrix<ValueType> result(rateMatrix);
  324. for (uint_fast64_t row = 0; row < result.getRowCount(); ++row) {
  325. for (auto& entry : result.getRow(row)) {
  326. entry.setValue(entry.getValue() / exitRates[row]);
  327. }
  328. }
  329. return result;
  330. }
  331. template <typename ValueType>
  332. storm::storage::SparseMatrix<ValueType> SparseCtmcCslHelper<ValueType>::computeGeneratorMatrix(storm::storage::SparseMatrix<ValueType> const& rateMatrix, std::vector<ValueType> const& exitRates) {
  333. storm::storage::SparseMatrix<ValueType> generatorMatrix(rateMatrix, true);
  334. // Place the negative exit rate on the diagonal.
  335. for (uint_fast64_t row = 0; row < generatorMatrix.getRowCount(); ++row) {
  336. for (auto& entry : generatorMatrix.getRow(row)) {
  337. if (entry.getColumn() == row) {
  338. entry.setValue(-exitRates[row]);
  339. }
  340. }
  341. }
  342. return generatorMatrix;
  343. }
  344. template <typename ValueType>
  345. std::vector<ValueType> SparseCtmcCslHelper<ValueType>::computeLongRunAverageProbabilities(storm::storage::SparseMatrix<ValueType> const& probabilityMatrix, storm::storage::BitVector const& psiStates, std::vector<ValueType> const* exitRateVector, bool qualitative, storm::solver::LinearEquationSolverFactory<ValueType> const& linearEquationSolverFactory) {
  346. // If there are no goal states, we avoid the computation and directly return zero.
  347. uint_fast64_t numberOfStates = probabilityMatrix.getRowCount();
  348. if (psiStates.empty()) {
  349. return std::vector<ValueType>(numberOfStates, storm::utility::zero<ValueType>());
  350. }
  351. // Likewise, if all bits are set, we can avoid the computation.
  352. if (psiStates.full()) {
  353. return std::vector<ValueType>(numberOfStates, storm::utility::one<ValueType>());
  354. }
  355. // Start by decomposing the DTMC into its BSCCs.
  356. storm::storage::StronglyConnectedComponentDecomposition<ValueType> bsccDecomposition(probabilityMatrix, storm::storage::BitVector(probabilityMatrix.getRowCount(), true), false, true);
  357. STORM_LOG_DEBUG("Found " << bsccDecomposition.size() << " BSCCs.");
  358. // Get some data members for convenience.
  359. ValueType one = storm::utility::one<ValueType>();
  360. ValueType zero = storm::utility::zero<ValueType>();
  361. // Prepare the vector holding the LRA values for each of the BSCCs.
  362. std::vector<ValueType> bsccLra(bsccDecomposition.size(), zero);
  363. // First we check which states are in BSCCs.
  364. storm::storage::BitVector statesInBsccs(numberOfStates);
  365. storm::storage::BitVector firstStatesInBsccs(numberOfStates);
  366. for (uint_fast64_t currentBsccIndex = 0; currentBsccIndex < bsccDecomposition.size(); ++currentBsccIndex) {
  367. storm::storage::StronglyConnectedComponent const& bscc = bsccDecomposition[currentBsccIndex];
  368. // Gather information for later use.
  369. bool first = true;
  370. for (auto const& state : bscc) {
  371. statesInBsccs.set(state);
  372. if (first) {
  373. firstStatesInBsccs.set(state);
  374. }
  375. first = false;
  376. }
  377. }
  378. storm::storage::BitVector statesNotInBsccs = ~statesInBsccs;
  379. STORM_LOG_DEBUG("Found " << statesInBsccs.getNumberOfSetBits() << " states in BSCCs.");
  380. // Prepare a vector holding the index within all states that are in BSCCs for every state.
  381. std::vector<uint_fast64_t> indexInStatesInBsccs;
  382. // Prepare a vector that maps the index within the set of all states in BSCCs to the index of the containing BSCC.
  383. std::vector<uint_fast64_t> stateToBsccIndexMap;
  384. if (!statesInBsccs.empty()) {
  385. firstStatesInBsccs = firstStatesInBsccs % statesInBsccs;
  386. // Then we construct an equation system that yields the steady state probabilities for all states in BSCCs.
  387. storm::storage::SparseMatrix<ValueType> bsccEquationSystem = probabilityMatrix.getSubmatrix(false, statesInBsccs, statesInBsccs, true);
  388. // Since in the fix point equation, we need to multiply the vector from the left, we convert this to a
  389. // multiplication from the right by transposing the system.
  390. bsccEquationSystem = bsccEquationSystem.transpose(false, true);
  391. // Create an auxiliary structure that makes it easy to look up the indices within the set of BSCC states.
  392. uint_fast64_t lastIndex = 0;
  393. uint_fast64_t currentNumberOfSetBits = 0;
  394. indexInStatesInBsccs.reserve(probabilityMatrix.getRowCount());
  395. for (auto index : statesInBsccs) {
  396. while (lastIndex <= index) {
  397. indexInStatesInBsccs.push_back(currentNumberOfSetBits);
  398. ++lastIndex;
  399. }
  400. ++currentNumberOfSetBits;
  401. }
  402. stateToBsccIndexMap.resize(statesInBsccs.getNumberOfSetBits());
  403. for (uint_fast64_t currentBsccIndex = 0; currentBsccIndex < bsccDecomposition.size(); ++currentBsccIndex) {
  404. storm::storage::StronglyConnectedComponent const& bscc = bsccDecomposition[currentBsccIndex];
  405. for (auto const& state : bscc) {
  406. stateToBsccIndexMap[indexInStatesInBsccs[state]] = currentBsccIndex;
  407. }
  408. }
  409. // Now build the final equation system matrix, the initial guess and the right-hand side in one go.
  410. std::vector<ValueType> bsccEquationSystemRightSide(bsccEquationSystem.getColumnCount(), zero);
  411. storm::storage::SparseMatrixBuilder<ValueType> builder;
  412. for (uint_fast64_t row = 0; row < bsccEquationSystem.getRowCount(); ++row) {
  413. // If the current row is the first one belonging to a BSCC, we substitute it by the constraint that the
  414. // values for states of this BSCC must sum to one. However, in order to have a non-zero value on the
  415. // diagonal, we add the constraint of the BSCC that produces a 1 on the diagonal.
  416. if (firstStatesInBsccs.get(row)) {
  417. uint_fast64_t requiredBscc = stateToBsccIndexMap[row];
  418. storm::storage::StronglyConnectedComponent const& bscc = bsccDecomposition[requiredBscc];
  419. for (auto const& state : bscc) {
  420. builder.addNextValue(row, indexInStatesInBsccs[state], one);
  421. }
  422. bsccEquationSystemRightSide[row] = one;
  423. } else {
  424. // Otherwise, we copy the row, and subtract 1 from the diagonal.
  425. for (auto& entry : bsccEquationSystem.getRow(row)) {
  426. if (entry.getColumn() == row) {
  427. builder.addNextValue(row, entry.getColumn(), entry.getValue() - one);
  428. } else {
  429. builder.addNextValue(row, entry.getColumn(), entry.getValue());
  430. }
  431. }
  432. }
  433. }
  434. // Create the initial guess for the LRAs. We take a uniform distribution over all states in a BSCC.
  435. std::vector<ValueType> bsccEquationSystemSolution(bsccEquationSystem.getColumnCount(), zero);
  436. for (uint_fast64_t bsccIndex = 0; bsccIndex < bsccDecomposition.size(); ++bsccIndex) {
  437. storm::storage::StronglyConnectedComponent const& bscc = bsccDecomposition[bsccIndex];
  438. for (auto const& state : bscc) {
  439. bsccEquationSystemSolution[indexInStatesInBsccs[state]] = one / bscc.size();
  440. }
  441. }
  442. bsccEquationSystem = builder.build();
  443. {
  444. std::unique_ptr<storm::solver::LinearEquationSolver<ValueType>> solver = linearEquationSolverFactory.create(std::move(bsccEquationSystem));
  445. solver->solveEquationSystem(bsccEquationSystemSolution, bsccEquationSystemRightSide);
  446. }
  447. // If exit rates were given, we need to 'fix' the results to also account for the timing behaviour.
  448. if (exitRateVector != nullptr) {
  449. std::vector<ValueType> bsccTotalValue(bsccDecomposition.size(), zero);
  450. for (auto stateIter = statesInBsccs.begin(); stateIter != statesInBsccs.end(); ++stateIter) {
  451. bsccTotalValue[stateToBsccIndexMap[indexInStatesInBsccs[*stateIter]]] += bsccEquationSystemSolution[indexInStatesInBsccs[*stateIter]] * (one / (*exitRateVector)[*stateIter]);
  452. }
  453. for (auto stateIter = statesInBsccs.begin(); stateIter != statesInBsccs.end(); ++stateIter) {
  454. bsccEquationSystemSolution[indexInStatesInBsccs[*stateIter]] = (bsccEquationSystemSolution[indexInStatesInBsccs[*stateIter]] * (one / (*exitRateVector)[*stateIter])) / bsccTotalValue[stateToBsccIndexMap[indexInStatesInBsccs[*stateIter]]];
  455. }
  456. }
  457. // Calculate LRA Value for each BSCC from steady state distribution in BSCCs.
  458. for (uint_fast64_t bsccIndex = 0; bsccIndex < bsccDecomposition.size(); ++bsccIndex) {
  459. storm::storage::StronglyConnectedComponent const& bscc = bsccDecomposition[bsccIndex];
  460. for (auto const& state : bscc) {
  461. if (psiStates.get(state)) {
  462. bsccLra[stateToBsccIndexMap[indexInStatesInBsccs[state]]] += bsccEquationSystemSolution[indexInStatesInBsccs[state]];
  463. }
  464. }
  465. }
  466. for (uint_fast64_t bsccIndex = 0; bsccIndex < bsccDecomposition.size(); ++bsccIndex) {
  467. STORM_LOG_DEBUG("Found LRA " << bsccLra[bsccIndex] << " for BSCC " << bsccIndex << ".");
  468. }
  469. } else {
  470. for (uint_fast64_t bsccIndex = 0; bsccIndex < bsccDecomposition.size(); ++bsccIndex) {
  471. storm::storage::StronglyConnectedComponent const& bscc = bsccDecomposition[bsccIndex];
  472. // At this point, all BSCCs are known to contain exactly one state, which is why we can set all values
  473. // directly (based on whether or not the contained state is a psi state).
  474. if (psiStates.get(*bscc.begin())) {
  475. bsccLra[bsccIndex] = storm::utility::one<ValueType>();
  476. }
  477. }
  478. for (uint_fast64_t bsccIndex = 0; bsccIndex < bsccDecomposition.size(); ++bsccIndex) {
  479. STORM_LOG_DEBUG("Found LRA " << bsccLra[bsccIndex] << " for BSCC " << bsccIndex << ".");
  480. }
  481. }
  482. std::vector<ValueType> rewardSolution;
  483. if (!statesNotInBsccs.empty()) {
  484. // Calculate LRA for states not in bsccs as expected reachability rewards.
  485. // Target states are states in bsccs, transition reward is the lra of the bscc for each transition into a
  486. // bscc and 0 otherwise. This corresponds to the sum of LRAs in BSCC weighted by the reachability probability
  487. // of the BSCC.
  488. std::vector<ValueType> rewardRightSide;
  489. rewardRightSide.reserve(statesNotInBsccs.getNumberOfSetBits());
  490. for (auto state : statesNotInBsccs) {
  491. ValueType reward = zero;
  492. for (auto entry : probabilityMatrix.getRow(state)) {
  493. if (statesInBsccs.get(entry.getColumn())) {
  494. reward += entry.getValue() * bsccLra[stateToBsccIndexMap[indexInStatesInBsccs[entry.getColumn()]]];
  495. }
  496. }
  497. rewardRightSide.push_back(reward);
  498. }
  499. storm::storage::SparseMatrix<ValueType> rewardEquationSystemMatrix = probabilityMatrix.getSubmatrix(false, statesNotInBsccs, statesNotInBsccs, true);
  500. rewardEquationSystemMatrix.convertToEquationSystem();
  501. rewardSolution = std::vector<ValueType>(rewardEquationSystemMatrix.getColumnCount(), one);
  502. {
  503. std::unique_ptr<storm::solver::LinearEquationSolver<ValueType>> solver = linearEquationSolverFactory.create(std::move(rewardEquationSystemMatrix));
  504. solver->solveEquationSystem(rewardSolution, rewardRightSide);
  505. }
  506. }
  507. // Fill the result vector.
  508. std::vector<ValueType> result(numberOfStates);
  509. auto rewardSolutionIter = rewardSolution.begin();
  510. for (uint_fast64_t bsccIndex = 0; bsccIndex < bsccDecomposition.size(); ++bsccIndex) {
  511. storm::storage::StronglyConnectedComponent const& bscc = bsccDecomposition[bsccIndex];
  512. for (auto const& state : bscc) {
  513. result[state] = bsccLra[bsccIndex];
  514. }
  515. }
  516. for (auto state : statesNotInBsccs) {
  517. STORM_LOG_ASSERT(rewardSolutionIter != rewardSolution.end(), "Too few elements in solution.");
  518. // Take the value from the reward computation. Since the n-th state not in any bscc is the n-th
  519. // entry in rewardSolution we can just take the next value from the iterator.
  520. result[state] = *rewardSolutionIter;
  521. ++rewardSolutionIter;
  522. }
  523. return result;
  524. }
  525. template <typename ValueType>
  526. std::vector<ValueType> SparseCtmcCslHelper<ValueType>::computeReachabilityTimes(storm::storage::SparseMatrix<ValueType> const& rateMatrix, storm::storage::SparseMatrix<ValueType> const& backwardTransitions, std::vector<ValueType> const& exitRateVector, storm::storage::BitVector const& initialStates, storm::storage::BitVector const& targetStates, bool qualitative, storm::solver::LinearEquationSolverFactory<ValueType> const& linearEquationSolverFactory) {
  527. // Compute expected time on CTMC by reduction to DTMC with rewards.
  528. storm::storage::SparseMatrix<ValueType> probabilityMatrix = computeProbabilityMatrix(rateMatrix, exitRateVector);
  529. // Initialize rewards.
  530. std::vector<ValueType> totalRewardVector;
  531. for (size_t i = 0; i < exitRateVector.size(); ++i) {
  532. if (targetStates[i] || storm::utility::isZero(exitRateVector[i])) {
  533. // Set reward for target states or states without outgoing transitions to 0.
  534. totalRewardVector.push_back(storm::utility::zero<ValueType>());
  535. } else {
  536. // Reward is (1 / exitRate).
  537. totalRewardVector.push_back(storm::utility::one<ValueType>() / exitRateVector[i]);
  538. }
  539. }
  540. return storm::modelchecker::helper::SparseDtmcPrctlHelper<ValueType>::computeReachabilityRewards(probabilityMatrix, backwardTransitions, totalRewardVector, targetStates, qualitative, linearEquationSolverFactory);
  541. }
  542. template class SparseCtmcCslHelper<double>;
  543. template std::vector<double> SparseCtmcCslHelper<double>::computeInstantaneousRewards(storm::storage::SparseMatrix<double> const& rateMatrix, std::vector<double> const& exitRateVector, storm::models::sparse::StandardRewardModel<double> const& rewardModel, double timeBound, storm::solver::LinearEquationSolverFactory<double> const& linearEquationSolverFactory);
  544. template std::vector<double> SparseCtmcCslHelper<double>::computeCumulativeRewards(storm::storage::SparseMatrix<double> const& rateMatrix, std::vector<double> const& exitRateVector, storm::models::sparse::StandardRewardModel<double> const& rewardModel, double timeBound, storm::solver::LinearEquationSolverFactory<double> const& linearEquationSolverFactory);
  545. template std::vector<double> SparseCtmcCslHelper<double>::computeReachabilityRewards(storm::storage::SparseMatrix<double> const& rateMatrix, storm::storage::SparseMatrix<double> const& backwardTransitions, std::vector<double> const& exitRateVector, storm::models::sparse::StandardRewardModel<double> const& rewardModel, storm::storage::BitVector const& targetStates, bool qualitative, storm::solver::LinearEquationSolverFactory<double> const& linearEquationSolverFactory);
  546. template std::vector<storm::RationalNumber> SparseCtmcCslHelper<storm::RationalNumber>::computeLongRunAverageProbabilities(storm::storage::SparseMatrix<storm::RationalNumber> const& probabilityMatrix, storm::storage::BitVector const& psiStates, std::vector<storm::RationalNumber> const* exitRateVector, bool qualitative, storm::solver::LinearEquationSolverFactory<storm::RationalNumber> const& linearEquationSolverFactory);
  547. template std::vector<storm::RationalFunction> SparseCtmcCslHelper<storm::RationalFunction>::computeLongRunAverageProbabilities(storm::storage::SparseMatrix<storm::RationalFunction> const& probabilityMatrix, storm::storage::BitVector const& psiStates, std::vector<storm::RationalFunction> const* exitRateVector, bool qualitative, storm::solver::LinearEquationSolverFactory<storm::RationalFunction> const& linearEquationSolverFactory);
  548. }
  549. }
  550. }