Browse Source

Added a pseudo model which can be constructed from only a matrix to look and behave like a model for use in Decomposition classes

Former-commit-id: f8fdc5a9b6
tempestpy_adaptions
PBerger 11 years ago
parent
commit
57b6208eee
  1. 2
      src/modelchecker/prctl/TopologicalValueIterationMdpPrctlModelChecker.h
  2. 108
      src/models/PseudoModel.cpp
  3. 90
      src/models/PseudoModel.h
  4. 105
      src/solver/TopologicalValueIterationNondeterministicLinearEquationSolver.cpp
  5. 3
      src/solver/TopologicalValueIterationNondeterministicLinearEquationSolver.h
  6. 1
      src/storage/Decomposition.h
  7. 33
      src/storage/StronglyConnectedComponentDecomposition.cpp
  8. 55
      src/storage/StronglyConnectedComponentDecomposition.h

2
src/modelchecker/prctl/TopologicalValueIterationMdpPrctlModelChecker.h

@ -68,7 +68,7 @@ private:
bool relative = s->getOptionByLongName("relative").getArgument(0).getValueAsBoolean();
// Now, we need to determine the SCCs of the MDP and a topological sort.
std::vector<std::vector<uint_fast64_t>> stronglyConnectedComponents = storm::utility::graph::performSccDecomposition(this->getModel(), stronglyConnectedComponents, stronglyConnectedComponentsDependencyGraph);
std::vector<std::vector<uint_fast64_t>> stronglyConnectedComponents = storm::utility::graph::performSccDecomposition(matrix, nondeterministicChoiceIndices, stronglyConnectedComponents, stronglyConnectedComponentsDependencyGraph);
storm::storage::SparseMatrix<T> stronglyConnectedComponentsDependencyGraph = this->getModel().extractSccDependencyGraph(stronglyConnectedComponents);
std::vector<uint_fast64_t> topologicalSort = storm::utility::graph::getTopologicalSort(stronglyConnectedComponentsDependencyGraph);

108
src/models/PseudoModel.cpp

