Browse Source

improved interface of solver environment

tempestpy_adaptions
TimQu 7 years ago
parent
commit
1cff0fcbbb
  1. 58
      src/storm/environment/solver/SolverEnvironment.cpp
  2. 9
      src/storm/environment/solver/SolverEnvironment.h
  3. 52
      src/storm/solver/IterativeMinMaxLinearEquationSolver.cpp
  4. 27
      src/storm/solver/StandardGameSolver.cpp
  5. 52
      src/storm/solver/SymbolicMinMaxLinearEquationSolver.cpp

58
src/storm/environment/solver/SolverEnvironment.cpp

@ -13,6 +13,7 @@
#include "storm/utility/macros.h"
#include "storm/exceptions/InvalidEnvironmentException.h"
#include "storm/exceptions/UnexpectedException.h"
namespace storm {
@ -87,8 +88,8 @@ namespace storm {
return linearEquationSolverType;
}
void SolverEnvironment::setLinearEquationSolverType(storm::solver::EquationSolverType const& value) {
linearEquationSolverTypeSetFromDefault = false;
void SolverEnvironment::setLinearEquationSolverType(storm::solver::EquationSolverType const& value, bool assumeSetFromDefault) {
linearEquationSolverTypeSetFromDefault = assumeSetFromDefault;
linearEquationSolverType = value;
}
@ -96,46 +97,49 @@ namespace storm {
return linearEquationSolverTypeSetFromDefault;
}
boost::optional<storm::RationalNumber> SolverEnvironment::getPrecisionOfCurrentLinearEquationSolver() const {
switch (getLinearEquationSolverType()) {
std::pair<boost::optional<storm::RationalNumber>, boost::optional<bool>> SolverEnvironment::getPrecisionOfLinearEquationSolver(storm::solver::EquationSolverType const& solverType) const {
std::pair<boost::optional<storm::RationalNumber>, boost::optional<bool>> result;
switch (solverType) {
case storm::solver::EquationSolverType::Gmmxx:
return gmmxx().getPrecision();
result.first = gmmxx().getPrecision();
break;
case storm::solver::EquationSolverType::Eigen:
return eigen().getPrecision();
result.first = eigen().getPrecision();
break;
case storm::solver::EquationSolverType::Native:
return native().getPrecision();
result.first = native().getPrecision();
result.second = native().getRelativeTerminationCriterion();
break;
case storm::solver::EquationSolverType::Elimination:
return boost::none;
break;
case storm::solver::EquationSolverType::Topological:
result = getPrecisionOfLinearEquationSolver(topological().getUnderlyingSolverType());
break;
default:
STORM_LOG_ASSERT(false, "The selected solver type is unknown.");
STORM_LOG_THROW(false, storm::exceptions::UnexpectedException, "The selected solver type is unknown.");
}
return result;
}
void SolverEnvironment::setLinearEquationSolverPrecision(storm::RationalNumber const& value) {
void SolverEnvironment::setLinearEquationSolverPrecision(boost::optional<storm::RationalNumber> const& newPrecision, boost::optional<bool> const& relativePrecision) {
// Assert that each solver type is handled in this method.
STORM_LOG_ASSERT(getLinearEquationSolverType() == storm::solver::EquationSolverType::Native ||
getLinearEquationSolverType() == storm::solver::EquationSolverType::Gmmxx ||
getLinearEquationSolverType() == storm::solver::EquationSolverType::Eigen ||
getLinearEquationSolverType() == storm::solver::EquationSolverType::Elimination ||
getLinearEquationSolverType() == storm::solver::EquationSolverType::Topological,
"The current solver type is not respected in this method.");
native().setPrecision(value);
gmmxx().setPrecision(value);
eigen().setPrecision(value);
// Elimination and Topological solver do not have a precision
}
void SolverEnvironment::setLinearEquationSolverRelativeTerminationCriterion(bool value) {
STORM_LOG_ASSERT(getLinearEquationSolverType() == storm::solver::EquationSolverType::Native ||
getLinearEquationSolverType() == storm::solver::EquationSolverType::Gmmxx ||
getLinearEquationSolverType() == storm::solver::EquationSolverType::Eigen ||
getLinearEquationSolverType() == storm::solver::EquationSolverType::Elimination ||
getLinearEquationSolverType() == storm::solver::EquationSolverType::Topological,
"The current solver type is not respected in this method.");
native().setRelativeTerminationCriterion(value);
// Elimination, gmm, eigen, and topological solver do not have an option for relative termination criterion
if (newPrecision) {
native().setPrecision(newPrecision.get());
gmmxx().setPrecision(newPrecision.get());
eigen().setPrecision(newPrecision.get());
// Elimination and Topological solver do not have a precision
}
if (relativePrecision) {
native().setRelativeTerminationCriterion(relativePrecision.get());
// gmm, eigen, elimination, and topological solvers do not have a precision
}
}
}

9
src/storm/environment/solver/SolverEnvironment.h

@ -41,12 +41,11 @@ namespace storm {
void setForceSoundness(bool value);
storm::solver::EquationSolverType const& getLinearEquationSolverType() const;
void setLinearEquationSolverType(storm::solver::EquationSolverType const& value);
void setLinearEquationSolverType(storm::solver::EquationSolverType const& value, bool assumeSetFromDefault = false);
bool isLinearEquationSolverTypeSetFromDefaultValue() const;
boost::optional<storm::RationalNumber> getPrecisionOfCurrentLinearEquationSolver() const;
void setLinearEquationSolverPrecision(storm::RationalNumber const& value);
void setLinearEquationSolverRelativeTerminationCriterion(bool value);
std::pair<boost::optional<storm::RationalNumber>, boost::optional<bool>> getPrecisionOfLinearEquationSolver(storm::solver::EquationSolverType const& solverType) const;
void setLinearEquationSolverPrecision(boost::optional<storm::RationalNumber> const& newPrecision, boost::optional<bool> const& relativePrecision = boost::none);
private:
SubEnvironment<EigenSolverEnvironment> eigenSolverEnvironment;

52
src/storm/solver/IterativeMinMaxLinearEquationSolver.cpp

@ -114,19 +114,32 @@ namespace storm {
// The solver that we will use throughout the procedure.
std::unique_ptr<storm::solver::LinearEquationSolver<ValueType>> solver;
// The linear equation solver should be at least as precise as this solver
std::unique_ptr<storm::Environment> environmentOfSolver;
boost::optional<storm::RationalNumber> precOfSolver = env.solver().getPrecisionOfCurrentLinearEquationSolver();
if (!storm::NumberTraits<ValueType>::IsExact && precOfSolver && precOfSolver.get() > env.solver().minMax().getPrecision()) {
environmentOfSolver = std::make_unique<storm::Environment>(env);
environmentOfSolver->solver().setLinearEquationSolverPrecision(env.solver().minMax().getPrecision());
std::unique_ptr<storm::Environment> environmentOfSolverStorage;
auto precOfSolver = env.solver().getPrecisionOfLinearEquationSolver(env.solver().getLinearEquationSolverType());
if (!storm::NumberTraits<ValueType>::IsExact) {
bool changePrecision = precOfSolver.first && precOfSolver.first.get() > env.solver().minMax().getPrecision();
bool changeRelative = precOfSolver.second && !precOfSolver.second.get() && env.solver().minMax().getRelativeTerminationCriterion();
if (changePrecision || changeRelative) {
environmentOfSolverStorage = std::make_unique<storm::Environment>(env);
boost::optional<storm::RationalNumber> newPrecision;
boost::optional<bool> newRelative;
if (changePrecision) {
newPrecision = env.solver().minMax().getPrecision();
}
if (changeRelative) {
newRelative = true;
}
environmentOfSolverStorage->solver().setLinearEquationSolverPrecision(newPrecision, newRelative);
}
}
storm::Environment const& environmentOfSolver = environmentOfSolverStorage ? *environmentOfSolverStorage : env;
SolverStatus status = SolverStatus::InProgress;
uint64_t iterations = 0;
this->startMeasureProgress();
do {
// Solve the equation system for the 'DTMC'.
solveInducedEquationSystem(environmentOfSolver ? *environmentOfSolver : env, solver, scheduler, x, subB, b);
solveInducedEquationSystem(environmentOfSolver, solver, scheduler, x, subB, b);
// Go through the multiplication result and see whether we can improve any of the choices.
bool schedulerImproved = false;
@ -311,14 +324,27 @@ namespace storm {
// Solve the equation system induced by the initial scheduler.
std::unique_ptr<storm::solver::LinearEquationSolver<ValueType>> linEqSolver;
// The linear equation solver should be at least as precise as this solver
std::unique_ptr<storm::Environment> environmentOfSolver;
boost::optional<storm::RationalNumber> precOfSolver = env.solver().getPrecisionOfCurrentLinearEquationSolver();
if (!storm::NumberTraits<ValueType>::IsExact && precOfSolver && precOfSolver.get() > env.solver().minMax().getPrecision()) {
environmentOfSolver = std::make_unique<storm::Environment>(env);
environmentOfSolver->solver().setLinearEquationSolverPrecision(env.solver().minMax().getPrecision());
std::unique_ptr<storm::Environment> environmentOfSolverStorage;
auto precOfSolver = env.solver().getPrecisionOfLinearEquationSolver(env.solver().getLinearEquationSolverType());
if (!storm::NumberTraits<ValueType>::IsExact) {
bool changePrecision = precOfSolver.first && precOfSolver.first.get() > env.solver().minMax().getPrecision();
bool changeRelative = precOfSolver.second && !precOfSolver.second.get() && env.solver().minMax().getRelativeTerminationCriterion();
if (changePrecision || changeRelative) {
environmentOfSolverStorage = std::make_unique<storm::Environment>(env);
boost::optional<storm::RationalNumber> newPrecision;
boost::optional<bool> newRelative;
if (changePrecision) {
newPrecision = env.solver().minMax().getPrecision();
}
if (changeRelative) {
newRelative = true;
}
environmentOfSolverStorage->solver().setLinearEquationSolverPrecision(newPrecision, newRelative);
}
}
storm::Environment const& environmentOfSolver = environmentOfSolverStorage ? *environmentOfSolverStorage : env;
solveInducedEquationSystem(environmentOfSolver ? *environmentOfSolver : env, linEqSolver, this->getInitialScheduler(), x, *auxiliaryRowGroupVector, b);
solveInducedEquationSystem(environmentOfSolver, linEqSolver, this->getInitialScheduler(), x, *auxiliaryRowGroupVector, b);
// If we were given an initial scheduler and are maximizing (minimizing), our current solution becomes
// always less-or-equal (greater-or-equal) than the actual solution.
guarantee = maximize(dir) ? SolverGuarantee::LessOrEqual : SolverGuarantee::GreaterOrEqual;

27
src/storm/solver/StandardGameSolver.cpp

@ -76,13 +76,34 @@ namespace storm {
uint64_t maxIter = env.solver().game().getMaximalNumberOfIterations();
// The linear equation solver should be at least as precise as this solver
std::unique_ptr<storm::Environment> environmentOfSolverStorage;
auto precOfSolver = env.solver().getPrecisionOfLinearEquationSolver(env.solver().getLinearEquationSolverType());
if (!storm::NumberTraits<ValueType>::IsExact) {
bool changePrecision = precOfSolver.first && precOfSolver.first.get() > env.solver().game().getPrecision();
bool changeRelative = precOfSolver.second && !precOfSolver.second.get() && env.solver().game().getRelativeTerminationCriterion();
if (changePrecision || changeRelative) {
environmentOfSolverStorage = std::make_unique<storm::Environment>(env);
boost::optional<storm::RationalNumber> newPrecision;
boost::optional<bool> newRelative;
if (changePrecision) {
newPrecision = env.solver().game().getPrecision();
}
if (changeRelative) {
newRelative = true;
}
environmentOfSolverStorage->solver().setLinearEquationSolverPrecision(newPrecision, newRelative);
}
}
storm::Environment const& environmentOfSolver = environmentOfSolverStorage ? *environmentOfSolverStorage : env;
// Solve the equation system induced by the two schedulers.
storm::storage::SparseMatrix<ValueType> submatrix;
getInducedMatrixVector(x, b, player1Choices, player2Choices, submatrix, subB);
if (this->linearEquationSolverFactory->getEquationProblemFormat(env) == LinearEquationSolverProblemFormat::EquationSystem) {
if (this->linearEquationSolverFactory->getEquationProblemFormat(environmentOfSolver) == LinearEquationSolverProblemFormat::EquationSystem) {
submatrix.convertToEquationSystem();
}
auto submatrixSolver = linearEquationSolverFactory->create(env, std::move(submatrix));
auto submatrixSolver = linearEquationSolverFactory->create(environmentOfSolver, std::move(submatrix));
if (this->lowerBound) { submatrixSolver->setLowerBound(this->lowerBound.get()); }
if (this->upperBound) { submatrixSolver->setUpperBound(this->upperBound.get()); }
submatrixSolver->setCachingEnabled(true);
@ -92,7 +113,7 @@ namespace storm {
do {
// Solve the equation system for the 'DTMC'.
// FIXME: we need to remove the 0- and 1- states to make the solution unique.
submatrixSolver->solveEquations(env, x, subB);
submatrixSolver->solveEquations(environmentOfSolver, x, subB);
bool schedulerImproved = extractChoices(player1Dir, player2Dir, x, b, *auxiliaryP2RowGroupVector, player1Choices, player2Choices);

52
src/storm/solver/SymbolicMinMaxLinearEquationSolver.cpp

@ -244,13 +244,26 @@ namespace storm {
// If we were given an initial scheduler, we take its solution as the starting point.
if (this->hasInitialScheduler()) {
// The linear equation solver should be at least as precise as this solver
std::unique_ptr<storm::Environment> environmentOfSolver;
boost::optional<storm::RationalNumber> precOfSolver = env.solver().getPrecisionOfCurrentLinearEquationSolver();
if (!storm::NumberTraits<ValueType>::IsExact && precOfSolver && precOfSolver.get() > env.solver().minMax().getPrecision()) {
environmentOfSolver = std::make_unique<storm::Environment>(env);
environmentOfSolver->solver().setLinearEquationSolverPrecision(env.solver().minMax().getPrecision());
std::unique_ptr<storm::Environment> environmentOfSolverStorage;
auto precOfSolver = env.solver().getPrecisionOfLinearEquationSolver(env.solver().getLinearEquationSolverType());
if (!storm::NumberTraits<ValueType>::IsExact) {
bool changePrecision = precOfSolver.first && precOfSolver.first.get() > env.solver().minMax().getPrecision();
bool changeRelative = precOfSolver.second && !precOfSolver.second.get() && env.solver().minMax().getRelativeTerminationCriterion();
if (changePrecision || changeRelative) {
environmentOfSolverStorage = std::make_unique<storm::Environment>(env);
boost::optional<storm::RationalNumber> newPrecision;
boost::optional<bool> newRelative;
if (changePrecision) {
newPrecision = env.solver().minMax().getPrecision();
}
if (changeRelative) {
newRelative = true;
}
environmentOfSolverStorage->solver().setLinearEquationSolverPrecision(newPrecision, newRelative);
}
}
localX = solveEquationsWithScheduler(environmentOfSolver ? *environmentOfSolver : env, this->getInitialScheduler(), x, b);
storm::Environment const& environmentOfSolver = environmentOfSolverStorage ? *environmentOfSolverStorage : env;
localX = solveEquationsWithScheduler(environmentOfSolver, this->getInitialScheduler(), x, b);
} else {
localX = this->getLowerBoundsVector();
}
@ -311,20 +324,33 @@ namespace storm {
// Initialize linear equation solver.
// It should be at least as precise as this solver.
std::unique_ptr<storm::Environment> environmentOfSolver;
boost::optional<storm::RationalNumber> precOfSolver = env.solver().getPrecisionOfCurrentLinearEquationSolver();
if (!storm::NumberTraits<ValueType>::IsExact && precOfSolver && precOfSolver.get() > env.solver().minMax().getPrecision()) {
environmentOfSolver = std::make_unique<storm::Environment>(env);
environmentOfSolver->solver().setLinearEquationSolverPrecision(env.solver().minMax().getPrecision());
std::unique_ptr<storm::Environment> environmentOfSolverStorage;
auto precOfSolver = env.solver().getPrecisionOfLinearEquationSolver(env.solver().getLinearEquationSolverType());
if (!storm::NumberTraits<ValueType>::IsExact) {
bool changePrecision = precOfSolver.first && precOfSolver.first.get() > env.solver().minMax().getPrecision();
bool changeRelative = precOfSolver.second && !precOfSolver.second.get() && env.solver().minMax().getRelativeTerminationCriterion();
if (changePrecision || changeRelative) {
environmentOfSolverStorage = std::make_unique<storm::Environment>(env);
boost::optional<storm::RationalNumber> newPrecision;
boost::optional<bool> newRelative;
if (changePrecision) {
newPrecision = env.solver().minMax().getPrecision();
}
if (changeRelative) {
newRelative = true;
}
environmentOfSolverStorage->solver().setLinearEquationSolverPrecision(newPrecision, newRelative);
}
}
storm::Environment const& environmentOfSolver = environmentOfSolverStorage ? *environmentOfSolverStorage : env;
std::unique_ptr<SymbolicLinearEquationSolver<DdType, ValueType>> linearEquationSolver = linearEquationSolverFactory->create(environmentOfSolver ? *environmentOfSolver : env, this->allRows, this->rowMetaVariables, this->columnMetaVariables, this->rowColumnMetaVariablePairs);
std::unique_ptr<SymbolicLinearEquationSolver<DdType, ValueType>> linearEquationSolver = linearEquationSolverFactory->create(environmentOfSolver, this->allRows, this->rowMetaVariables, this->columnMetaVariables, this->rowColumnMetaVariablePairs);
this->forwardBounds(*linearEquationSolver);
// Iteratively solve and improve the scheduler.
uint64_t maxIter = env.solver().minMax().getMaximalNumberOfIterations();
while (!converged && iterations < maxIter) {
storm::dd::Add<DdType, ValueType> schedulerX = solveEquationsWithScheduler(environmentOfSolver ? *environmentOfSolver : env, *linearEquationSolver, scheduler, currentSolution, b, diagonal);
storm::dd::Add<DdType, ValueType> schedulerX = solveEquationsWithScheduler(environmentOfSolver, *linearEquationSolver, scheduler, currentSolution, b, diagonal);
// Policy improvement step.
storm::dd::Add<DdType, ValueType> choiceValues = this->A.multiplyMatrix(schedulerX.swapVariables(this->rowColumnMetaVariablePairs), this->columnMetaVariables) + b;

Loading…
Cancel
Save