Browse Source

Some minor fixes to GraphAnalyzer and model checkers.

tempestpy_adaptions
dehnert 12 years ago
parent
commit
5d849018de
  1. 19
      src/modelchecker/EigenDtmcPrctlModelChecker.h
  2. 21
      src/modelchecker/GmmxxDtmcPrctlModelChecker.h
  3. 120
      src/utility/GraphAnalyzer.h

19
src/modelchecker/EigenDtmcPrctlModelChecker.h

@ -126,17 +126,16 @@ public:
// Then, we need to identify the states which have to be taken out of the matrix, i.e. // Then, we need to identify the states which have to be taken out of the matrix, i.e.
// all states that have probability 0 and 1 of satisfying the until-formula. // all states that have probability 0 and 1 of satisfying the until-formula.
storm::storage::BitVector notExistsPhiUntilPsiStates(this->getModel().getNumberOfStates());
storm::storage::BitVector alwaysPhiUntilPsiStates(this->getModel().getNumberOfStates());
storm::utility::GraphAnalyzer::getPhiUntilPsiStates<double>(this->getModel(), *leftStates, *rightStates, &notExistsPhiUntilPsiStates, &alwaysPhiUntilPsiStates);
notExistsPhiUntilPsiStates.complement();
storm::storage::BitVector statesWithProbability0(this->getModel().getNumberOfStates());
storm::storage::BitVector statesWithProbability1(this->getModel().getNumberOfStates());
storm::utility::GraphAnalyzer::performProb01(this->getModel(), *leftStates, *rightStates, &statesWithProbability0, &statesWithProbability1);
delete leftStates; delete leftStates;
delete rightStates; delete rightStates;
LOG4CPLUS_INFO(logger, "Found " << notExistsPhiUntilPsiStates.getNumberOfSetBits() << " 'no' states.");
LOG4CPLUS_INFO(logger, "Found " << alwaysPhiUntilPsiStates.getNumberOfSetBits() << " 'yes' states.");
storm::storage::BitVector maybeStates = ~(notExistsPhiUntilPsiStates | alwaysPhiUntilPsiStates);
LOG4CPLUS_INFO(logger, "Found " << statesWithProbability0.getNumberOfSetBits() << " 'no' states.");
LOG4CPLUS_INFO(logger, "Found " << statesWithProbability1.getNumberOfSetBits() << " 'yes' states.");
storm::storage::BitVector maybeStates = ~(statesWithProbability0 | statesWithProbability1);
LOG4CPLUS_INFO(logger, "Found " << maybeStates.getNumberOfSetBits() << " 'maybe' states."); LOG4CPLUS_INFO(logger, "Found " << maybeStates.getNumberOfSetBits() << " 'maybe' states.");
// Create resulting vector and set values accordingly. // Create resulting vector and set values accordingly.
@ -175,7 +174,7 @@ public:
Type *pb = &(b[0]); // get the address storing the data for b Type *pb = &(b[0]); // get the address storing the data for b
MapType vectorB(pb, b.size()); // vectorB shares data MapType vectorB(pb, b.size()); // vectorB shares data
this->getModel().getTransitionProbabilityMatrix()->getConstrainedRowCountVector(maybeStates, alwaysPhiUntilPsiStates, &x);
this->getModel().getTransitionProbabilityMatrix()->getConstrainedRowCountVector(maybeStates, statesWithProbability1, &x);
Eigen::BiCGSTAB<Eigen::SparseMatrix<Type, 1, int_fast32_t>> solver; Eigen::BiCGSTAB<Eigen::SparseMatrix<Type, 1, int_fast32_t>> solver;
solver.compute(*eigenSubMatrix); solver.compute(*eigenSubMatrix);
@ -212,8 +211,8 @@ public:
delete eigenSubMatrix; delete eigenSubMatrix;
} }
storm::utility::setVectorValues<Type>(result, notExistsPhiUntilPsiStates, storm::utility::constGetZero<Type>());
storm::utility::setVectorValues<Type>(result, alwaysPhiUntilPsiStates, storm::utility::constGetOne<Type>());
storm::utility::setVectorValues<Type>(result, statesWithProbability0, storm::utility::constGetZero<Type>());
storm::utility::setVectorValues<Type>(result, statesWithProbability1, storm::utility::constGetOne<Type>());
return result; return result;
} }

21
src/modelchecker/GmmxxDtmcPrctlModelChecker.h