@ -0,0 +1,108 @@
#include "src/models/PseudoModel.h"
#include "src/utility/constants.h"
#include "src/models/AbstractModel.h"
namespace storm {
namespace models {
template<typename ValueType>
ModelBasedPseudoModel<ValueType>::ModelBasedPseudoModel(storm::models::AbstractModel<ValueType> const& model) : _model(model) {
// Intentionally left empty.
}
template<typename ValueType>
NonDeterministicMatrixBasedPseudoModel<ValueType>::NonDeterministicMatrixBasedPseudoModel(storm::storage::SparseMatrix<ValueType> const& matrix, std::vector<uint_fast64_t> const& nondeterministicChoiceIndices) : _matrix(matrix), _nondeterministicChoiceIndices(nondeterministicChoiceIndices) {
// Intentionally left empty.
}
template<typename ValueType>
DeterministicMatrixBasedPseudoModel<ValueType>::DeterministicMatrixBasedPseudoModel(storm::storage::SparseMatrix<ValueType> const& matrix) : _matrix(matrix) {
// Intentionally left empty.
}
template<typename ValueType>
typename storm::storage::SparseMatrix<ValueType>::const_rows
ModelBasedPseudoModel<ValueType>::getRows(uint_fast64_t state) const {
return this->_model.getRows(state);
}
template<typename ValueType>
typename storm::storage::SparseMatrix<ValueType>::const_rows
NonDeterministicMatrixBasedPseudoModel<ValueType>::getRows(uint_fast64_t state) const {
return this->_matrix.getRows(this->_nondeterministicChoiceIndices[state], this->_nondeterministicChoiceIndices[state + 1] - 1);
}
template<typename ValueType>
typename storm::storage::SparseMatrix<ValueType>::const_rows
DeterministicMatrixBasedPseudoModel<ValueType>::getRows(uint_fast64_t state) const {
return this->_matrix.getRows(state, state);
}
template<typename ValueType>
uint_fast64_t
ModelBasedPseudoModel<ValueType>::getNumberOfStates() const {
return this->_model.getNumberOfStates();
}
template<typename ValueType>
uint_fast64_t
NonDeterministicMatrixBasedPseudoModel<ValueType>::getNumberOfStates() const {
return this->_matrix.getColumnCount();
}
template<typename ValueType>
uint_fast64_t
DeterministicMatrixBasedPseudoModel<ValueType>::getNumberOfStates() const {
return this->_matrix.getColumnCount();
}
template<typename ValueType>
storm::storage::SparseMatrix<ValueType>
AbstractPseudoModel<ValueType>::extractPartitionDependencyGraph(storm::storage::Decomposition<storm::storage::StateBlock> const& decomposition) const {
uint_fast64_t numberOfStates = decomposition.size();
// First, we need to create a mapping of states to their SCC index, to ease the computation of dependency transitions later.
std::vector<uint_fast64_t> stateToBlockMap(this->getNumberOfStates());
for (uint_fast64_t i = 0; i < decomposition.size(); ++i) {
for (auto state : decomposition[i]) {
stateToBlockMap[state] = i;
}
}
// The resulting sparse matrix will have as many rows/columns as there are blocks in the partition.
storm::storage::SparseMatrixBuilder<ValueType> dependencyGraphBuilder(numberOfStates, numberOfStates);
for (uint_fast64_t currentBlockIndex = 0; currentBlockIndex < decomposition.size(); ++currentBlockIndex) {
// Get the next block.
typename storm::storage::StateBlock const& block = decomposition[currentBlockIndex];
// Now, we determine the blocks which are reachable (in one step) from the current block.
boost::container::flat_set<uint_fast64_t> allTargetBlocks;
for (auto state : block) {
for (auto const& transitionEntry : this->getRows(state)) {
uint_fast64_t targetBlock = stateToBlockMap[transitionEntry.first];
// We only need to consider transitions that are actually leaving the SCC.
if (targetBlock != currentBlockIndex) {
allTargetBlocks.insert(targetBlock);
}
}
}
// Now we can just enumerate all the target SCCs and insert the corresponding transitions.
for (auto targetBlock : allTargetBlocks) {
dependencyGraphBuilder.addNextValue(currentBlockIndex, targetBlock, storm::utility::constantOne<ValueType>());
}
}
return dependencyGraphBuilder.build();
}
template class ModelBasedPseudoModel<double>;
template class NonDeterministicMatrixBasedPseudoModel<double>;
template class DeterministicMatrixBasedPseudoModel<double>;
template class ModelBasedPseudoModel<int>;
template class NonDeterministicMatrixBasedPseudoModel<int>;
template class DeterministicMatrixBasedPseudoModel<int>;
} // namespace models
} // namespace storm

90
src/models/PseudoModel.h

