Browse Source
Added AcyclicMinMaxLinearEquationSolver and AcyclicLinearEquationSolver which are optimized for many calls on an acyclic model.
tempestpy_adaptions
Added AcyclicMinMaxLinearEquationSolver and AcyclicLinearEquationSolver which are optimized for many calls on an acyclic model.
tempestpy_adaptions
Tim Quatmann
5 years ago
15 changed files with 498 additions and 28 deletions
-
4src/storm/settings/modules/CoreSettings.cpp
-
4src/storm/settings/modules/MinMaxEquationSolverSettings.cpp
-
117src/storm/solver/AcyclicLinearEquationSolver.cpp
-
40src/storm/solver/AcyclicLinearEquationSolver.h
-
100src/storm/solver/AcyclicMinMaxLinearEquationSolver.cpp
-
50src/storm/solver/AcyclicMinMaxLinearEquationSolver.h
-
10src/storm/solver/LinearEquationSolver.cpp
-
30src/storm/solver/LinearEquationSolverRequirements.cpp
-
7src/storm/solver/LinearEquationSolverRequirements.h
-
5src/storm/solver/MinMaxLinearEquationSolver.cpp
-
25src/storm/solver/MinMaxLinearEquationSolverRequirements.cpp
-
9src/storm/solver/MinMaxLinearEquationSolverRequirements.h
-
4src/storm/solver/SolverSelectionOptions.cpp
-
4src/storm/solver/SolverSelectionOptions.h
-
117src/storm/solver/helper/AcyclicSolverHelper.cpp
@ -0,0 +1,117 @@ |
|||
#include "storm/solver/AcyclicLinearEquationSolver.h"
|
|||
|
|||
#include "storm/solver/helper/AcyclicSolverHelper.cpp"
|
|||
|
|||
#include "storm/utility/vector.h"
|
|||
|
|||
namespace storm { |
|||
namespace solver { |
|||
|
|||
template<typename ValueType> |
|||
AcyclicLinearEquationSolver<ValueType>::AcyclicLinearEquationSolver() { |
|||
// Intentionally left empty.
|
|||
} |
|||
|
|||
template<typename ValueType> |
|||
AcyclicLinearEquationSolver<ValueType>::AcyclicLinearEquationSolver(storm::storage::SparseMatrix<ValueType> const& A) { |
|||
this->setMatrix(A); |
|||
} |
|||
|
|||
template<typename ValueType> |
|||
AcyclicLinearEquationSolver<ValueType>::AcyclicLinearEquationSolver(storm::storage::SparseMatrix<ValueType>&& A) { |
|||
this->setMatrix(std::move(A)); |
|||
} |
|||
|
|||
template<typename ValueType> |
|||
void AcyclicLinearEquationSolver<ValueType>::setMatrix(storm::storage::SparseMatrix<ValueType> const& A) { |
|||
localA.reset(); |
|||
this->A = &A; |
|||
clearCache(); |
|||
} |
|||
|
|||
template<typename ValueType> |
|||
void AcyclicLinearEquationSolver<ValueType>::setMatrix(storm::storage::SparseMatrix<ValueType>&& A) { |
|||
localA = std::make_unique<storm::storage::SparseMatrix<ValueType>>(std::move(A)); |
|||
this->A = localA.get(); |
|||
clearCache(); |
|||
} |
|||
|
|||
template<typename ValueType> |
|||
uint64_t AcyclicLinearEquationSolver<ValueType>::getMatrixRowCount() const { |
|||
return this->A->getRowCount(); |
|||
} |
|||
|
|||
template<typename ValueType> |
|||
uint64_t AcyclicLinearEquationSolver<ValueType>::getMatrixColumnCount() const { |
|||
return this->A->getColumnCount(); |
|||
} |
|||
|
|||
template<typename ValueType> |
|||
bool AcyclicLinearEquationSolver<ValueType>::internalSolveEquations(Environment const& env, std::vector<ValueType>& x, std::vector<ValueType> const& b) const { |
|||
STORM_LOG_ASSERT(x.size() == this->A->getRowGroupCount(), "Provided x-vector has invalid size."); |
|||
STORM_LOG_ASSERT(b.size() == this->A->getRowCount(), "Provided b-vector has invalid size."); |
|||
|
|||
if (!multiplier) { |
|||
// We have not allocated cache memory, yet
|
|||
rowOrdering = helper::computeTopologicalGroupOrdering(*this->A); |
|||
if (!rowOrdering) { |
|||
// It is not required to reorder the elements.
|
|||
this->multiplier = storm::solver::MultiplierFactory<ValueType>().create(env, *this->A); |
|||
} else { |
|||
bFactors.clear(); |
|||
orderedMatrix = helper::createReorderedMatrix(*this->A, *rowOrdering, bFactors); |
|||
this->multiplier = storm::solver::MultiplierFactory<ValueType>().create(env, *orderedMatrix); |
|||
} |
|||
auxiliaryRowVector = std::vector<ValueType>(); |
|||
} |
|||
|
|||
std::vector<ValueType> const* bPtr = &b; |
|||
if (rowOrdering) { |
|||
STORM_LOG_ASSERT(rowOrdering->size() == b.size(), "b-vector has unexpected size."); |
|||
auxiliaryRowVector->resize(b.size()); |
|||
storm::utility::vector::selectVectorValues(*auxiliaryRowVector, *rowOrdering, b); |
|||
for (auto const& bFactor : bFactors) { |
|||
(*auxiliaryRowVector)[bFactor.first] *= bFactor.second; |
|||
} |
|||
bPtr = &auxiliaryRowVector.get(); |
|||
} |
|||
|
|||
this->multiplier->multiplyGaussSeidel(env, x, bPtr, true); |
|||
|
|||
if (!this->isCachingEnabled()) { |
|||
this->clearCache(); |
|||
} |
|||
return true; |
|||
} |
|||
|
|||
template<typename ValueType> |
|||
LinearEquationSolverProblemFormat AcyclicLinearEquationSolver<ValueType>::getEquationProblemFormat(storm::Environment const& env) const { |
|||
return LinearEquationSolverProblemFormat::FixedPointSystem; |
|||
} |
|||
|
|||
template<typename ValueType> |
|||
LinearEquationSolverRequirements AcyclicLinearEquationSolver<ValueType>::getRequirements(Environment const& env) const { |
|||
// Return the requirements of the underlying solver
|
|||
LinearEquationSolverRequirements requirements; |
|||
requirements.requireAcyclic(); |
|||
return requirements; |
|||
} |
|||
|
|||
template<typename ValueType> |
|||
void AcyclicLinearEquationSolver<ValueType>::clearCache() const { |
|||
multiplier.reset(); |
|||
orderedMatrix = boost::none; |
|||
rowOrdering = boost::none; |
|||
auxiliaryRowVector = boost::none; |
|||
bFactors.clear(); |
|||
} |
|||
|
|||
// Explicitly instantiate the min max linear equation solver.
|
|||
template class AcyclicLinearEquationSolver<double>; |
|||
|
|||
#ifdef STORM_HAVE_CARL
|
|||
template class AcyclicLinearEquationSolver<storm::RationalNumber>; |
|||
template class AcyclicLinearEquationSolver<storm::RationalFunction>; |
|||
#endif
|
|||
} |
|||
} |
@ -0,0 +1,100 @@ |
|||
#include "storm/solver/AcyclicMinMaxLinearEquationSolver.h"
|
|||
|
|||
#include "storm/solver/helper/AcyclicSolverHelper.cpp"
|
|||
|
|||
#include "storm/utility/vector.h"
|
|||
|
|||
namespace storm { |
|||
namespace solver { |
|||
|
|||
template<typename ValueType> |
|||
AcyclicMinMaxLinearEquationSolver<ValueType>::AcyclicMinMaxLinearEquationSolver() { |
|||
// Intentionally left empty.
|
|||
} |
|||
|
|||
template<typename ValueType> |
|||
AcyclicMinMaxLinearEquationSolver<ValueType>::AcyclicMinMaxLinearEquationSolver(storm::storage::SparseMatrix<ValueType> const& A) : StandardMinMaxLinearEquationSolver<ValueType>(A) { |
|||
// Intentionally left empty.
|
|||
} |
|||
|
|||
template<typename ValueType> |
|||
AcyclicMinMaxLinearEquationSolver<ValueType>::AcyclicMinMaxLinearEquationSolver(storm::storage::SparseMatrix<ValueType>&& A) : StandardMinMaxLinearEquationSolver<ValueType>(std::move(A)) { |
|||
// Intentionally left empty.
|
|||
} |
|||
|
|||
|
|||
template<typename ValueType> |
|||
bool AcyclicMinMaxLinearEquationSolver<ValueType>::internalSolveEquations(Environment const& env, OptimizationDirection dir, std::vector<ValueType>& x, std::vector<ValueType> const& b) const { |
|||
STORM_LOG_ASSERT(x.size() == this->A->getRowGroupCount(), "Provided x-vector has invalid size."); |
|||
STORM_LOG_ASSERT(b.size() == this->A->getRowCount(), "Provided b-vector has invalid size."); |
|||
// Allocate memory for the scheduler (if required)
|
|||
if (this->isTrackSchedulerSet()) { |
|||
if (this->schedulerChoices) { |
|||
this->schedulerChoices->resize(this->A->getRowGroupCount()); |
|||
} else { |
|||
this->schedulerChoices = std::vector<uint_fast64_t>(this->A->getRowGroupCount()); |
|||
} |
|||
} |
|||
|
|||
if (!multiplier) { |
|||
// We have not allocated cache memory, yet
|
|||
rowGroupOrdering = helper::computeTopologicalGroupOrdering(*this->A); |
|||
if (!rowGroupOrdering) { |
|||
// It is not required to reorder the elements.
|
|||
this->multiplier = storm::solver::MultiplierFactory<ValueType>().create(env, *this->A); |
|||
} else { |
|||
bFactors.clear(); |
|||
orderedMatrix = helper::createReorderedMatrix(*this->A, *rowGroupOrdering, bFactors); |
|||
this->multiplier = storm::solver::MultiplierFactory<ValueType>().create(env, *orderedMatrix); |
|||
} |
|||
auxiliaryRowVector = std::vector<ValueType>(); |
|||
} |
|||
|
|||
std::vector<ValueType> const* bPtr = &b; |
|||
if (rowGroupOrdering) { |
|||
STORM_LOG_ASSERT(rowGroupOrdering->size() == b.size(), "b-vector has unexpected size."); |
|||
auxiliaryRowVector->resize(b.size()); |
|||
storm::utility::vector::selectVectorValues(*auxiliaryRowVector, *rowGroupOrdering, b); |
|||
for (auto const& bFactor : bFactors) { |
|||
(*auxiliaryRowVector)[bFactor.first] *= bFactor.second; |
|||
} |
|||
bPtr = &auxiliaryRowVector.get(); |
|||
} |
|||
|
|||
if (this->isTrackSchedulerSet()) { |
|||
this->multiplier->multiplyAndReduceGaussSeidel(env, dir, x, bPtr, &this->schedulerChoices.get(), true); |
|||
} else { |
|||
this->multiplier->multiplyAndReduceGaussSeidel(env, dir, x, bPtr, nullptr, true); |
|||
} |
|||
|
|||
if (!this->isCachingEnabled()) { |
|||
this->clearCache(); |
|||
} |
|||
return true; |
|||
} |
|||
|
|||
template<typename ValueType> |
|||
MinMaxLinearEquationSolverRequirements AcyclicMinMaxLinearEquationSolver<ValueType>::getRequirements(Environment const& env, boost::optional<storm::solver::OptimizationDirection> const& direction, bool const& hasInitialScheduler) const { |
|||
// Return the requirements of the underlying solver
|
|||
MinMaxLinearEquationSolverRequirements requirements; |
|||
requirements.requireAcyclic(); |
|||
return requirements; |
|||
} |
|||
|
|||
template<typename ValueType> |
|||
void AcyclicMinMaxLinearEquationSolver<ValueType>::clearCache() const { |
|||
multiplier.reset(); |
|||
orderedMatrix = boost::none; |
|||
rowGroupOrdering = boost::none; |
|||
auxiliaryRowVector = boost::none; |
|||
bFactors.clear(); |
|||
} |
|||
|
|||
// Explicitly instantiate the min max linear equation solver.
|
|||
template class AcyclicMinMaxLinearEquationSolver<double>; |
|||
|
|||
#ifdef STORM_HAVE_CARL
|
|||
template class AcyclicMinMaxLinearEquationSolver<storm::RationalNumber>; |
|||
#endif
|
|||
} |
|||
} |
@ -0,0 +1,50 @@ |
|||
#pragma once |
|||
|
|||
#include <memory> |
|||
|
|||
#include "storm/solver/Multiplier.h" |
|||
#include "storm/solver/StandardMinMaxLinearEquationSolver.h" |
|||
|
|||
namespace storm { |
|||
|
|||
class Environment; |
|||
|
|||
namespace solver { |
|||
|
|||
/*! |
|||
* This solver can be used on equation systems that are known to be acyclic. |
|||
* It is optimized for solving many instances of the equation system with the same underlying matrix. |
|||
*/ |
|||
template<typename ValueType> |
|||
class AcyclicMinMaxLinearEquationSolver : public StandardMinMaxLinearEquationSolver<ValueType> { |
|||
public: |
|||
AcyclicMinMaxLinearEquationSolver(); |
|||
AcyclicMinMaxLinearEquationSolver(storm::storage::SparseMatrix<ValueType> const& A); |
|||
AcyclicMinMaxLinearEquationSolver(storm::storage::SparseMatrix<ValueType>&& A); |
|||
|
|||
virtual ~AcyclicMinMaxLinearEquationSolver() { |
|||
} |
|||
|
|||
virtual void clearCache() const override; |
|||
|
|||
virtual MinMaxLinearEquationSolverRequirements getRequirements(Environment const& env, boost::optional<storm::solver::OptimizationDirection> const& direction = boost::none, bool const& hasInitialScheduler = false) const override ; |
|||
|
|||
protected: |
|||
|
|||
virtual bool internalSolveEquations(storm::Environment const& env, OptimizationDirection d, std::vector<ValueType>& x, std::vector<ValueType> const& b) const override; |
|||
|
|||
private: |
|||
// cached multiplier either with original matrix or ordered matrix |
|||
mutable std::unique_ptr<storm::solver::Multiplier<ValueType>> multiplier; |
|||
// cached matrix for the multiplier (only if different from original matrix) |
|||
mutable boost::optional<storm::storage::SparseMatrix<ValueType>> orderedMatrix; |
|||
// cached row group ordering (only if not identity) |
|||
mutable boost::optional<std::vector<uint64_t>> rowGroupOrdering; // A.rowGroupCount() entries |
|||
// can be used if the entries in 'b' need to be reordered |
|||
mutable boost::optional<std::vector<ValueType>> auxiliaryRowVector; // A.rowCount() entries |
|||
// contains factors applied to scale the entries of the 'b' vector |
|||
mutable std::vector<std::pair<uint64_t, ValueType>> bFactors; |
|||
|
|||
}; |
|||
} |
|||
} |
@ -0,0 +1,117 @@ |
|||
#include "storm/storage/SparseMatrix.h"
|
|||
#include "storm/utility/constants.h"
|
|||
#include "storm/utility/vector.h"
|
|||
#include "storm/exceptions/InvalidStateException.h"
|
|||
#include "storm/exceptions/InvalidEnvironmentException.h"
|
|||
#include "storm/exceptions/UnexpectedException.h"
|
|||
#include "storm/exceptions/UnmetRequirementException.h"
|
|||
|
|||
namespace storm { |
|||
namespace solver { |
|||
|
|||
namespace helper { |
|||
|
|||
|
|||
/*!
|
|||
* Returns a reordering of the matrix row(groups) and columns such that we can solve the (minmax or linear) equation system in one go. |
|||
* More precisely, let x be the result and i an arbitrary rowgroup index. Solving for rowgroup x[i] only requires knowledge of the result at rowgroups x[i+1], x[i+2], ... |
|||
*/ |
|||
template<typename ValueType> |
|||
boost::optional<std::vector<uint64_t>> computeTopologicalGroupOrdering(storm::storage::SparseMatrix<ValueType> const& matrix) { |
|||
uint64_t numGroups = matrix.getRowGroupCount(); |
|||
bool orderedMatrixRequired = false; |
|||
std::vector<uint64_t> result; |
|||
result.reserve(numGroups); |
|||
storm::storage::BitVector processed(numGroups, false); |
|||
storm::storage::BitVector visited(numGroups, false); |
|||
std::vector<uint64_t> stack; |
|||
|
|||
// It's more likely that states without a successor are located at the end (due to the way we build the model).
|
|||
// We therefore process the matrix in reverse order.
|
|||
uint64_t startState = numGroups; |
|||
while (startState > 0) { |
|||
--startState; |
|||
// skip already processed states
|
|||
if (processed.get(startState)) continue; |
|||
|
|||
// Now do a dfs from start state.
|
|||
stack.push_back(startState); |
|||
while (!stack.empty()) { |
|||
uint64_t current = stack.back(); |
|||
if (visited.get(current)) { |
|||
// We are backtracking, so add this state now
|
|||
result.push_back(current); |
|||
processed.set(current); |
|||
stack.pop_back(); |
|||
} else { |
|||
visited.set(current); |
|||
for (auto const& entry : matrix.getRowGroup(current)) { |
|||
if (!processed.get(entry.getColumn()) && !storm::utility::isZero(entry.getValue())) { |
|||
orderedMatrixRequired = true; |
|||
STORM_LOG_THROW(!visited.get(entry.getColumn()), storm::exceptions::UnmetRequirementException, "The model is not acyclic."); |
|||
stack.push_back(entry.getColumn()); |
|||
} |
|||
} |
|||
// If there are no successors to process, we will add the current state to the result in the next iteration.
|
|||
} |
|||
} |
|||
} |
|||
// we will do backwards iterations, so the order has to be reversed
|
|||
if (orderedMatrixRequired) { |
|||
std::reverse(result.begin(), result.end()); |
|||
return result; |
|||
} else { |
|||
return boost::none; |
|||
} |
|||
} |
|||
|
|||
/// reorders the row group such that the i'th row of the new matrix corresponds to the order[i]'th row of the source matrix.
|
|||
/// Also eliminates selfloops p>0 and inserts 1/p into the bFactors
|
|||
template<typename ValueType> |
|||
storm::storage::SparseMatrix<ValueType> createReorderedMatrix(storm::storage::SparseMatrix<ValueType> const& matrix, std::vector<uint64_t> const& newToOrigIndexMap, std::vector<std::pair<uint64_t, ValueType>>& bFactors) { |
|||
std::vector<uint64_t> origToNewMap(newToOrigIndexMap.size(), std::numeric_limits<uint64>::max()); |
|||
for (uint64_t i = 0; i < newToOrigIndexMap.size(); ++i) { |
|||
origToNewMap[newToOrigIndexMap[i]] = i; |
|||
} |
|||
|
|||
bool hasRowGrouping = !matrix.hasTrivialRowGrouping(); |
|||
storm::storage::SparseMatrixBuilder<ValueType> builder(matrix.getRowCount(), matrix.getColumnCount(), matrix.getNonzeroEntryCount(), false, hasRowGrouping, hasRowGrouping ? matrix.getRowGroupCount() : static_cast<uint64_t>(0)); |
|||
uint64_t newRow = 0; |
|||
for (uint64_t newRowGroup = 0; newRowGroup < newToOrigIndexMap.size(); ++newRowGroup) { |
|||
auto const& origRowGroup = newToOrigIndexMap[newRowGroup]; |
|||
if (hasRowGrouping) { |
|||
builder.newRowGroup(newRowGroup); |
|||
} |
|||
for (uint64_t origRow = matrix.getRowGroupIndices()[origRowGroup]; origRow < matrix.getRowGroupIndices()[origRowGroup + 1]; ++origRow) { |
|||
for (auto const& entry : matrix.getRow(origRow)) { |
|||
if (storm::utility::isZero(entry.getValue())) { |
|||
continue; |
|||
} |
|||
if (entry.getColumn() == origRowGroup) { |
|||
if (storm::utility::isOne(entry.getValue())) { |
|||
// a one selfloop can only mean that there is never a non-zero value at the b vector for the current row.
|
|||
// Let's "assert" this by considering infinity. (This is necessary to avoid division by zero)
|
|||
bFactors.emplace_back(newRow, storm::utility::infinity<ValueType>()); |
|||
} |
|||
ValueType factor = storm::utility::one<ValueType>() / (storm::utility::one<ValueType>() - entry.getValue()); |
|||
bFactors.emplace_back(newRow, factor); |
|||
} |
|||
builder.addNextValue(newRow, origToNewMap[entry.getColumn()], entry.getValue()); |
|||
} |
|||
++newRow; |
|||
} |
|||
} |
|||
auto result = builder.build(matrix.getRowCount(), matrix.getColumnCount(), matrix.getRowGroupCount()); |
|||
// apply the bFactors to the relevant rows
|
|||
for (auto const& bFactor : bFactors) { |
|||
STORM_LOG_ASSERT(!storm::utility::isInfinity(bFactor.second) || storm::utility::isZero(result.getRowSum(bFactor.first)), "The input matrix does not seem to be probabilistic."); |
|||
for (auto& entry : result.getRow(bFactor.first)) { |
|||
entry.setValue(entry.getValue() * bFactor.second); |
|||
} |
|||
} |
|||
STORM_LOG_DEBUG("Reordered " << matrix.getDimensionsAsString() << " with " << bFactors.size() << " selfloop entries for acyclic solving."); |
|||
return result; |
|||
} |
|||
} |
|||
} |
|||
} |
Write
Preview
Loading…
Cancel
Save
Reference in new issue