Browse Source

fixed issue related to row groups in sparse matrix and adapted the affected calling sites

Former-commit-id: 96c6fd7e59
tempestpy_adaptions
dehnert 9 years ago
parent
commit
6a80348150
  1. 21
      src/parser/NondeterministicSparseTransitionParser.cpp
  2. 33
      src/solver/GameSolver.cpp
  3. 8
      src/solver/SymbolicGameSolver.cpp
  4. 32
      src/storage/SparseMatrix.cpp
  5. 1
      src/storage/prism/Program.cpp
  6. 4
      src/utility/ConstantsComparator.cpp
  7. 42
      src/utility/constants.cpp
  8. 13
      src/utility/solver.cpp
  9. 16
      src/utility/solver.h

21
src/parser/NondeterministicSparseTransitionParser.cpp

@ -117,15 +117,21 @@ namespace storm {
// If we have switched the source state, we possibly need to insert the rows of the last
// source state.
if (source != lastSource) {
curRow += ((modelInformation.getRowGroupIndices())[lastSource + 1] - (modelInformation.getRowGroupIndices())[lastSource]) -(lastChoice + 1);
curRow += ((modelInformation.getRowGroupIndices())[lastSource + 1] - (modelInformation.getRowGroupIndices())[lastSource]) - (lastChoice + 1);
}
// If we skipped some states, we need to reserve empty rows for all their nondeterministic
// choices.
// choices and create the row groups.
for (uint_fast64_t i = lastSource + 1; i < source; ++i) {
matrixBuilder.newRowGroup(modelInformation.getRowGroupIndices()[i]);
curRow += ((modelInformation.getRowGroupIndices())[i + 1] - (modelInformation.getRowGroupIndices())[i]);
}
// If we moved to the next source, we need to open the next row group.
if (source != lastSource) {
matrixBuilder.newRowGroup(modelInformation.getRowGroupIndices()[source]);
}
// If we advanced to the next state, but skipped some choices, we have to reserve rows
// for them
if (source != lastSource) {
@ -155,7 +161,7 @@ namespace storm {
}
}
if (source != lastSource) {
// Add create a new rowGroup for the source, if this is the first choice we encounter for this state.
// Create a new rowGroup for the source, if this is the first choice we encounter for this state.
matrixBuilder.newRowGroup(curRow);
}
}
@ -178,13 +184,8 @@ namespace storm {
// Since we assume the transition rewards are for the transitions of the model, we copy the rowGroupIndices.
if (isRewardFile) {
// We already have rowGroup 0.
for (uint_fast64_t index = 1; index < modelInformation.getRowGroupIndices().size() - 1; index++) {
matrixBuilder.newRowGroup(modelInformation.getRowGroupIndices()[index]);
}
} else {
for (uint_fast64_t node = lastSource + 1; node <= firstPass.highestStateIndex; node++) {
matrixBuilder.newRowGroup(curRow + 1);
for (uint_fast64_t node = lastSource + 1; node < modelInformation.getRowGroupCount(); node++) {
matrixBuilder.newRowGroup(modelInformation.getRowGroupIndices()[node]);
}
}

33
src/solver/GameSolver.cpp

@ -39,26 +39,27 @@ namespace storm {
}
for (uint_fast64_t pl1State = 0; pl1State < numberOfPlayer1States; ++pl1State) {
uint_fast64_t startRow = player1Matrix.getRowGroupIndices()[pl1State];
uint_fast64_t endRow = player1Matrix.getRowGroupIndices()[pl1State + 1];
storm::storage::SparseMatrix<storm::storage::sparse::state_type>::const_rows relevantRows = player1Matrix.getRowGroup(pl1State);
storm::storage::SparseMatrix<storm::storage::sparse::state_type>::const_iterator it = relevantRows.begin();
storm::storage::SparseMatrix<storm::storage::sparse::state_type>::const_iterator ite = relevantRows.end();
if (relevantRows.getNumberOfEntries() > 0) {
storm::storage::SparseMatrix<storm::storage::sparse::state_type>::const_iterator it = relevantRows.begin();
storm::storage::SparseMatrix<storm::storage::sparse::state_type>::const_iterator ite = relevantRows.end();
// Set the first value.
tmpResult[pl1State] = player2Result[it->getColumn()];
++it;
// Set the first value.
tmpResult[pl1State] = player2Result[it->getColumn()];
++it;
// Now iterate through the different values and pick the extremal one.
if (player1Goal == OptimizationDirection::Minimize) {
for (; it != ite; ++it) {
tmpResult[pl1State] = std::min(tmpResult[pl1State], player2Result[it->getColumn()]);
// Now iterate through the different values and pick the extremal one.
if (player1Goal == OptimizationDirection::Minimize) {
for (; it != ite; ++it) {
tmpResult[pl1State] = std::min(tmpResult[pl1State], player2Result[it->getColumn()]);
}
} else {
for (; it != ite; ++it) {
tmpResult[pl1State] = std::max(tmpResult[pl1State], player2Result[it->getColumn()]);
}
}
} else {
for (; it != ite; ++it) {
tmpResult[pl1State] = std::max(tmpResult[pl1State], player2Result[it->getColumn()]);
}
tmpResult[pl1State] = storm::utility::zero<ValueType>();
}
}
@ -67,7 +68,7 @@ namespace storm {
std::swap(x, tmpResult);
++iterations;
} while (!converged);
} while (!converged && iterations < maximalNumberOfIterations);
}
template class GameSolver<double>;

8
src/solver/SymbolicGameSolver.cpp

@ -26,7 +26,7 @@ namespace storm {
uint_fast64_t iterations = 0;
bool converged = false;
while (!converged && iterations < maximalNumberOfIterations) {
do {
// Compute tmp = A * x + b
storm::dd::Add<Type> xCopyAsColumn = xCopy.swapVariables(this->rowColumnMetaVariablePairs);
storm::dd::Add<Type> tmp = this->gameMatrix.multiplyMatrix(xCopyAsColumn, this->columnMetaVariables);
@ -34,12 +34,12 @@ namespace storm {
// Now abstract from player 2 and player 1 variables.
switch (player2Goal) {
case OptimizationDirection::Minimize: tmp.minAbstract(this->player2Variables); break;
case OptimizationDirection::Minimize: tmp = tmp.minAbstract(this->player2Variables); break;
case OptimizationDirection::Maximize: tmp = tmp.maxAbstract(this->player2Variables); break;
}
switch (player1Goal) {
case OptimizationDirection::Minimize: tmp.minAbstract(this->player1Variables); break;
case OptimizationDirection::Minimize: tmp = tmp.minAbstract(this->player1Variables); break;
case OptimizationDirection::Maximize: tmp = tmp.maxAbstract(this->player1Variables); break;
}
@ -52,7 +52,7 @@ namespace storm {
}
++iterations;
}
} while (!converged && iterations < maximalNumberOfIterations);
return xCopy;
}

32
src/storage/SparseMatrix.cpp

@ -1,4 +1,5 @@
#include <boost/functional/hash.hpp>
#include <src/storage/sparse/StateType.h>
// To detect whether the usage of TBB is possible, this include is neccessary
#include "storm-config.h"
@ -143,9 +144,18 @@ namespace storm {
template<typename ValueType>
void SparseMatrixBuilder<ValueType>::newRowGroup(index_type startingRow) {
STORM_LOG_THROW(hasCustomRowGrouping, storm::exceptions::InvalidStateException, "Matrix was not created to have a custom row grouping.");
STORM_LOG_THROW(rowGroupIndices.empty() || startingRow >= rowGroupIndices.back(), storm::exceptions::InvalidStateException, "Illegal row group with negative size.");
STORM_LOG_THROW(startingRow >= lastRow, storm::exceptions::InvalidStateException, "Illegal row group with negative size.");
rowGroupIndices.push_back(startingRow);
++currentRowGroup;
// Close all rows from the most recent one to the starting row.
for (index_type i = lastRow + 1; i <= startingRow; ++i) {
rowIndications.push_back(currentEntryCount);
}
// Reset the most recently seen row/column to allow for proper insertion of the following elements.
lastRow = startingRow;
lastColumn = 0;
}
template<typename ValueType>
@ -156,7 +166,7 @@ namespace storm {
rowCount = std::max(rowCount, initialRowCount);
}
rowCount = std::max(rowCount, overriddenRowCount);
// If the current row count was overridden, we may need to add empty rows.
for (index_type i = lastRow + 1; i < rowCount; ++i) {
rowIndications.push_back(currentEntryCount);
@ -196,7 +206,7 @@ namespace storm {
rowGroupIndices.push_back(rowCount);
}
}
return SparseMatrix<ValueType>(columnCount, std::move(rowIndications), std::move(columnsAndValues), std::move(rowGroupIndices), hasCustomRowGrouping);
}
@ -1163,16 +1173,16 @@ namespace storm {
// Explicitly instantiate the entry, builder and the matrix.
// double
template class MatrixEntry<typename SparseMatrix<double>::index_type, double>;
template std::ostream& operator<<(std::ostream& out, MatrixEntry<uint_fast64_t, double> const& entry);
template std::ostream& operator<<(std::ostream& out, MatrixEntry<typename SparseMatrix<double>::index_type, double> const& entry);
template class SparseMatrixBuilder<double>;
template class SparseMatrix<double>;
template std::ostream& operator<<(std::ostream& out, SparseMatrix<double> const& matrix);
template std::vector<double> SparseMatrix<double>::getPointwiseProductRowSumVector(storm::storage::SparseMatrix<double> const& otherMatrix) const;
template bool SparseMatrix<double>::isSubmatrixOf(SparseMatrix<double> const& matrix) const;
// float
// float
template class MatrixEntry<typename SparseMatrix<float>::index_type, float>;
template std::ostream& operator<<(std::ostream& out, MatrixEntry<uint_fast64_t, float> const& entry);
template std::ostream& operator<<(std::ostream& out, MatrixEntry<typename SparseMatrix<float>::index_type, float> const& entry);
template class SparseMatrixBuilder<float>;
template class SparseMatrix<float>;
template std::ostream& operator<<(std::ostream& out, SparseMatrix<float> const& matrix);
@ -1181,12 +1191,20 @@ namespace storm {
// int
template class MatrixEntry<typename SparseMatrix<int>::index_type, int>;
template std::ostream& operator<<(std::ostream& out, MatrixEntry<uint_fast64_t, int> const& entry);
template std::ostream& operator<<(std::ostream& out, MatrixEntry<typename SparseMatrix<int>::index_type, int> const& entry);
template class SparseMatrixBuilder<int>;
template class SparseMatrix<int>;
template std::ostream& operator<<(std::ostream& out, SparseMatrix<int> const& matrix);
template bool SparseMatrix<int>::isSubmatrixOf(SparseMatrix<int> const& matrix) const;
// state_type
template class MatrixEntry<typename SparseMatrix<storm::storage::sparse::state_type>::index_type, storm::storage::sparse::state_type>;
template std::ostream& operator<<(std::ostream& out, MatrixEntry<typename SparseMatrix<storm::storage::sparse::state_type>::index_type, storm::storage::sparse::state_type> const& entry);
template class SparseMatrixBuilder<storm::storage::sparse::state_type>;
template class SparseMatrix<storm::storage::sparse::state_type>;
template std::ostream& operator<<(std::ostream& out, SparseMatrix<storm::storage::sparse::state_type> const& matrix);
template bool SparseMatrix<int>::isSubmatrixOf(SparseMatrix<storm::storage::sparse::state_type> const& matrix) const;
#ifdef STORM_HAVE_CARL
// Rat Function
template class MatrixEntry<typename SparseMatrix<RationalFunction>::index_type, RationalFunction>;

1
src/storage/prism/Program.cpp

@ -1011,7 +1011,6 @@ namespace storm {
}
// Now we are in a position to start the enumeration over all command variables.
uint_fast64_t modelCount = 0;
solver->allSat(allCommandVariables, [&] (storm::solver::SmtSolver::ModelReference& modelReference) -> bool {
// Now we need to reconstruct the chosen commands from the valuation of the command variables.
std::vector<std::vector<std::reference_wrapper<Command const>>> chosenCommands(possibleCommands.size());

4
src/utility/ConstantsComparator.cpp

@ -2,6 +2,7 @@
#include <cstdlib>
#include <cmath>
#include <src/storage/sparse/StateType.h>
#include "src/utility/constants.h"
#include "src/settings/SettingsManager.h"
@ -94,7 +95,8 @@ namespace storm {
template class ConstantsComparator<double>;
template class ConstantsComparator<float>;
template class ConstantsComparator<int>;
template class ConstantsComparator<storm::storage::sparse::state_type>;
#ifdef STORM_HAVE_CARL
template class ConstantsComparator<RationalFunction>;
template class ConstantsComparator<Polynomial>;

42
src/utility/constants.cpp

@ -1,7 +1,7 @@
#include "src/utility/constants.h"
#include "src/storage/SparseMatrix.h"
#include "src/storage/sparse/StateType.h"
#include "src/storage/SparseMatrix.h"
#include "src/settings/SettingsManager.h"
#include "src/settings/modules/GeneralSettings.h"
@ -83,27 +83,13 @@ namespace storm {
return std::pow(value, exponent);
}
template<>
double simplify(double value) {
// In the general case, we don't to anything here, but merely return the value. If something else is
// supposed to happen here, the templated function can be specialized for this particular type.
return value;
}
template<>
float simplify(float value) {
// In the general case, we don't to anything here, but merely return the value. If something else is
// supposed to happen here, the templated function can be specialized for this particular type.
template<typename ValueType>
ValueType simplify(ValueType value) {
// In the general case, we don't to anything here, but merely return the value. If something else is
// supposed to happen here, the templated function can be specialized for this particular type.
return value;
}
template<>
int simplify(int value) {
// In the general case, we don't to anything here, but merely return the value. If something else is
// supposed to happen here, the templated function can be specialized for this particular type.
return value;
}
#ifdef STORM_HAVE_CARL
template<>
RationalFunction& simplify(RationalFunction& value);
@ -201,7 +187,23 @@ namespace storm {
template storm::storage::MatrixEntry<storm::storage::sparse::state_type, int> simplify(storm::storage::MatrixEntry<storm::storage::sparse::state_type, int> matrixEntry);
template storm::storage::MatrixEntry<storm::storage::sparse::state_type, int>& simplify(storm::storage::MatrixEntry<storm::storage::sparse::state_type, int>& matrixEntry);
template storm::storage::MatrixEntry<storm::storage::sparse::state_type, int>&& simplify(storm::storage::MatrixEntry<storm::storage::sparse::state_type, int>&& matrixEntry);
template bool isOne(storm::storage::sparse::state_type const& value);
template bool isZero(storm::storage::sparse::state_type const& value);
template bool isConstant(storm::storage::sparse::state_type const& value);
template storm::storage::sparse::state_type one();
template storm::storage::sparse::state_type zero();
template storm::storage::sparse::state_type infinity();
template storm::storage::sparse::state_type pow(storm::storage::sparse::state_type const& value, uint_fast64_t exponent);
template storm::storage::sparse::state_type simplify(storm::storage::sparse::state_type value);
template storm::storage::MatrixEntry<storm::storage::sparse::state_type, storm::storage::sparse::state_type> simplify(storm::storage::MatrixEntry<storm::storage::sparse::state_type, storm::storage::sparse::state_type> matrixEntry);
template storm::storage::MatrixEntry<storm::storage::sparse::state_type, storm::storage::sparse::state_type>& simplify(storm::storage::MatrixEntry<storm::storage::sparse::state_type, storm::storage::sparse::state_type>& matrixEntry);
template storm::storage::MatrixEntry<storm::storage::sparse::state_type, storm::storage::sparse::state_type>&& simplify(storm::storage::MatrixEntry<storm::storage::sparse::state_type, storm::storage::sparse::state_type>&& matrixEntry);
#ifdef STORM_HAVE_CARL
template bool isOne(RationalFunction const& value);
template bool isZero(RationalFunction const& value);

13
src/utility/solver.cpp

@ -1,13 +1,13 @@
#include "src/utility/solver.h"
#include "src/solver/SymbolicGameSolver.h"
#include <vector>
#include "src/solver/SymbolicLinearEquationSolver.h"
#include "src/solver/SymbolicMinMaxLinearEquationSolver.h"
#include "src/solver/SymbolicGameSolver.h"
#include "src/solver/NativeLinearEquationSolver.h"
#include "src/solver/GmmxxLinearEquationSolver.h"
#include "src/solver/GameSolver.h"
#include "src/solver/NativeMinMaxLinearEquationSolver.h"
#include "src/solver/GmmxxMinMaxLinearEquationSolver.h"
@ -102,8 +102,7 @@ namespace storm {
template<typename ValueType>
void MinMaxLinearEquationSolverFactory<ValueType>::setPreferredTechnique(storm::solver::MinMaxTechniqueSelection preferredTech) {
this->prefTech = preferredTech;
}
}
template<typename ValueType>
std::unique_ptr<storm::solver::MinMaxLinearEquationSolver<ValueType>> MinMaxLinearEquationSolverFactory<ValueType>::create(storm::storage::SparseMatrix<ValueType> const& matrix, bool trackPolicy) const {
@ -133,6 +132,10 @@ namespace storm {
}
template<typename ValueType>
std::unique_ptr<storm::solver::GameSolver<ValueType>> GameSolverFactory<ValueType>::create(storm::storage::SparseMatrix<storm::storage::sparse::state_type> const& player1Matrix, storm::storage::SparseMatrix<ValueType> const& player2Matrix) const {
return std::unique_ptr<storm::solver::GameSolver<ValueType>>(new storm::solver::GameSolver<ValueType>(player1Matrix, player2Matrix));
}
std::unique_ptr<storm::solver::LpSolver> LpSolverFactory::create(std::string const& name, storm::solver::LpSolverTypeSelection solvT) const {
storm::solver::LpSolverType t;
@ -192,7 +195,7 @@ namespace storm {
template class GmmxxLinearEquationSolverFactory<double>;
template class NativeLinearEquationSolverFactory<double>;
template class MinMaxLinearEquationSolverFactory<double>;
template class GameSolverFactory<double>;
}
}
}

16
src/utility/solver.h

@ -4,6 +4,7 @@
#include <set>
#include <vector>
#include <memory>
#include <src/storage/sparse/StateType.h>
#include "src/storage/dd/DdType.h"
#include "src/solver/SolverSelectionOptions.h"
@ -13,13 +14,15 @@ namespace storm {
template<storm::dd::DdType T> class SymbolicGameSolver;
template<storm::dd::DdType T, typename V> class SymbolicLinearEquationSolver;
template<storm::dd::DdType T, typename V> class SymbolicMinMaxLinearEquationSolver;
template<typename V> class GameSolver;
template<typename V> class LinearEquationSolver;
template<typename V> class MinMaxLinearEquationSolver;
class LpSolver;
class SmtSolver;
template<typename ValueType> class NativeLinearEquationSolver;
enum class NativeLinearEquationSolverSolutionMethod;
enum class
NativeLinearEquationSolverSolutionMethod;
}
namespace storage {
@ -93,7 +96,7 @@ namespace storm {
public:
MinMaxLinearEquationSolverFactory(storm::solver::EquationSolverTypeSelection solverType = storm::solver::EquationSolverTypeSelection::FROMSETTINGS);
/*!
* Creates a new nondeterministic linear equation solver instance with the given matrix.
* Creates a new min/max linear equation solver instance with the given matrix.
*/
virtual std::unique_ptr<storm::solver::MinMaxLinearEquationSolver<ValueType>> create(storm::storage::SparseMatrix<ValueType> const& matrix, bool trackPolicy = false) const;
void setSolverType(storm::solver::EquationSolverTypeSelection solverType);
@ -106,6 +109,15 @@ namespace storm {
/// Notice that we save the selection enum here, which allows different solvers to use different techniques.
storm::solver::MinMaxTechniqueSelection prefTech;
};
template<typename ValueType>
class GameSolverFactory {
public:
/*!
* Creates a new game solver instance with the given matrices.
*/
virtual std::unique_ptr<storm::solver::GameSolver<ValueType>> create(storm::storage::SparseMatrix<storm::storage::sparse::state_type> const& player1Matrix, storm::storage::SparseMatrix<ValueType> const& player2Matrix) const;
};
class LpSolverFactory {
public:

Loading…
Cancel
Save