@ -0,0 +1,90 @@
#ifndef STORM_MODELS_PSEUDOMODEL_H_
#define STORM_MODELS_PSEUDOMODEL_H_
#include <cstdint>
#include "src/storage/SparseMatrix.h"
#include "src/storage/Decomposition.h"
namespace storm {
namespace models {
// Forward declare the abstract model class.
template <typename ValueType> class AbstractModel;
/*!
* This classes encapsulate the model/transitionmatrix on which the SCC decomposition is performed.
* The Abstract Base class is specialized by the two possible representations:
* - For a model the implementation ModelBasedPseudoModel hands the call to getRows() through to the model
* - For a matrix of a nondeterministic model the implementation NonDeterministicMatrixBasedPseudoModel emulates the call
* on the matrix itself like the model function would
* - For a matrix of a deterministic model the implementation DeterministicMatrixBasedPseudoModel emulates the call
* on the matrix itself like the model function would
*/
template <typename ValueType>
class AbstractPseudoModel {
public:
AbstractPseudoModel() {}
virtual ~AbstractPseudoModel() {}
virtual typename storm::storage::SparseMatrix<ValueType>::const_rows getRows(uint_fast64_t state) const = 0;
/*!
* Calculates the number of states in the represented system.
* @return The Number of States in the underlying model/transition matrix
*/
virtual uint_fast64_t getNumberOfStates() const = 0;
/*!
* Extracts the dependency graph from a (pseudo-) model according to the given partition.
*
* @param decomposition A decomposition containing the blocks of the partition of the system.
* @return A sparse matrix with bool entries that represents the dependency graph of the blocks of the partition.
*/
virtual storm::storage::SparseMatrix<ValueType> extractPartitionDependencyGraph(storm::storage::Decomposition<storm::storage::StateBlock> const& decomposition) const;
};
template <typename ValueType>
class ModelBasedPseudoModel : public AbstractPseudoModel<ValueType> {
public:
/*!
* Creates an encapsulation for the SCC decomposition based on a model
* @param model The Model on which the decomposition is to be performed
*/
ModelBasedPseudoModel(storm::models::AbstractModel<ValueType> const& model);
virtual typename storm::storage::SparseMatrix<ValueType>::const_rows getRows(uint_fast64_t state) const override;
virtual uint_fast64_t getNumberOfStates() const override;
private:
storm::models::AbstractModel<ValueType> const& _model;
};
template <typename ValueType>
class NonDeterministicMatrixBasedPseudoModel : public AbstractPseudoModel<ValueType> {
public:
/*!
* Creates an encapsulation for the SCC decomposition based on a matrix
* @param matrix The Matrix on which the decomposition is to be performed
*/
NonDeterministicMatrixBasedPseudoModel(storm::storage::SparseMatrix<ValueType> const& matrix, std::vector<uint_fast64_t> const& nondeterministicChoiceIndices);
virtual typename storm::storage::SparseMatrix<ValueType>::const_rows getRows(uint_fast64_t state) const override;
virtual uint_fast64_t getNumberOfStates() const override;
private:
storm::storage::SparseMatrix<ValueType> const& _matrix;
std::vector<uint_fast64_t> const& _nondeterministicChoiceIndices;
};
template <typename ValueType>
class DeterministicMatrixBasedPseudoModel : public AbstractPseudoModel<ValueType> {
public:
/*!
* Creates an encapsulation for the SCC decomposition based on a matrix
* @param matrix The Matrix on which the decomposition is to be performed
*/
DeterministicMatrixBasedPseudoModel(storm::storage::SparseMatrix<ValueType> const& matrix);
virtual typename storm::storage::SparseMatrix<ValueType>::const_rows getRows(uint_fast64_t state) const override;
virtual uint_fast64_t getNumberOfStates() const override;
private:
storm::storage::SparseMatrix<ValueType> const& _matrix;
};
}
}
#endif // STORM_MODELS_PSEUDOMODEL_H_

105
src/solver/TopologicalValueIterationNondeterministicLinearEquationSolver.cpp