@ -111,18 +111,17 @@ public:
// Then, we need to identify the states which have to be taken out of the matrix, i.e. // Then, we need to identify the states which have to be taken out of the matrix, i.e.
// all states that have probability 0 and 1 of satisfying the until-formula. // all states that have probability 0 and 1 of satisfying the until-formula.
storm::storage::BitVector notExistsPhiUntilPsiStates(this->getModel().getNumberOfStates());
storm::storage::BitVector alwaysPhiUntilPsiStates(this->getModel().getNumberOfStates());
storm::utility::GraphAnalyzer::getPhiUntilPsiStates(this->getModel(), *leftStates, *rightStates, &notExistsPhiUntilPsiStates, &alwaysPhiUntilPsiStates);
notExistsPhiUntilPsiStates.complement();
storm::storage::BitVector statesWithProbability0(this->getModel().getNumberOfStates());
storm::storage::BitVector statesWithProbability1(this->getModel().getNumberOfStates());
storm::utility::GraphAnalyzer::performProb01(this->getModel(), *leftStates, *rightStates, &statesWithProbability0, &statesWithProbability1);
// Delete sub-results that are obsolete now. // Delete sub-results that are obsolete now.
delete leftStates; delete leftStates;
delete rightStates; delete rightStates;
LOG4CPLUS_INFO(logger, "Found " << notExistsPhiUntilPsiStates.getNumberOfSetBits() << " 'no' states.");
LOG4CPLUS_INFO(logger, "Found " << alwaysPhiUntilPsiStates.getNumberOfSetBits() << " 'yes' states.");
storm::storage::BitVector maybeStates = ~(notExistsPhiUntilPsiStates | alwaysPhiUntilPsiStates);
LOG4CPLUS_INFO(logger, "Found " << statesWithProbability0.getNumberOfSetBits() << " 'no' states.");
LOG4CPLUS_INFO(logger, "Found " << statesWithProbability1.getNumberOfSetBits() << " 'yes' states.");
storm::storage::BitVector maybeStates = ~(statesWithProbability0 | statesWithProbability1);
LOG4CPLUS_INFO(logger, "Found " << maybeStates.getNumberOfSetBits() << " 'maybe' states."); LOG4CPLUS_INFO(logger, "Found " << maybeStates.getNumberOfSetBits() << " 'maybe' states.");
// Create resulting vector. // Create resulting vector.
@ -149,7 +148,7 @@ public:
// Prepare the right-hand side of the equation system. For entry i this corresponds to // Prepare the right-hand side of the equation system. For entry i this corresponds to
// the accumulated probability of going from state i to some 'yes' state. // the accumulated probability of going from state i to some 'yes' state.
std::vector<Type> b(mayBeStatesSetBitCount); std::vector<Type> b(mayBeStatesSetBitCount);
this->getModel().getTransitionMatrix()->getConstrainedRowCountVector(maybeStates, alwaysPhiUntilPsiStates, &b);
this->getModel().getTransitionMatrix()->getConstrainedRowCountVector(maybeStates, statesWithProbability1, &b);
// Solve the corresponding system of linear equations. // Solve the corresponding system of linear equations.
this->solveLinearEquationSystem(*gmmxxMatrix, x, b); this->solveLinearEquationSystem(*gmmxxMatrix, x, b);
@ -162,8 +161,8 @@ public:
} }
// Set values of resulting vector that are known exactly. // Set values of resulting vector that are known exactly.
storm::utility::setVectorValues<Type>(result, notExistsPhiUntilPsiStates, storm::utility::constGetZero<Type>());
storm::utility::setVectorValues<Type>(result, alwaysPhiUntilPsiStates, storm::utility::constGetOne<Type>());
storm::utility::setVectorValues<Type>(result, statesWithProbability0, storm::utility::constGetZero<Type>());
storm::utility::setVectorValues<Type>(result, statesWithProbability1, storm::utility::constGetOne<Type>());
return result; return result;
} }
@ -253,7 +252,7 @@ public:
// Determine which states have a reward of infinity by definition. // Determine which states have a reward of infinity by definition.
storm::storage::BitVector infinityStates(this->getModel().getNumberOfStates()); storm::storage::BitVector infinityStates(this->getModel().getNumberOfStates());
storm::storage::BitVector trueStates(this->getModel().getNumberOfStates(), true); storm::storage::BitVector trueStates(this->getModel().getNumberOfStates(), true);
storm::utility::GraphAnalyzer::getAlwaysPhiUntilPsiStates(this->getModel(), trueStates, *targetStates, &infinityStates);
storm::utility::GraphAnalyzer::performProb1(this->getModel(), trueStates, *targetStates, &infinityStates);
infinityStates.complement(); infinityStates.complement();
// Create resulting vector. // Create resulting vector.

