Browse Source

Moved creation of SCC-dependency graph into abstract model class. Added functionality to sparse matrix class to not give the number of nonzeros upfront, but to to grow on demand.

tempestpy_adaptions
dehnert 12 years ago
parent
commit
69395face2
  1. 41
      src/models/AbstractDeterministicModel.h
  2. 10
      src/models/AbstractModel.h
  3. 44
      src/models/AbstractNondeterministicModel.h
  4. 32
      src/models/GraphTransitions.h
  5. 79
      src/storage/SparseMatrix.h
  6. 114
      src/utility/GraphAnalyzer.h

41
src/models/AbstractDeterministicModel.h

@ -45,6 +45,47 @@ class AbstractDeterministicModel: public AbstractModel<T> {
AbstractDeterministicModel(AbstractDeterministicModel const& other) : AbstractModel<T>(other) {
// Intentionally left empty.
}
/*!
* Extracts the SCC dependency graph from the model according to the given SCC decomposition.
*
* @param stronglyConnectedComponents A vector containing the SCCs of the system.
* @param stateToSccMap A mapping from state indices to
*/
virtual storm::storage::SparseMatrix<bool> extractSccDependencyGraph(std::vector<std::vector<uint_fast64_t>> const& stronglyConnectedComponents, std::map<uint_fast64_t, uint_fast64_t> const& stateToSccMap) {
// The resulting sparse matrix will have as many rows/columns as there are SCCs.
uint_fast64_t numberOfStates = stronglyConnectedComponents.size();
storm::storage::SparseMatrix<bool> sccDependencyGraph(numberOfStates);
sccDependencyGraph.initialize();
for (uint_fast64_t currentSccIndex = 0; currentSccIndex < stronglyConnectedComponents.size(); ++currentSccIndex) {
// Get the actual SCC.
std::vector<uint_fast64_t> const& scc = stronglyConnectedComponents[currentSccIndex];
// Now, we determine the SCCs which are reachable (in one step) from the current SCC.
std::set<uint_fast64_t> allTargetSccs;
for (auto state : scc) {
for (typename storm::storage::SparseMatrix<T>::ConstIndexIterator succIt = this->getTransitionMatrix()->constColumnIteratorBegin(state), succIte = this->getTransitionMatrix()->constColumnIteratorEnd(state); succIt != succIte; ++succIt) {
uint_fast64_t targetScc = stateToSccMap.find(*succIt)->second;
// We only need to consider transitions that are actually leaving the SCC.
if (targetScc != currentSccIndex) {
allTargetSccs.insert(targetScc);
}
}
}
// Now we can just enumerate all the target SCCs and insert the corresponding transitions.
for (auto targetScc : allTargetSccs) {
sccDependencyGraph.insertNextValue(currentSccIndex, targetScc, true);
}
}
// Finalize the matrix.
sccDependencyGraph.finalize(true);
return sccDependencyGraph;
}
};
} // namespace models

10
src/models/AbstractModel.h