@ -1,34 +1,49 @@
#include "src/solver/TopologicalValueIterationNativeNondeterministicLinearEquationSolver.h"
#include "src/solver/TopologicalValueIterationNondeterministicLinearEquationSolver.h"
#include <utility>
#include "src/settings/Settings.h"
#include "src/utility/vector.h"
#include "src/utility/graph.h"
#include "src/models/PseudoModel.h"
#include "src/storage/StronglyConnectedComponentDecomposition.h"
namespace storm {
namespace solver {
template<typename ValueType>
TopologicalValueIterationNondeterministicLinearEquationSolver<ValueType>::TopologicalValueIterationNondeterministicLinearEquationSolver() {
// // Intentionally left empty.
// Get the settings object to customize solving.
storm::settings::Settings* settings = storm::settings::Settings::getInstance();
// Get appropriate settings.
maximalNumberOfIterations = settings->getOptionByLongName("maxiter").getArgument(0).getValueAsUnsignedInteger();
precision = settings->getOptionByLongName("precision").getArgument(0).getValueAsDouble();
relative = !settings->isSet("absolute");
}
template<typename ValueType>
TopologicalValueIterationNondeterministicLinearEquationSolver<ValueType>::TopologicalValueIterationNondeterministicLinearEquationSolver(double precision, uint_fast64_t maximalNumberOfIterations, bool relative) : NativeNondeterministicLinearEquationSolver<ValueType>(precision, maximalNumberOfIterations, relative) {
TopologicalValueIterationNondeterministicLinearEquationSolver<ValueType>::TopologicalValueIterationNondeterministicLinearEquationSolver(double precision, uint_fast64_t maximalNumberOfIterations, bool relative) : NativeNondeterministicLinearEquationSolver(precision, maximalNumberOfIterations, relative) {
// Intentionally left empty.
}
template<typename ValueType>
TopologicalValueIterationNondeterministicLinearEquationSolver<ValueType>* TopologicalValueIterationNondeterministicLinearEquationSolver<ValueType>::clone() const {
return new NativeNondeterministicLinearEquationSolver<ValueType>(*this);
NondeterministicLinearEquationSolver<ValueType>* TopologicalValueIterationNondeterministicLinearEquationSolver<ValueType>::clone() const {
return new TopologicalValueIterationNondeterministicLinearEquationSolver<ValueType>(*this);
}
template<typename ValueType>
void TopologicalValueIterationNondeterministicLinearEquationSolver<ValueType>::solveEquationSystem(bool minimize, storm::storage::SparseMatrix<ValueType> const& A, std::vector<ValueType>& x, std::vector<ValueType> const& b, std::vector<uint_fast64_t> const& nondeterministicChoiceIndices, std::vector<ValueType>* multiplyResult, std::vector<ValueType>* newX) const {
// Now, we need to determine the SCCs of the MDP and a topological sort.
std::vector<std::vector<uint_fast64_t>> stronglyConnectedComponents = storm::utility::graph::performSccDecomposition(this->getModel(), stronglyConnectedComponents, stronglyConnectedComponentsDependencyGraph);
storm::storage::SparseMatrix<T> stronglyConnectedComponentsDependencyGraph = this->getModel().extractSccDependencyGraph(stronglyConnectedComponents);
//std::vector<std::vector<uint_fast64_t>> stronglyConnectedComponents = storm::utility::graph::performSccDecomposition(this->getModel(), stronglyConnectedComponents, stronglyConnectedComponentsDependencyGraph);
//storm::storage::SparseMatrix<T> stronglyConnectedComponentsDependencyGraph = this->getModel().extractSccDependencyGraph(stronglyConnectedComponents);
storm::models::NonDeterministicMatrixBasedPseudoModel<ValueType> pseudoModel(A, nondeterministicChoiceIndices);
storm::storage::StronglyConnectedComponentDecomposition<ValueType> sccDecomposition(*static_cast<storm::models::AbstractPseudoModel<ValueType>*>(&pseudoModel), false, false);
storm::storage::SparseMatrix<ValueType> stronglyConnectedComponentsDependencyGraph = pseudoModel.extractPartitionDependencyGraph(sccDecomposition);
std::vector<uint_fast64_t> topologicalSort = storm::utility::graph::getTopologicalSort(stronglyConnectedComponentsDependencyGraph);
// Set up the environment for the power method.
@ -52,14 +67,14 @@ namespace storm {
// Iterate over all SCCs of the MDP as specified by the topological sort. This guarantees that an SCC is only
// solved after all SCCs it depends on have been solved.
for (auto sccIndexIt = topologicalSort.begin(); sccIndexIt != topologicalSort.end() && converged; ++sccIndexIt) {
std::vector<uint_fast64_t> const& scc = stronglyConnectedComponents[*sccIndexIt];
std::vector<uint_fast64_t> const& scc = sccDecomposition[*sccIndexIt];
// For the current SCC, we need to perform value iteration until convergence.
localIterations = 0;
converged = false;
while (!converged && localIterations < maxIterations) {
while (!converged && localIterations < maximalNumberOfIterations) {
// Compute x' = A*x + b.
matrix.multiplyWithVector(scc, nondeterministicChoiceIndices, *currentX, multiplyResult);
A.multiplyWithVector(scc, nondeterministicChoiceIndices, *currentX, multiplyResult);
storm::utility::addVectors(scc, nondeterministicChoiceIndices, multiplyResult, b);
/*
@ -119,76 +134,6 @@ namespace storm {
else {
LOG4CPLUS_WARN(logger, "Iterative solver did not converged after " << currentMaxLocalIterations << " iterations.");
}
/*
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
*/
// Set up the environment for the power method. If scratch memory was not provided, we need to create it.
bool multiplyResultMemoryProvided = true;
if (multiplyResult == nullptr) {
multiplyResult = new std::vector<ValueType>(A.getRowCount());
multiplyResultMemoryProvided = false;
}
std::vector<ValueType>* currentX = &x;
bool xMemoryProvided = true;
if (newX == nullptr) {
newX = new std::vector<ValueType>(x.size());
xMemoryProvided = false;
}
std::vector<ValueType>* swap = nullptr;
uint_fast64_t iterations = 0;
bool converged = false;
// Keep track of which of the vectors for x is the auxiliary copy.
std::vector<ValueType>* copyX = newX;
// Proceed with the iterations as long as the method did not converge or reach the
// user-specified maximum number of iterations.
while (!converged && iterations < maximalNumberOfIterations) {
// Compute x' = A*x + b.
A.multiplyWithVector(*currentX, *multiplyResult);
storm::utility::vector::addVectorsInPlace(*multiplyResult, b);
// Reduce the vector x' by applying min/max for all non-deterministic choices as given by the topmost
// element of the min/max operator stack.
if (minimize) {
storm::utility::vector::reduceVectorMin(*multiplyResult, *newX, nondeterministicChoiceIndices);
} else {
storm::utility::vector::reduceVectorMax(*multiplyResult, *newX, nondeterministicChoiceIndices);
}
// Determine whether the method converged.
converged = storm::utility::vector::equalModuloPrecision(*currentX, *newX, precision, relative);
// Update environment variables.
std::swap(currentX, newX);
++iterations;
}
// Check if the solver converged and issue a warning otherwise.
if (converged) {
LOG4CPLUS_INFO(logger, "Iterative solver converged after " << iterations << " iterations.");
} else {
LOG4CPLUS_WARN(logger, "Iterative solver did not converge after " << iterations << " iterations.");
}
// If we performed an odd number of iterations, we need to swap the x and currentX, because the newest result
// is currently stored in currentX, but x is the output vector.
if (currentX == copyX) {
std::swap(x, *currentX);
}
if (!xMemoryProvided) {
delete copyX;
}
if (!multiplyResultMemoryProvided) {
delete multiplyResult;
}
}
// Explicitly instantiate the solver.

3
src/solver/TopologicalValueIterationNativeNondeterministicLinearEquationSolver.h → src/solver/TopologicalValueIterationNondeterministicLinearEquationSolver.h

@ -1,7 +1,6 @@
#ifndef STORM_SOLVER_TOPOLOGICALVALUEITERATIONNONDETERMINISTICLINEAREQUATIONSOLVER_H_
#define STORM_SOLVER_TOPOLOGICALVALUEITERATIONNONDETERMINISTICLINEAREQUATIONSOLVER_H_
#include "src/solver/NondeterministicLinearEquationSolver.h"
#include "src/solver/NativeNondeterministicLinearEquationSolver.h"
namespace storm {
@ -11,7 +10,7 @@ namespace storm {
* A class that uses SCC Decompositions to solve a linear equation system
*/
template<class ValueType>
class TopologicalValueIterationNondeterministicLinearEquationSolver : public NondeterministicLinearEquationSolver<ValueType> {
class TopologicalValueIterationNondeterministicLinearEquationSolver : public NativeNondeterministicLinearEquationSolver<ValueType> {
public:
/*!
* Constructs a nondeterministic linear equation solver with parameters being set according to the settings

1
src/storage/Decomposition.h

@ -133,7 +133,6 @@ namespace storm {
// Declare the streaming operator as a friend function to enable output of decompositions.
template<typename BlockTimePrime>
friend std::ostream& operator<<(std::ostream& out, Decomposition<BlockTimePrime> const& decomposition);
protected:
// The blocks of the decomposition.
std::vector<Block> blocks;

33
src/storage/StronglyConnectedComponentDecomposition.cpp

@ -1,5 +1,6 @@
#include "src/storage/StronglyConnectedComponentDecomposition.h"
#include "src/models/AbstractModel.h"
#include "src/models/PseudoModel.h"
namespace storm {
namespace storage {
@ -10,18 +11,36 @@ namespace storm {
template <typename ValueType>
StronglyConnectedComponentDecomposition<ValueType>::StronglyConnectedComponentDecomposition(storm::models::AbstractModel<ValueType> const& model, bool dropNaiveSccs, bool onlyBottomSccs) : Decomposition() {
performSccDecomposition(model, dropNaiveSccs, onlyBottomSccs);
performSccDecomposition(storm::models::ModelBasedPseudoModel<ValueType>(model), dropNaiveSccs, onlyBottomSccs);
}
template <typename ValueType>
StronglyConnectedComponentDecomposition<ValueType>::StronglyConnectedComponentDecomposition(storm::models::AbstractPseudoModel<ValueType> const& pseudoModel, bool dropNaiveSccs, bool onlyBottomSccs) : Decomposition() {
performSccDecomposition(pseudoModel, dropNaiveSccs, onlyBottomSccs);
}
template <typename ValueType>
StronglyConnectedComponentDecomposition<ValueType>::StronglyConnectedComponentDecomposition(storm::models::AbstractModel<ValueType> const& model, StateBlock const& block, bool dropNaiveSccs, bool onlyBottomSccs) {
storm::storage::BitVector subsystem(model.getNumberOfStates(), block.begin(), block.end());
performSccDecomposition(model, subsystem, dropNaiveSccs, onlyBottomSccs);
storm::models::ModelBasedPseudoModel<ValueType> encapsulatedModel(model);
storm::storage::BitVector subsystem(encapsulatedModel.getNumberOfStates(), block.begin(), block.end());
performSccDecomposition(*static_cast<storm::models::AbstractPseudoModel<ValueType>*>(&encapsulatedModel), subsystem, dropNaiveSccs, onlyBottomSccs);
}
template <typename ValueType>
StronglyConnectedComponentDecomposition<ValueType>::StronglyConnectedComponentDecomposition(storm::models::AbstractPseudoModel<ValueType> const& pseudoModel, StateBlock const& block, bool dropNaiveSccs, bool onlyBottomSccs) {
storm::storage::BitVector subsystem(pseudoModel.getNumberOfStates(), block.begin(), block.end());
performSccDecomposition(pseudoModel, subsystem, dropNaiveSccs, onlyBottomSccs);
}
template <typename ValueType>
StronglyConnectedComponentDecomposition<ValueType>::StronglyConnectedComponentDecomposition(storm::models::AbstractModel<ValueType> const& model, storm::storage::BitVector const& subsystem, bool dropNaiveSccs, bool onlyBottomSccs) {
performSccDecomposition(model, subsystem, dropNaiveSccs, onlyBottomSccs);
storm::models::ModelBasedPseudoModel<ValueType> encapsulatedModel(model);
performSccDecomposition(*static_cast<storm::models::AbstractPseudoModel<ValueType>*>(&encapsulatedModel), subsystem, dropNaiveSccs, onlyBottomSccs);
}
template <typename ValueType>
StronglyConnectedComponentDecomposition<ValueType>::StronglyConnectedComponentDecomposition(storm::models::AbstractPseudoModel<ValueType> const& pseudoModel, storm::storage::BitVector const& subsystem, bool dropNaiveSccs, bool onlyBottomSccs) {
performSccDecomposition(pseudoModel, subsystem, dropNaiveSccs, onlyBottomSccs);
}
template <typename ValueType>
@ -47,7 +66,7 @@ namespace storm {
}
template <typename ValueType>
void StronglyConnectedComponentDecomposition<ValueType>::performSccDecomposition(storm::models::AbstractModel<ValueType> const& model, storm::storage::BitVector const& subsystem, bool dropNaiveSccs, bool onlyBottomSccs) {
void StronglyConnectedComponentDecomposition<ValueType>::performSccDecomposition(storm::models::AbstractPseudoModel<ValueType> const& model, storm::storage::BitVector const& subsystem, bool dropNaiveSccs, bool onlyBottomSccs) {
// Set up the environment of Tarjan's algorithm.
uint_fast64_t numberOfStates = model.getNumberOfStates();
std::vector<uint_fast64_t> tarjanStack;
@ -68,7 +87,7 @@ namespace storm {
template <typename ValueType>
void StronglyConnectedComponentDecomposition<ValueType>::performSccDecomposition(storm::models::AbstractModel<ValueType> const& model, bool dropNaiveSccs, bool onlyBottomSccs) {
void StronglyConnectedComponentDecomposition<ValueType>::performSccDecomposition(storm::models::AbstractPseudoModel<ValueType> const& model, bool dropNaiveSccs, bool onlyBottomSccs) {
// Prepare a block that contains all states for a call to the other overload of this function.
storm::storage::BitVector fullSystem(model.getNumberOfStates(), true);
@ -77,7 +96,7 @@ namespace storm {
}
template <typename ValueType>
void StronglyConnectedComponentDecomposition<ValueType>::performSccDecompositionHelper(storm::models::AbstractModel<ValueType> const& model, uint_fast64_t startState, storm::storage::BitVector const& subsystem, uint_fast64_t& currentIndex, std::vector<uint_fast64_t>& stateIndices, std::vector<uint_fast64_t>& lowlinks, std::vector<uint_fast64_t>& tarjanStack, storm::storage::BitVector& tarjanStackStates, storm::storage::BitVector& visitedStates, bool dropNaiveSccs, bool onlyBottomSccs) {
void StronglyConnectedComponentDecomposition<ValueType>::performSccDecompositionHelper(storm::models::AbstractPseudoModel<ValueType> const& model, uint_fast64_t startState, storm::storage::BitVector const& subsystem, uint_fast64_t& currentIndex, std::vector<uint_fast64_t>& stateIndices, std::vector<uint_fast64_t>& lowlinks, std::vector<uint_fast64_t>& tarjanStack, storm::storage::BitVector& tarjanStackStates, storm::storage::BitVector& visitedStates, bool dropNaiveSccs, bool onlyBottomSccs) {
// Create the stacks needed for turning the recursive formulation of Tarjan's algorithm into an iterative
// version. In particular, we keep one stack for states and one stack for the iterators. The last one is not
// strictly needed, but reduces iteration work when all successors of a particular state are considered.

55
src/storage/StronglyConnectedComponentDecomposition.h

@ -3,11 +3,15 @@
#include "src/storage/Decomposition.h"
#include "src/storage/BitVector.h"
#include "src/storage/SparseMatrix.h"
namespace storm {
namespace models {
// Forward declare the abstract model class.
template <typename ValueType> class AbstractModel;
// Forward declare the abstract pseudo-model class.
template <typename ValueType> class AbstractPseudoModel;
}
namespace storage {
@ -34,6 +38,17 @@ namespace storm {
*/
StronglyConnectedComponentDecomposition(storm::models::AbstractModel<ValueType> const& model, bool dropNaiveSccs = false, bool onlyBottomSccs = false);
/*
* Creates an SCC decomposition of the given encapsulated model.
*
* @param pseudoModel The encapsulated model to decompose into SCCs.
* @param dropNaiveSccs A flag that indicates whether trivial SCCs (i.e. SCCs consisting of just one state
* without a self-loop) are to be kept in the decomposition.
* @param onlyBottomSccs If set to true, only bottom SCCs, i.e. SCCs in which all states have no way of
* leaving the SCC), are kept.
*/
StronglyConnectedComponentDecomposition(storm::models::AbstractPseudoModel<ValueType> const& pseudoModel, bool dropNaiveSccs = false, bool onlyBottomSccs = false);
/*
* Creates an SCC decomposition of the given block in the given model.
*
@ -46,6 +61,18 @@ namespace storm {
*/
StronglyConnectedComponentDecomposition(storm::models::AbstractModel<ValueType> const& model, StateBlock const& block, bool dropNaiveSccs = false, bool onlyBottomSccs = false);
/*
* Creates an SCC decomposition of the given block in the given encapsulated model.
*
* @param pseudoModel The encapsulated model whose block to decompose.
* @param block The block to decompose into SCCs.
* @param dropNaiveSccs A flag that indicates whether trivial SCCs (i.e. SCCs consisting of just one state
* without a self-loop) are to be kept in the decomposition.
* @param onlyBottomSccs If set to true, only bottom SCCs, i.e. SCCs in which all states have no way of
* leaving the SCC), are kept.
*/
StronglyConnectedComponentDecomposition(storm::models::AbstractPseudoModel<ValueType> const& pseudoModel, StateBlock const& block, bool dropNaiveSccs = false, bool onlyBottomSccs = false);
/*
* Creates an SCC decomposition of the given subsystem in the given model.
*
@ -58,35 +85,49 @@ namespace storm {
*/
StronglyConnectedComponentDecomposition(storm::models::AbstractModel<ValueType> const& model, storm::storage::BitVector const& subsystem, bool dropNaiveSccs = false, bool onlyBottomSccs = false);
/*
* Creates an SCC decomposition of the given subsystem in the given encapsulated model.
*
* @param pseudoModel The encapsulated model that contains the block.
* @param subsystem A bit vector indicating which subsystem to consider for the decomposition into SCCs.
* @param dropNaiveSccs A flag that indicates whether trivial SCCs (i.e. SCCs consisting of just one state
* without a self-loop) are to be kept in the decomposition.
* @param onlyBottomSccs If set to true, only bottom SCCs, i.e. SCCs in which all states have no way of
* leaving the SCC), are kept.
*/
StronglyConnectedComponentDecomposition(storm::models::AbstractPseudoModel<ValueType> const& pseudoModel, storm::storage::BitVector const& subsystem, bool dropNaiveSccs = false, bool onlyBottomSccs = false);
/*!
* Creates an SCC decomposition by copying the given SCC decomposition.
*
* @oaram other The SCC decomposition to copy.
* @param other The SCC decomposition to copy.
*/
StronglyConnectedComponentDecomposition(StronglyConnectedComponentDecomposition const& other);
/*!
* Assigns the contents of the given SCC decomposition to the current one by copying its contents.
*
* @oaram other The SCC decomposition from which to copy-assign.
* @param other The SCC decomposition from which to copy-assign.
*/
StronglyConnectedComponentDecomposition& operator=(StronglyConnectedComponentDecomposition const& other);
/*!
* Creates an SCC decomposition by moving the given SCC decomposition.
*
* @oaram other The SCC decomposition to move.
* @param other The SCC decomposition to move.
*/
StronglyConnectedComponentDecomposition(StronglyConnectedComponentDecomposition&& other);
/*!
* Assigns the contents of the given SCC decomposition to the current one by moving its contents.
*
* @oaram other The SCC decomposition from which to copy-assign.
* @param other The SCC decomposition from which to copy-assign.
*/
StronglyConnectedComponentDecomposition& operator=(StronglyConnectedComponentDecomposition&& other);
private:
/*!
* Performs the SCC decomposition of the given model. As a side-effect this fills the vector of
* blocks of the decomposition.
@ -97,7 +138,7 @@ namespace storm {
* @param onlyBottomSccs If set to true, only bottom SCCs, i.e. SCCs in which all states have no way of
* leaving the SCC), are kept.
*/
void performSccDecomposition(storm::models::AbstractModel<ValueType> const& model, bool dropNaiveSccs, bool onlyBottomSccs);
void performSccDecomposition(storm::models::AbstractPseudoModel<ValueType> const& model, bool dropNaiveSccs, bool onlyBottomSccs);
/*
* Performs the SCC decomposition of the given block in the given model. As a side-effect this fills
@ -110,7 +151,7 @@ namespace storm {
* @param onlyBottomSccs If set to true, only bottom SCCs, i.e. SCCs in which all states have no way of
* leaving the SCC), are kept.
*/
void performSccDecomposition(storm::models::AbstractModel<ValueType> const& model, storm::storage::BitVector const& subsystem, bool dropNaiveSccs, bool onlyBottomSccs);
void performSccDecomposition(storm::models::AbstractPseudoModel<ValueType> const& model, storm::storage::BitVector const& subsystem, bool dropNaiveSccs, bool onlyBottomSccs);
/*!
* A helper function that performs the SCC decomposition given all auxiliary data structures. As a
@ -130,7 +171,7 @@ namespace storm {
* @param onlyBottomSccs If set to true, only bottom SCCs, i.e. SCCs in which all states have no way of
* leaving the SCC), are kept.
*/
void performSccDecompositionHelper(storm::models::AbstractModel<ValueType> const& model, uint_fast64_t startState, storm::storage::BitVector const& subsystem, uint_fast64_t& currentIndex, std::vector<uint_fast64_t>& stateIndices, std::vector<uint_fast64_t>& lowlinks, std::vector<uint_fast64_t>& tarjanStack, storm::storage::BitVector& tarjanStackStates, storm::storage::BitVector& visitedStates, bool dropNaiveSccs, bool onlyBottomSccs);
void performSccDecompositionHelper(storm::models::AbstractPseudoModel<ValueType> const& model, uint_fast64_t startState, storm::storage::BitVector const& subsystem, uint_fast64_t& currentIndex, std::vector<uint_fast64_t>& stateIndices, std::vector<uint_fast64_t>& lowlinks, std::vector<uint_fast64_t>& tarjanStack, storm::storage::BitVector& tarjanStackStates, storm::storage::BitVector& visitedStates, bool dropNaiveSccs, bool onlyBottomSccs);
};
}
}

Loading…
Cancel
Save