120
src/utility/GraphAnalyzer.h

@ -22,30 +22,58 @@ namespace utility {
class GraphAnalyzer { class GraphAnalyzer {
public: public:
/*!
* Computes the sets of states that have probability 0 or 1, respectively, of satisfying phi until psi in a
* deterministic model.
* @param model The model whose graph structure to search.
* @param phiStates The set of all states satisfying phi.
* @param psiStates The set of all states satisfying psi.
* @param statesWithProbability0 A pointer to a bit vector that is initially empty and will contain all states with
* probability 0 after the invocation of the function.
* @param statesWithProbability1 A pointer to a bit vector that is initially empty and will contain all states with
* probability 1 after the invocation of the function.
*/
template <class T>
static void performProb01(storm::models::AbstractDeterministicModel<T>& model, const storm::storage::BitVector& phiStates, const storm::storage::BitVector& psiStates, storm::storage::BitVector* statesWithProbability0, storm::storage::BitVector* statesWithProbability1) {
// Check for valid parameters.
if (statesWithProbability0 == nullptr) {
LOG4CPLUS_ERROR(logger, "Parameter 'statesWithProbability0' must not be null.");
throw storm::exceptions::InvalidArgumentException("Parameter 'statesWithProbability0' must not be null.");
}
if (statesWithProbability1 == nullptr) {
LOG4CPLUS_ERROR(logger, "Parameter 'statesWithProbability1' must not be null.");
throw storm::exceptions::InvalidArgumentException("Parameter 'statesWithProbability1' must not be null.");
}
// Perform the actual search.
GraphAnalyzer::performProbGreater0(model, phiStates, psiStates, statesWithProbability0);
GraphAnalyzer::performProb1(model, phiStates, psiStates, *statesWithProbability0, statesWithProbability1);
statesWithProbability0->complement();
}
/*! /*!
* Performs a backwards depth-first search trough the underlying graph structure * Performs a backwards depth-first search trough the underlying graph structure
* of the given model to determine which states of the model can reach one
* of the given target states whilst always staying in the set of filter states
* before. The resulting states are written to the given bit vector.
* of the given model to determine which states of the model have a positive probability
* of satisfying phi until psi. The resulting states are written to the given bit vector.
* @param model The model whose graph structure to search. * @param model The model whose graph structure to search.
* @param phiStates A bit vector of all states satisfying phi. * @param phiStates A bit vector of all states satisfying phi.
* @param psiStates A bit vector of all states satisfying psi. * @param psiStates A bit vector of all states satisfying psi.
* @param existsPhiUntilPsiStates A pointer to the result of the search for states that possess
* a paths satisfying phi until psi.
* @param statesWithProbabilityGreater0 A pointer to the result of the search for states that possess
* a positive probability of satisfying phi until psi.
*/ */
template <class T> template <class T>
static void getExistsPhiUntilPsiStates(storm::models::AbstractDeterministicModel<T>& model, const storm::storage::BitVector& phiStates, const storm::storage::BitVector& psiStates, storm::storage::BitVector* existsPhiUntilPsiStates) {
static void performProbGreater0(storm::models::AbstractDeterministicModel<T>& model, const storm::storage::BitVector& phiStates, const storm::storage::BitVector& psiStates, storm::storage::BitVector* statesWithProbabilityGreater0) {
// Check for valid parameter. // Check for valid parameter.
if (existsPhiUntilPsiStates == nullptr) {
LOG4CPLUS_ERROR(logger, "Parameter 'existsPhiUntilPhiStates' must not be null.");
throw storm::exceptions::InvalidArgumentException("Parameter 'existsPhiUntilPhiStates' must not be null.");
if (statesWithProbabilityGreater0 == nullptr) {
LOG4CPLUS_ERROR(logger, "Parameter 'statesWithProbabilityGreater0' must not be null.");
throw storm::exceptions::InvalidArgumentException("Parameter 'statesWithProbabilityGreater0' must not be null.");
} }
// Get the backwards transition relation from the model to ease the search. // Get the backwards transition relation from the model to ease the search.
storm::models::GraphTransitions<T>& backwardTransitions = model.getBackwardTransitions(); storm::models::GraphTransitions<T>& backwardTransitions = model.getBackwardTransitions();
// Add all psi states as the already satisfy the condition. // Add all psi states as the already satisfy the condition.
*existsPhiUntilPsiStates |= psiStates;
*statesWithProbabilityGreater0 |= psiStates;
// Initialize the stack used for the DFS with the states // Initialize the stack used for the DFS with the states
std::vector<uint_fast64_t> stack; std::vector<uint_fast64_t> stack;
@ -58,8 +86,8 @@ public:
stack.pop_back(); stack.pop_back();
for(auto it = backwardTransitions.beginStateSuccessorsIterator(currentState); it != backwardTransitions.endStateSuccessorsIterator(currentState); ++it) { for(auto it = backwardTransitions.beginStateSuccessorsIterator(currentState); it != backwardTransitions.endStateSuccessorsIterator(currentState); ++it) {
if (phiStates.get(*it) && !existsPhiUntilPsiStates->get(*it)) {
existsPhiUntilPsiStates->set(*it, true);
if (phiStates.get(*it) && !statesWithProbabilityGreater0->get(*it)) {
statesWithProbabilityGreater0->set(*it, true);
stack.push_back(*it); stack.push_back(*it);
} }
} }
@ -75,27 +103,29 @@ public:
* @param model The model whose graph structure to search. * @param model The model whose graph structure to search.
* @param phiStates A bit vector of all states satisfying phi. * @param phiStates A bit vector of all states satisfying phi.
* @param psiStates A bit vector of all states satisfying psi. * @param psiStates A bit vector of all states satisfying psi.
* @param existsPhiUntilPsiStates A reference to a bit vector of states that possess a path
* satisfying phi until psi.
* @param statesWithProbabilityGreater0 A reference to a bit vector of states that possess a positive
* probability mass of satisfying phi until psi.
* @param alwaysPhiUntilPsiStates A pointer to the result of the search for states that only * @param alwaysPhiUntilPsiStates A pointer to the result of the search for states that only
* have paths satisfying phi until psi. * have paths satisfying phi until psi.
*/ */
template <class T> template <class T>
static void getAlwaysPhiUntilPsiStates(storm::models::AbstractDeterministicModel<T>& model, const storm::storage::BitVector& phiStates, const storm::storage::BitVector& psiStates, const storm::storage::BitVector& existsPhiUntilPsiStates, storm::storage::BitVector* alwaysPhiUntilPsiStates) {
static void performProb1(storm::models::AbstractDeterministicModel<T>& model, const storm::storage::BitVector& phiStates, const storm::storage::BitVector& psiStates, const storm::storage::BitVector& statesWithProbabilityGreater0, storm::storage::BitVector* statesWithProbability1) {
// Check for valid parameter. // Check for valid parameter.
if (alwaysPhiUntilPsiStates == nullptr) {
LOG4CPLUS_ERROR(logger, "Parameter 'alwaysPhiUntilPhiStates' must not be null.");
throw storm::exceptions::InvalidArgumentException("Parameter 'alwaysPhiUntilPhiStates' must not be null.");
if (statesWithProbability1 == nullptr) {
LOG4CPLUS_ERROR(logger, "Parameter 'statesWithProbability1' must not be null.");
throw storm::exceptions::InvalidArgumentException("Parameter 'statesWithProbability1' must not be null.");
} }
GraphAnalyzer::getExistsPhiUntilPsiStates(model, ~psiStates, ~existsPhiUntilPsiStates, alwaysPhiUntilPsiStates);
alwaysPhiUntilPsiStates->complement();
GraphAnalyzer::performProbGreater0(model, ~psiStates, ~statesWithProbabilityGreater0, statesWithProbability1);
statesWithProbability1->complement();
} }
/*! /*!
* Computes the set of states of the given model for which all paths lead to * Computes the set of states of the given model for which all paths lead to
* the given set of target states and only visit states from the filter set * the given set of target states and only visit states from the filter set
* before. The results are written to the given bit vector.
* before. In order to do this, it uses the given set of states that
* characterizes the states that possess at least one path to a target state.
* The results are written to the given bit vector.
* @param model The model whose graph structure to search. * @param model The model whose graph structure to search.
* @param phiStates A bit vector of all states satisfying phi. * @param phiStates A bit vector of all states satisfying phi.
* @param psiStates A bit vector of all states satisfying psi. * @param psiStates A bit vector of all states satisfying psi.
@ -103,49 +133,19 @@ public:
* have paths satisfying phi until psi. * have paths satisfying phi until psi.
*/ */
template <class T> template <class T>
static void getAlwaysPhiUntilPsiStates(storm::models::AbstractDeterministicModel<T>& model, const storm::storage::BitVector& phiStates, const storm::storage::BitVector& psiStates, storm::storage::BitVector* alwaysPhiUntilPsiStates) {
static void performProb1(storm::models::AbstractDeterministicModel<T>& model, const storm::storage::BitVector& phiStates, const storm::storage::BitVector& psiStates, storm::storage::BitVector* statesWithProbability1) {
// Check for valid parameter. // Check for valid parameter.
if (alwaysPhiUntilPsiStates == nullptr) {
LOG4CPLUS_ERROR(logger, "Parameter 'alwaysPhiUntilPhiStates' must not be null.");
throw storm::exceptions::InvalidArgumentException("Parameter 'alwaysPhiUntilPhiStates' must not be null.");
if (statesWithProbability1 == nullptr) {
LOG4CPLUS_ERROR(logger, "Parameter 'statesWithProbability1' must not be null.");
throw storm::exceptions::InvalidArgumentException("Parameter 'statesWithProbability1' must not be null.");
} }
storm::storage::BitVector existsPhiUntilPsiStates(model.getNumberOfStates());
GraphAnalyzer::getExistsPhiUntilPsiStates(model, phiStates, psiStates, &existsPhiUntilPsiStates);
GraphAnalyzer::getExistsPhiUntilPsiStates(model, ~psiStates, ~existsPhiUntilPsiStates, alwaysPhiUntilPsiStates);
alwaysPhiUntilPsiStates->complement();
storm::storage::BitVector* statesWithProbabilityGreater0 = new storm::storage::BitVector(model.getNumberOfStates());
GraphAnalyzer::performProbGreater0(model, phiStates, psiStates, statesWithProbabilityGreater0);
GraphAnalyzer::performProbGreater0(model, ~psiStates, ~(*statesWithProbabilityGreater0), statesWithProbability1);
delete statesWithProbabilityGreater0;
statesWithProbability1->complement();
} }
/*!
* Computes the set of states of the given model for which all paths lead to
* the given set of target states and only visit states from the filter set
* before.
* @param model The model whose graph structure to search.
* @param phiStates A bit vector of all states satisfying phi.
* @param psiStates A bit vector of all states satisfying psi.
* @param existsPhiUntilPsiStates A pointer to the result of the search for states that possess
* a path satisfying phi until psi.
* @param alwaysPhiUntilPsiStates A pointer to the result of the search for states that only
* have paths satisfying phi until psi.
*/
template <class T>
static void getPhiUntilPsiStates(storm::models::AbstractDeterministicModel<T>& model, const storm::storage::BitVector& phiStates, const storm::storage::BitVector& psiStates, storm::storage::BitVector* existsPhiUntilPsiStates, storm::storage::BitVector* alwaysPhiUntilPsiStates) {
// Check for valid parameters.
if (existsPhiUntilPsiStates == nullptr) {
LOG4CPLUS_ERROR(logger, "Parameter 'existsPhiUntilPhiStates' must not be null.");
throw storm::exceptions::InvalidArgumentException("Parameter 'existsPhiUntilPhiStates' must not be null.");
}
if (alwaysPhiUntilPsiStates == nullptr) {
LOG4CPLUS_ERROR(logger, "Parameter 'alwaysPhiUntilPhiStates' must not be null.");
throw storm::exceptions::InvalidArgumentException("Parameter 'alwaysPhiUntilPhiStates' must not be null.");
}
// Perform search.
GraphAnalyzer::getExistsPhiUntilPsiStates(model, phiStates, psiStates, existsPhiUntilPsiStates);
GraphAnalyzer::getAlwaysPhiUntilPsiStates(model, phiStates, psiStates, *existsPhiUntilPsiStates, alwaysPhiUntilPsiStates);
}
}; };
} // namespace utility } // namespace utility

Loading…
Cancel
Save