@ -85,6 +85,14 @@ class AbstractModel: public std::enable_shared_from_this<AbstractModel<T>> {
*/
virtual ModelType getType() const = 0;
/*!
* Extracts the SCC dependency graph from the model according to the given SCC decomposition.
*
* @param stronglyConnectedComponents A vector containing the SCCs of the system.
* @param stateToSccMap A mapping from state indices to
*/
virtual storm::storage::SparseMatrix<bool> extractSccDependencyGraph(std::vector<std::vector<uint_fast64_t>> const& stronglyConnectedComponents, std::map<uint_fast64_t, uint_fast64_t> const& stateToSccMap) = 0;
/*!
* Returns the state space size of the model.
* @return The size of the state space of the model.
@ -151,7 +159,7 @@ class AbstractModel: public std::enable_shared_from_this<AbstractModel<T>> {
* Returns the set of states with which the given state is labeled.
* @return The set of states with which the given state is labeled.
*/
std::set<std::string> const getPropositionsForState(uint_fast64_t const &state) const {
std::set<std::string> const getPropositionsForState(uint_fast64_t const& state) const {
return stateLabeling->getPropositionsForState(state);
}

44
src/models/AbstractNondeterministicModel.h

@ -51,7 +51,6 @@ class AbstractNondeterministicModel: public AbstractModel<T> {
// Intentionally left empty.
}
/*!
* Returns the number of choices for all states of the MDP.
* @return The number of choices for all states of the MDP.
@ -60,6 +59,49 @@ class AbstractNondeterministicModel: public AbstractModel<T> {
return this->getTransitionMatrix()->getRowCount();
}
/*!
* Extracts the SCC dependency graph from the model according to the given SCC decomposition.
*
* @param stronglyConnectedComponents A vector containing the SCCs of the system.
* @param stateToSccMap A mapping from state indices to
*/
virtual storm::storage::SparseMatrix<bool> extractSccDependencyGraph(std::vector<std::vector<uint_fast64_t>> const& stronglyConnectedComponents, std::map<uint_fast64_t, uint_fast64_t> const& stateToSccMap) {
// The resulting sparse matrix will have as many rows/columns as there are SCCs.
uint_fast64_t numberOfStates = stronglyConnectedComponents.size();
storm::storage::SparseMatrix<bool> sccDependencyGraph(numberOfStates);
sccDependencyGraph.initialize();
for (uint_fast64_t currentSccIndex = 0; currentSccIndex < stronglyConnectedComponents.size(); ++currentSccIndex) {
// Get the actual SCC.
std::vector<uint_fast64_t> const& scc = stronglyConnectedComponents[currentSccIndex];
// Now, we determine the SCCs which are reachable (in one step) from the current SCC.
std::set<uint_fast64_t> allTargetSccs;
for (auto state : scc) {
for (uint_fast64_t rowIndex = (*nondeterministicChoiceIndices)[state]; rowIndex < (*nondeterministicChoiceIndices)[state + 1]; ++rowIndex) {
for (typename storm::storage::SparseMatrix<T>::ConstIndexIterator succIt = this->getTransitionMatrix()->constColumnIteratorBegin(rowIndex), succIte = this->getTransitionMatrix()->constColumnIteratorEnd(rowIndex); succIt != succIte; ++succIt) {
uint_fast64_t targetScc = stateToSccMap.find(*succIt)->second;
// We only need to consider transitions that are actually leaving the SCC.
if (targetScc != currentSccIndex) {
allTargetSccs.insert(targetScc);
}
}
}
}
// Now we can just enumerate all the target SCCs and insert the corresponding transitions.
for (auto targetScc : allTargetSccs) {
sccDependencyGraph.insertNextValue(currentSccIndex, targetScc, true);
}
}
// Finalize the matrix.
sccDependencyGraph.finalize(true);
return sccDependencyGraph;
}
/*!
* Retrieves the size of the internal representation of the model in memory.
* @return the size of the internal representation of the model in memory

32
src/models/GraphTransitions.h

@ -119,38 +119,6 @@ public:
}
private:
void initializeFromSccDecomposition(GraphTransitions<T> const& transitions, std::vector<std::vector<uint_fast64_t>> const& stronglyConnectedComponents, std::map<uint_fast64_t, uint_fast64_t> const& stateToSccMap) {
for (uint_fast64_t currentSccIndex = 0; currentSccIndex < stronglyConnectedComponents.size(); ++currentSccIndex) {
// Mark beginning of current SCC.
stateIndications[currentSccIndex] = successorList.size();
// Get the actual SCC.
std::vector<uint_fast64_t> const& scc = stronglyConnectedComponents[currentSccIndex];
// Now, we determine the SCCs which are reachable (in one step) from the current SCC.
std::set<uint_fast64_t> allTargetSccs;
for (auto state : scc) {
for (stateSuccessorIterator succIt = transitions.beginStateSuccessorsIterator(state), succIte = transitions.endStateSuccessorsIterator(state); succIt != succIte; ++succIt) {
uint_fast64_t targetScc = stateToSccMap.find(*succIt)->second;
// We only need to consider transitions that are actually leaving the SCC.
if (targetScc != currentSccIndex) {
allTargetSccs.insert(targetScc);
}
}
}
// Now we can just enumerate all the target SCCs and insert the corresponding transitions.
for (auto targetScc : allTargetSccs) {
successorList.push_back(targetScc);
}
}
// Put the sentinel element at the end and initialize the number of transitions.
stateIndications[numberOfStates] = successorList.size();
numberOfTransitions = successorList.size();
}
/*!
* Initializes this graph transitions object using the forward transition
* relation given by means of a sparse matrix.

79
src/storage/SparseMatrix.h

@ -6,6 +6,7 @@
#include <algorithm>
#include <iostream>
#include <iterator>
#include <set>
#include "boost/integer/integer_mask.hpp"
#include "src/exceptions/InvalidStateException.h"
@ -55,8 +56,8 @@ public:
* If we only want to iterate over the columns or values of the non-zero entries of
* a row, we can simply iterate over the array (part) itself.
*/
typedef const uint_fast64_t* const ConstIndexIterator;
typedef const T* const ConstValueIterator;
typedef uint_fast64_t const* ConstIndexIterator;
typedef T const* ConstValueIterator;
/*!
* Iterator class that is able to iterate over the non-zero elements of a matrix and return the position
@ -232,7 +233,7 @@ public:
*
* @param nonZeroEntries The number of non-zero entries this matrix is going to hold.
*/
void initialize(uint_fast64_t nonZeroEntries) {
void initialize(uint_fast64_t nonZeroEntries = 0) {
// Check whether initializing the matrix is safe.
if (internalStatus != MatrixStatus::UnInitialized) {
triggerErrorState();
@ -260,11 +261,15 @@ public:
}
}
}
/*!
* Sets the matrix element at the given row and column to the given value.
* Sets the matrix element at the given row and column to the given value. After all elements have been added,
* a call to finalize() is mandatory.
* NOTE: This is a linear setter. It must be called consecutively for each element,
* row by row *and* column by column.
* NOTE: This method is different from insertNextValue(...) in that the number of nonzero elements must be known
* in advance (and passed to initialize()), because adding elements will not automatically increase the size of the
* underlying storage.
*
* @param row The row in which the matrix element is to be set.
* @param col The column in which the matrix element is to be set.
@ -276,7 +281,7 @@ public:
if ((row > rowCount) || (col > colCount)) {
triggerErrorState();
LOG4CPLUS_ERROR(logger, "Trying to add a value at illegal position (" << row << ", " << col << ") in matrix of size (" << rowCount << ", " << colCount << ").");
throw storm::exceptions::OutOfRangeException("Trying to add a value at illegal position.");
throw storm::exceptions::OutOfRangeException() << "Trying to add a value at illegal position (" << row << ", " << col << ") in matrix of size (" << rowCount << ", " << colCount << ").";
}
// If we switched to another row, we have to adjust the missing
@ -294,14 +299,53 @@ public:
++currentSize;
}
/*!
* Inserts a value at the given row and column with the given value. After all elements have been added,
* a call to finalize() is mandatory.
* NOTE: This is a linear inserter. It must be called consecutively for each element, row by row *and* column by
* column.
* NOTE: This method is different from addNextValue(...) in that the number of nonzero elements need not be known
* in advance, because inserting elements will automatically increase the size of the underlying storage.
*
* @param row The row in which the matrix element is to be set.
* @param col The column in which the matrix element is to be set.
* @param value The value that is to be set.
*/
void insertNextValue(const uint_fast64_t row, const uint_fast64_t col, T const& value) {
// Check whether the given row and column positions are valid and throw
// error otherwise.
if ((row > rowCount) || (col > colCount)) {
triggerErrorState();
LOG4CPLUS_ERROR(logger, "Trying to insert a value at illegal position (" << row << ", " << col << ") in matrix of size (" << rowCount << ", " << colCount << ").");
throw storm::exceptions::OutOfRangeException() << "Trying to insert a value at illegal position (" << row << ", " << col << ") in matrix of size (" << rowCount << ", " << colCount << ").";
}
// If we switched to another row, we have to adjust the missing entries in the rowIndications array.
if (row != lastRow) {
for (uint_fast64_t i = lastRow + 1; i <= row; ++i) {
rowIndications[i] = currentSize;
}
lastRow = row;
}
// Finally, set the element and increase the current size.
valueStorage.push_back(value);
columnIndications.push_back(col);
++nonZeroEntryCount;
++currentSize;
}
/*
* Finalizes the sparse matrix to indicate that initialization has been
* completed and the matrix may now be used.
* Finalizes the sparse matrix to indicate that initialization has been completed and the matrix may now be used.
*
* @param pushSentinelElement A boolean flag that indicates whether the sentinel element is to be pushed or inserted
* at a fixed location. If the elements have been added to the matrix via insertNextElement, this needs to be true
* and false otherwise.
*/
void finalize() {
// Check whether it's safe to finalize the matrix and throw error
// otherwise.
void finalize(bool pushSentinelElement = false) {
// Check whether it's safe to finalize the matrix and throw error otherwise.
if (!isInitialized()) {
triggerErrorState();
LOG4CPLUS_ERROR(logger, "Trying to finalize an uninitialized matrix.");
@ -319,11 +363,14 @@ public:
}
}
// Set a sentinel element at the last position of the row_indications
// array. This eases iteration work, as now the indices of row i
// are always between row_indications[i] and row_indications[i + 1],
// also for the first and last row.
rowIndications[rowCount] = nonZeroEntryCount;
// Set a sentinel element at the last position of the row_indications array. This eases iteration work, as
// now the indices of row i are always between rowIndications[i] and rowIndications[i + 1], also for the
// first and last row.
if (pushSentinelElement) {
rowIndications.push_back(nonZeroEntryCount);
} else {
rowIndications[rowCount] = nonZeroEntryCount;
}
setState(MatrixStatus::ReadReady);
}

114
src/utility/GraphAnalyzer.h

@ -421,7 +421,7 @@ public:
}
/*!
* Performs a decomposition of the given nondetermistic model into its SCCs.
* Performs a decomposition of the given nondeterministic model into its SCCs.
*
* @param model The nondeterminstic model to decompose.
* @return A pair whose first component represents the SCCs and whose second component represents the dependency
@ -430,51 +430,75 @@ public:
template <typename T>
static std::pair<std::vector<std::vector<uint_fast64_t>>, storm::models::GraphTransitions<T>> performSccDecomposition(storm::models::AbstractNondeterministicModel<T> const& model) {
LOG4CPLUS_INFO(logger, "Computing SCC decomposition.");
std::pair<std::vector<std::vector<uint_fast64_t>>, storm::models::GraphTransitions<T>> sccDecomposition;
uint_fast64_t numberOfStates = model.getNumberOfStates();
// Set up the environment of Tarjan's algorithm.
std::vector<uint_fast64_t> tarjanStack;
tarjanStack.reserve(numberOfStates);
storm::storage::BitVector tarjanStackStates(numberOfStates);
std::vector<uint_fast64_t> stateIndices(numberOfStates);
std::vector<uint_fast64_t> lowlinks(numberOfStates);
storm::storage::BitVector visitedStates(numberOfStates);
std::map<uint_fast64_t, uint_fast64_t> stateToSccMap;
// Start the search for SCCs from every vertex in the graph structure, because there is.
uint_fast64_t currentIndex = 0;
for (uint_fast64_t state = 0; state < numberOfStates; ++state) {
if (!visitedStates.get(state)) {
performSccDecompositionHelper(state, currentIndex, stateIndices, lowlinks, tarjanStack, tarjanStackStates, visitedStates, model.getTransitionMatrix(), sccDecomposition.first, stateToSccMap);
}
}
// Finally, determine the dependency graph over the SCCs and return result.
sccDecomposition.second = model.extractSccDependencyGraph(sccDecomposition.first, stateToSccMap);
// Get the forward transition relation from the model to ease the search.
std::vector<uint_fast64_t> const& nondeterministicChoiceIndices = model.getNondeterministicChoiceIndices();
storm::models::GraphTransitions<T> forwardTransitions(model.getTransitionMatrix(), nondeterministicChoiceIndices, true);
// Perform the actual SCC decomposition based on the graph-transitions of the system.
std::pair<std::vector<std::vector<uint_fast64_t>>, storm::models::GraphTransitions<T>> result = performSccDecomposition(nondeterministicChoiceIndices.size() - 1, forwardTransitions);
LOG4CPLUS_INFO(logger, "Done computing SCC decomposition.");
return result;
return sccDecomposition;
}
/*!
* Performs a topological sort of the states of the system according to the given transitions.
*
* @param transitions The transitions of the graph structure.
* @param matrix A square matrix representing the transition relation of the system.
* @return A vector of indices that is a topological sort of the states.
*/
template <typename T>
static std::vector<uint_fast64_t> getTopologicalSort(storm::models::GraphTransitions<T> const& transitions) {
static std::vector<uint_fast64_t> getTopologicalSort(storm::storage::SparseMatrix<T> const& matrix) {
if (matrix.getRowCount() != matrix.getColumnCount()) {
LOG4CPLUS_ERROR(logger, "Provided matrix is required to be square.");
throw storm::exceptions::InvalidArgumentException() << "Provided matrix is required to be square.";
}
uint_fast64_t numberOfStates = matrix.getRowCount();
// Prepare the result. This relies on the matrix being square.
std::vector<uint_fast64_t> topologicalSort;
topologicalSort.reserve(transitions.getNumberOfStates());
topologicalSort.reserve(numberOfStates);
// Prepare the stacks needed for recursion.
std::vector<uint_fast64_t> recursionStack;
recursionStack.reserve(transitions.getNumberOfStates());
std::vector<typename storm::models::GraphTransitions<T>::stateSuccessorIterator> iteratorRecursionStack;
iteratorRecursionStack.reserve(transitions.getNumberOfStates());
recursionStack.reserve(matrix.getRowCount());
std::vector<typename storm::storage::SparseMatrix<T>::ConstIndexIterator> iteratorRecursionStack;
iteratorRecursionStack.reserve(numberOfStates);
// Perform a depth-first search over the given transitions and record states in the reverse order they were visited.
storm::storage::BitVector visitedStates(transitions.getNumberOfStates());
for (uint_fast64_t state = 0; state < transitions.getNumberOfStates(); ++state) {
storm::storage::BitVector visitedStates(numberOfStates);
for (uint_fast64_t state = 0; state < numberOfStates; ++state) {
if (!visitedStates.get(state)) {
recursionStack.push_back(state);
iteratorRecursionStack.push_back(transitions.beginStateSuccessorsIterator(state));
iteratorRecursionStack.push_back(matrix.constColumnIteratorBegin(state));
recursionStepForward:
while (!recursionStack.empty()) {
uint_fast64_t currentState = recursionStack.back();
typename storm::models::GraphTransitions<T>::stateSuccessorIterator currentIt = iteratorRecursionStack.back();
typename storm::storage::SparseMatrix<T>::ConstIndexIterator currentIt = iteratorRecursionStack.back();
visitedStates.set(currentState, true);
recursionStepBackward:
for (; currentIt != transitions.endStateSuccessorsIterator(currentState); ++currentIt) {
for (; currentIt != matrix.constColumnIteratorEnd(currentState); ++currentIt) {
if (!visitedStates.get(*currentIt)) {
// Put unvisited successor on top of our recursion stack and remember that.
recursionStack.push_back(*currentIt);
@ -484,7 +508,7 @@ public:
iteratorRecursionStack.push_back(currentIt + 1);
// Also, put initial value for iterator on corresponding recursion stack.
iteratorRecursionStack.push_back(transitions.beginStateSuccessorsIterator(*currentIt));
iteratorRecursionStack.push_back(matrix.constColumnIteratorBegin(*currentIt));
goto recursionStepForward;
}
@ -514,40 +538,6 @@ public:
}
private:
/*!
* Performs an SCC decomposition of the system given by its forward transitions.
*
* @param forwardTransitions The (forward) transition relation of the model to decompose.
* @return A pair whose first component represents the SCCs and whose second component represents the dependency
* graph of the SCCs.
*/
template <typename T>
static std::pair<std::vector<std::vector<uint_fast64_t>>, storm::models::GraphTransitions<T>> performSccDecomposition(storm::models::GraphTransitions<T> const& forwardTransitions) {
std::pair<std::vector<std::vector<uint_fast64_t>>, storm::models::GraphTransitions<T>> sccDecomposition;
uint_fast64_t numberOfStates = forwardTransitions.getNumberOfStates();
// Set up the environment of Tarjan's algorithm.
std::vector<uint_fast64_t> tarjanStack;
tarjanStack.reserve(numberOfStates);
storm::storage::BitVector tarjanStackStates(numberOfStates);
std::vector<uint_fast64_t> stateIndices(numberOfStates);
std::vector<uint_fast64_t> lowlinks(numberOfStates);
storm::storage::BitVector visitedStates(numberOfStates);
std::map<uint_fast64_t, uint_fast64_t> stateToSccMap;
// Start the search for SCCs from every vertex in the graph structure, because there is.
uint_fast64_t currentIndex = 0;
for (uint_fast64_t state = 0; state < numberOfStates; ++state) {
if (!visitedStates.get(state)) {
performSccDecompositionHelper(state, currentIndex, stateIndices, lowlinks, tarjanStack, tarjanStackStates, visitedStates, forwardTransitions, sccDecomposition.first, stateToSccMap);
}
}
// Finally, determine the dependency graph over the SCCs and return result.
sccDecomposition.second = storm::models::GraphTransitions<T>(forwardTransitions, sccDecomposition.first, stateToSccMap);
return sccDecomposition;
}
/*!
* Performs an SCC decomposition using Tarjan's algorithm.
*
@ -559,30 +549,30 @@ private:
* @param tarjanStack A stack used for Tarjan's algorithm.
* @param tarjanStackStates A bit vector that represents all states that are currently contained in tarjanStack.
* @param visitedStates A bit vector that stores all states that have already been visited.
* @param forwardTransitions The (forward) transition relation of the graph structure.
* @param matrix The transition matrix representing the graph structure.
* @param stronglyConnectedComponents A vector of strongly connected components to which newly found SCCs are added.
* @param stateToSccMap A mapping from state indices to SCC indices that maps each state to its SCC.
*/
template <typename T>
static void performSccDecompositionHelper(uint_fast64_t startState, 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, storm::models::GraphTransitions<T> const& forwardTransitions, std::vector<std::vector<uint_fast64_t>>& stronglyConnectedComponents, std::map<uint_fast64_t, uint_fast64_t>& stateToSccMap) {
static void performSccDecompositionHelper(uint_fast64_t startState, 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, storm::storage::SparseMatrix<T> const& matrix, std::vector<std::vector<uint_fast64_t>>& stronglyConnectedComponents, std::map<uint_fast64_t, uint_fast64_t>& stateToSccMap) {
// 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.
std::vector<uint_fast64_t> recursionStateStack;
recursionStateStack.reserve(lowlinks.size());
std::vector<typename storm::models::GraphTransitions<T>::stateSuccessorIterator> recursionIteratorStack;
std::vector<typename storm::storage::SparseMatrix<T>::ConstIndexIterator> recursionIteratorStack;
recursionIteratorStack.reserve(lowlinks.size());
std::vector<bool> statesInStack(lowlinks.size());
// Initialize the recursion stacks with the given initial state (and its successor iterator).
recursionStateStack.push_back(startState);
recursionIteratorStack.push_back(forwardTransitions.beginStateSuccessorsIterator(startState));
recursionIteratorStack.push_back(matrix.constColumnIteratorBegin(startState));
recursionStepForward:
while (!recursionStateStack.empty()) {
uint_fast64_t currentState = recursionStateStack.back();
typename storm::models::GraphTransitions<T>::stateSuccessorIterator currentIt = recursionIteratorStack.back();
typename storm::storage::SparseMatrix<T>::ConstIndexIterator currentIt = recursionIteratorStack.back();
// Perform the treatment of newly discovered state as defined by Tarjan's algorithm
visitedStates.set(currentState, true);
@ -593,7 +583,7 @@ private:
tarjanStackStates.set(currentState, true);
// Now, traverse all successors of the current state.
for(; currentIt != forwardTransitions.endStateSuccessorsIterator(currentState); ++currentIt) {
for(; currentIt != matrix.constColumnIteratorEnd(currentState); ++currentIt) {
// If we have not visited the successor already, we need to perform the procedure
// recursively on the newly found state.
if (!visitedStates.get(*currentIt)) {
@ -606,7 +596,7 @@ private:
statesInStack[*currentIt] = true;
// Also, put initial value for iterator on corresponding recursion stack.
recursionIteratorStack.push_back(forwardTransitions.beginStateSuccessorsIterator(*currentIt));
recursionIteratorStack.push_back(matrix.constColumnIteratorBegin(*currentIt));
// Perform the actual recursion step in an iterative way.
goto recursionStepForward;

Loading…
Cancel
Save