Browse Source

Caching of solvers can now be enabled/disabled

Former-commit-id: 498622f45c [formerly 1713f554a0]
Former-commit-id: 1015c7fef8
tempestpy_adaptions
TimQu 8 years ago
parent
commit
53ff482947
  1. 7
      src/modelchecker/csl/helper/SparseCtmcCslHelper.cpp
  2. 1
      src/modelchecker/csl/helper/SparseMarkovAutomatonCslHelper.cpp
  3. 2
      src/modelchecker/multiobjective/pcaa/SparseMaPcaaWeightVectorChecker.cpp
  4. 4
      src/solver/EigenLinearEquationSolver.cpp
  5. 4
      src/solver/EliminationLinearEquationSolver.cpp
  6. 28
      src/solver/GmmxxLinearEquationSolver.cpp
  7. 7
      src/solver/GmmxxLinearEquationSolver.h
  8. 42
      src/solver/LinearEquationSolver.cpp
  9. 22
      src/solver/LinearEquationSolver.h
  10. 18
      src/solver/MinMaxLinearEquationSolver.cpp
  11. 20
      src/solver/MinMaxLinearEquationSolver.h
  12. 43
      src/solver/NativeLinearEquationSolver.cpp
  13. 6
      src/solver/NativeLinearEquationSolver.h
  14. 19
      src/solver/StandardMinMaxLinearEquationSolver.cpp
  15. 3
      src/solver/StandardMinMaxLinearEquationSolver.h
  16. 2
      src/utility/policyguessing.cpp

7
src/modelchecker/csl/helper/SparseCtmcCslHelper.cpp

@ -629,6 +629,7 @@ namespace storm {
}
std::unique_ptr<storm::solver::LinearEquationSolver<ValueType>> solver = linearEquationSolverFactory.create(std::move(uniformizedMatrix));
solver->setCachingEnabled(true);
if (!useMixedPoissonProbabilities && std::get<0>(foxGlynnResult) > 1) {
// Perform the matrix-vector multiplications (without adding).
@ -735,12 +736,6 @@ namespace storm {
template storm::storage::SparseMatrix<storm::RationalNumber> SparseCtmcCslHelper::computeProbabilityMatrix(storm::storage::SparseMatrix<storm::RationalNumber> const& rateMatrix, std::vector<storm::RationalNumber> const& exitRates);
template storm::storage::SparseMatrix<storm::RationalFunction> SparseCtmcCslHelper::computeProbabilityMatrix(storm::storage::SparseMatrix<storm::RationalFunction> const& rateMatrix, std::vector<storm::RationalFunction> const& exitRates);
#endif
}
}

1
src/modelchecker/csl/helper/SparseMarkovAutomatonCslHelper.cpp

@ -86,6 +86,7 @@ namespace storm {
}
std::unique_ptr<storm::solver::MinMaxLinearEquationSolver<ValueType>> solver = minMaxLinearEquationSolverFactory.create(aProbabilistic);
solver->setCachingEnabled(true);
// Perform the actual value iteration
// * loop until the step bound has been reached

2
src/modelchecker/multiobjective/pcaa/SparseMaPcaaWeightVectorChecker.cpp

@ -307,6 +307,7 @@ namespace storm {
result->solver = minMaxSolverFactory.create(PS.toPS);
result->solver->setOptimizationDirection(storm::solver::OptimizationDirection::Maximize);
result->solver->setTrackScheduler(true);
result->solver->setCachingEnabled(true);
result->b.resize(PS.getNumberOfChoices());
@ -365,6 +366,7 @@ namespace storm {
storm::storage::SparseMatrix<ValueType> linEqMatrix = PS.toPS.selectRowsFromRowGroups(optimalChoicesAtCurrentEpoch, true);
linEqMatrix.convertToEquationSystem();
linEq.solver = linEq.factory.create(std::move(linEqMatrix));
linEq.solver->setCachingEnabled(true);
}
// Get the results for the individual objectives.

4
src/solver/EigenLinearEquationSolver.cpp

@ -117,7 +117,7 @@ namespace storm {
template<typename ValueType>
void EigenLinearEquationSolver<ValueType>::setMatrix(storm::storage::SparseMatrix<ValueType> const& A) {
eigenA = storm::adapters::EigenAdapter::toEigenSparseMatrix<ValueType>(A);
this->resetAuxiliaryData();
this->clearCache();
}
template<typename ValueType>
@ -125,7 +125,7 @@ namespace storm {
// Take ownership of the matrix so it is destroyed after we have translated it to Eigen's format.
storm::storage::SparseMatrix<ValueType> localA(std::move(A));
this->setMatrix(localA);
this->resetAuxiliaryData();
this->clearCache();
}
template<typename ValueType>

4
src/solver/EliminationLinearEquationSolver.cpp

@ -47,14 +47,14 @@ namespace storm {
void EliminationLinearEquationSolver<ValueType>::setMatrix(storm::storage::SparseMatrix<ValueType> const& A) {
this->A = &A;
localA.reset();
this->resetAuxiliaryData();
this->clearCache();
}
template<typename ValueType>
void EliminationLinearEquationSolver<ValueType>::setMatrix(storm::storage::SparseMatrix<ValueType>&& A) {
localA = std::make_unique<storm::storage::SparseMatrix<ValueType>>(std::move(A));
this->A = localA.get();
this->resetAuxiliaryData();
this->clearCache();
}
template<typename ValueType>

28
src/solver/GmmxxLinearEquationSolver.cpp

@ -124,14 +124,14 @@ namespace storm {
void GmmxxLinearEquationSolver<ValueType>::setMatrix(storm::storage::SparseMatrix<ValueType> const& A) {
localA.reset();
this->A = &A;
resetAuxiliaryData();
clearCache();
}
template<typename ValueType>
void GmmxxLinearEquationSolver<ValueType>::setMatrix(storm::storage::SparseMatrix<ValueType>&& A) {
localA = std::make_unique<storm::storage::SparseMatrix<ValueType>>(std::move(A));
this->A = localA.get();
resetAuxiliaryData();
clearCache();
}
template<typename ValueType>
@ -187,6 +187,10 @@ namespace storm {
}
}
if(!this->isCachingEnabled()) {
clearCache();
}
// Check if the solver converged and issue a warning otherwise.
if (iter.converged()) {
STORM_LOG_INFO("Iterative solver converged after " << iter.get_iteration() << " iterations.");
@ -221,6 +225,10 @@ namespace storm {
} else {
gmm::mult(*gmmxxA, x, result);
}
if(!this->isCachingEnabled()) {
clearCache();
}
}
template<typename ValueType>
@ -238,10 +246,10 @@ namespace storm {
std::vector<ValueType>* currentX = &x;
if(!this->auxiliaryRowVector) {
this->auxiliaryRowVector = std::make_unique<std::vector<ValueType>>(getMatrixRowCount());
if(!this->cachedRowVector) {
this->cachedRowVector = std::make_unique<std::vector<ValueType>>(getMatrixRowCount());
}
std::vector<ValueType>* nextX = this->auxiliaryRowVector.get();
std::vector<ValueType>* nextX = this->cachedRowVector.get();
// Set up additional environment variables.
uint_fast64_t iterationCount = 0;
@ -264,10 +272,14 @@ namespace storm {
// If the last iteration did not write to the original x we have to swap the contents, because the
// output has to be written to the input parameter x.
if (currentX == this->auxiliaryRowVector.get()) {
if (currentX == this->cachedRowVector.get()) {
std::swap(x, *currentX);
}
if(!this->isCachingEnabled()) {
clearCache();
}
return iterationCount;
}
@ -282,12 +294,12 @@ namespace storm {
}
template<typename ValueType>
void GmmxxLinearEquationSolver<ValueType>::resetAuxiliaryData() const {
void GmmxxLinearEquationSolver<ValueType>::clearCache() const {
gmmxxA.reset();
iluPreconditioner.reset();
diagonalPreconditioner.reset();
jacobiDecomposition.reset();
LinearEquationSolver<ValueType>::resetAuxiliaryData();
LinearEquationSolver<ValueType>::clearCache();
}
template<typename ValueType>

7
src/solver/GmmxxLinearEquationSolver.h

@ -98,11 +98,8 @@ namespace storm {
void setSettings(GmmxxLinearEquationSolverSettings<ValueType> const& newSettings);
GmmxxLinearEquationSolverSettings<ValueType> const& getSettings() const;
/*
* Clears auxiliary data that has possibly been stored during previous calls of the solver.
*/
virtual void resetAuxiliaryData() const override;
virtual void clearCache() const override;
private:
/*!
@ -130,7 +127,7 @@ namespace storm {
// The settings used by the solver.
GmmxxLinearEquationSolverSettings<ValueType> settings;
// Auxiliary data obtained during solving
// cached data obtained during solving
mutable std::unique_ptr<gmm::csr_matrix<ValueType>> gmmxxA;
mutable std::unique_ptr<gmm::ilu_precond<gmm::csr_matrix<ValueType>>> iluPreconditioner;
mutable std::unique_ptr<gmm::diagonal_precond<gmm::csr_matrix<ValueType>>> diagonalPreconditioner;

42
src/solver/LinearEquationSolver.cpp

@ -14,21 +14,26 @@ namespace storm {
namespace solver {
template<typename ValueType>
LinearEquationSolver<ValueType>::LinearEquationSolver() {
LinearEquationSolver<ValueType>::LinearEquationSolver() : cachingEnabled(false) {
// Intentionally left empty.
}
template<typename ValueType>
void LinearEquationSolver<ValueType>::repeatedMultiply(std::vector<ValueType>& x, std::vector<ValueType> const* b, uint_fast64_t n) const {
if(!auxiliaryRowVector) {
auxiliaryRowVector = std::make_unique<std::vector<ValueType>>(getMatrixRowCount());
if(!cachedRowVector) {
cachedRowVector = std::make_unique<std::vector<ValueType>>(getMatrixRowCount());
}
// We enable caching for this. But remember how the old setting was
bool cachingWasEnabled = isCachingEnabled();
setCachingEnabled(true);
// Set up some temporary variables so that we can just swap pointers instead of copying the result after
// each iteration.
std::vector<ValueType>* currentX = &x;
std::vector<ValueType>* nextX = auxiliaryRowVector.get();
std::vector<ValueType>* nextX = cachedRowVector.get();
// Now perform matrix-vector multiplication as long as we meet the bound.
for (uint_fast64_t i = 0; i < n; ++i) {
@ -38,16 +43,37 @@ namespace storm {
// If we performed an odd number of repetitions, we need to swap the contents of currentVector and x,
// because the output is supposed to be stored in the input vector x.
if (currentX == auxiliaryRowVector.get()) {
if (currentX == cachedRowVector.get()) {
std::swap(x, *currentX);
}
// restore the old caching setting
setCachingEnabled(cachingWasEnabled);
if(!isCachingEnabled()) {
clearCache();
}
}
template<typename ValueType>
void LinearEquationSolver<ValueType>::resetAuxiliaryData() const {
auxiliaryRowVector.reset();
void LinearEquationSolver<ValueType>::setCachingEnabled(bool value) const {
if(cachingEnabled && !value) {
// caching will be turned off. Hence we clear the cache at this point
clearCache();
}
cachingEnabled = value;
}
template<typename ValueType>
bool LinearEquationSolver<ValueType>::isCachingEnabled() const {
return cachingEnabled;
}
template<typename ValueType>
void LinearEquationSolver<ValueType>::clearCache() const {
cachedRowVector.reset();
}
template<typename ValueType>
std::unique_ptr<LinearEquationSolver<ValueType>> LinearEquationSolverFactory<ValueType>::create(storm::storage::SparseMatrix<ValueType>&& matrix) const {
return create(matrix);

22
src/solver/LinearEquationSolver.h

@ -67,15 +67,26 @@ namespace storm {
* @param n The number of times to perform the multiplication.
*/
void repeatedMultiply(std::vector<ValueType>& x, std::vector<ValueType> const* b, uint_fast64_t n) const;
/*!
* Sets whether some of the generated data during solver calls should be cached.
* This possibly decreases the runtime of subsequent calls but also increases memory consumption.
*/
void setCachingEnabled(bool value) const;
/*!
* Retrieves whether some of the generated data during solver calls should be cached.
*/
bool isCachingEnabled() const;
/*
* Clears auxiliary data that has possibly been stored during previous calls of the solver.
* Clears the currently cached data that has been stored during previous calls of the solver.
*/
virtual void resetAuxiliaryData() const;
virtual void clearCache() const;
protected:
// auxiliary storage. If set, this vector has getMatrixRowCount() entries.
mutable std::unique_ptr<std::vector<ValueType>> auxiliaryRowVector;
mutable std::unique_ptr<std::vector<ValueType>> cachedRowVector;
private:
/*!
@ -87,6 +98,9 @@ namespace storm {
* Retrieves the column count of the matrix associated with this solver.
*/
virtual uint64_t getMatrixColumnCount() const = 0;
/// Whether some of the generated data during solver calls should be cached.
mutable bool cachingEnabled;
};

18
src/solver/MinMaxLinearEquationSolver.cpp

@ -18,7 +18,7 @@ namespace storm {
namespace solver {
template<typename ValueType>
MinMaxLinearEquationSolver<ValueType>::MinMaxLinearEquationSolver(OptimizationDirectionSetting direction) : direction(direction), trackScheduler(false) {
MinMaxLinearEquationSolver<ValueType>::MinMaxLinearEquationSolver(OptimizationDirectionSetting direction) : direction(direction), trackScheduler(false), cachingEnabled(false) {
// Intentionally left empty.
}
@ -80,7 +80,21 @@ namespace storm {
}
template<typename ValueType>
void MinMaxLinearEquationSolver<ValueType>::resetAuxiliaryData() const {
void MinMaxLinearEquationSolver<ValueType>::setCachingEnabled(bool value) {
if(cachingEnabled && !value) {
// caching will be turned off. Hence we clear the cache at this point
clearCache();
}
cachingEnabled = value;
}
template<typename ValueType>
bool MinMaxLinearEquationSolver<ValueType>::isCachingEnabled() const {
return cachingEnabled;
}
template<typename ValueType>
void MinMaxLinearEquationSolver<ValueType>::clearCache() const {
// Intentionally left empty.
}

20
src/solver/MinMaxLinearEquationSolver.h

@ -128,11 +128,22 @@ namespace storm {
* Gets whether the precision is taken to be absolute or relative
*/
virtual bool getRelative() const = 0;
/*!
* Sets whether some of the generated data during solver calls should be cached.
* This possibly decreases the runtime of subsequent calls but also increases memory consumption.
*/
void setCachingEnabled(bool value);
/*!
* Retrieves whether some of the generated data during solver calls should be cached.
*/
bool isCachingEnabled() const;
/*
* Clears auxiliary data that has possibly been stored during previous calls of the solver.
* Clears the currently cached data that has been stored during previous calls of the solver.
*/
virtual void resetAuxiliaryData() const;
virtual void clearCache() const;
protected:
@ -144,6 +155,11 @@ namespace storm {
/// The scheduler (if it could be successfully generated).
mutable boost::optional<std::unique_ptr<storm::storage::TotalScheduler>> scheduler;
private:
/// Whether some of the generated data during solver calls should be cached.
bool cachingEnabled;
};
template<typename ValueType>

43
src/solver/NativeLinearEquationSolver.cpp

@ -97,21 +97,21 @@ namespace storm {
void NativeLinearEquationSolver<ValueType>::setMatrix(storm::storage::SparseMatrix<ValueType> const& A) {
localA.reset();
this->A = &A;
resetAuxiliaryData();
clearCache();
}
template<typename ValueType>
void NativeLinearEquationSolver<ValueType>::setMatrix(storm::storage::SparseMatrix<ValueType>&& A) {
localA = std::make_unique<storm::storage::SparseMatrix<ValueType>>(std::move(A));
this->A = localA.get();
resetAuxiliaryData();
clearCache();
}
template<typename ValueType>
bool NativeLinearEquationSolver<ValueType>::solveEquations(std::vector<ValueType>& x, std::vector<ValueType> const& b) const {
if(!this->auxiliaryRowVector) {
this->auxiliaryRowVector = std::make_unique<std::vector<ValueType>>(getMatrixRowCount());
if(!this->cachedRowVector) {
this->cachedRowVector = std::make_unique<std::vector<ValueType>>(getMatrixRowCount());
}
if (this->getSettings().getSolutionMethod() == NativeLinearEquationSolverSettings<ValueType>::SolutionMethod::SOR || this->getSettings().getSolutionMethod() == NativeLinearEquationSolverSettings<ValueType>::SolutionMethod::GaussSeidel) {
@ -126,16 +126,20 @@ namespace storm {
A->performSuccessiveOverRelaxationStep(omega, x, b);
// Now check if the process already converged within our precision.
converged = storm::utility::vector::equalModuloPrecision<ValueType>(*this->auxiliaryRowVector, x, static_cast<ValueType>(this->getSettings().getPrecision()), this->getSettings().getRelativeTerminationCriterion()) || (this->hasCustomTerminationCondition() && this->getTerminationCondition().terminateNow(x));
converged = storm::utility::vector::equalModuloPrecision<ValueType>(*this->cachedRowVector, x, static_cast<ValueType>(this->getSettings().getPrecision()), this->getSettings().getRelativeTerminationCriterion()) || (this->hasCustomTerminationCondition() && this->getTerminationCondition().terminateNow(x));
// If we did not yet converge, we need to backup the contents of x.
if (!converged) {
*this->auxiliaryRowVector = x;
*this->cachedRowVector = x;
}
// Increase iteration count so we can abort if convergence is too slow.
++iterationCount;
}
if(!this->isCachingEnabled()) {
clearCache();
}
return converged;
@ -148,7 +152,7 @@ namespace storm {
std::vector<ValueType> const& jacobiD = jacobiDecomposition->second;
std::vector<ValueType>* currentX = &x;
std::vector<ValueType>* nextX = this->auxiliaryRowVector.get();
std::vector<ValueType>* nextX = this->cachedRowVector.get();
// Set up additional environment variables.
uint_fast64_t iterationCount = 0;
@ -172,10 +176,14 @@ namespace storm {
// If the last iteration did not write to the original x we have to swap the contents, because the
// output has to be written to the input parameter x.
if (currentX == this->auxiliaryRowVector.get()) {
if (currentX == this->cachedRowVector.get()) {
std::swap(x, *currentX);
}
if(!this->isCachingEnabled()) {
clearCache();
}
return iterationCount < this->getSettings().getMaximalNumberOfIterations();
}
}
@ -189,14 +197,19 @@ namespace storm {
}
} else {
// If the two vectors are aliases, we need to create a temporary.
if(!this->auxiliaryRowVector) {
this->auxiliaryRowVector = std::make_unique<std::vector<ValueType>>(getMatrixRowCount());
if(!this->cachedRowVector) {
this->cachedRowVector = std::make_unique<std::vector<ValueType>>(getMatrixRowCount());
}
A->multiplyWithVector(x, *this->auxiliaryRowVector);
A->multiplyWithVector(x, *this->cachedRowVector);
if (b != nullptr) {
storm::utility::vector::addVectors(*this->auxiliaryRowVector, *b, result);
storm::utility::vector::addVectors(*this->cachedRowVector, *b, result);
} else {
result.swap(*this->auxiliaryRowVector);
result.swap(*this->cachedRowVector);
}
if(!this->isCachingEnabled()) {
clearCache();
}
}
}
@ -212,9 +225,9 @@ namespace storm {
}
template<typename ValueType>
void NativeLinearEquationSolver<ValueType>::resetAuxiliaryData() const {
void NativeLinearEquationSolver<ValueType>::clearCache() const {
jacobiDecomposition.reset();
LinearEquationSolver<ValueType>::resetAuxiliaryData();
LinearEquationSolver<ValueType>::clearCache();
}
template<typename ValueType>

6
src/solver/NativeLinearEquationSolver.h

@ -55,10 +55,7 @@ namespace storm {
void setSettings(NativeLinearEquationSolverSettings<ValueType> const& newSettings);
NativeLinearEquationSolverSettings<ValueType> const& getSettings() const;
/*
* Clears auxiliary data that has possibly been stored during previous calls of the solver.
*/
virtual void resetAuxiliaryData() const override;
virtual void clearCache() const override;
private:
virtual uint64_t getMatrixRowCount() const override;
@ -75,6 +72,7 @@ namespace storm {
// The settings used by the solver.
NativeLinearEquationSolverSettings<ValueType> settings;
// cached auxiliary data
mutable std::unique_ptr<std::pair<storm::storage::SparseMatrix<ValueType>, std::vector<ValueType>>> jacobiDecomposition;
};

19
src/solver/StandardMinMaxLinearEquationSolver.cpp

@ -111,6 +111,7 @@ namespace storm {
// Create a solver that we will use throughout the procedure. We will modify the matrix in each iteration.
auto solver = linearEquationSolverFactory->create(std::move(submatrix));
solver->setCachingEnabled(true);
Status status = Status::InProgress;
uint64_t iterations = 0;
@ -170,6 +171,10 @@ namespace storm {
this->scheduler = std::make_unique<storm::storage::TotalScheduler>(std::move(scheduler));
}
if(!this->isCachingEnabled()) {
clearCache();
}
if(status == Status::Converged || status == Status::TerminatedEarly) {
return true;
} else{
@ -206,6 +211,7 @@ namespace storm {
bool StandardMinMaxLinearEquationSolver<ValueType>::solveEquationsValueIteration(OptimizationDirection dir, std::vector<ValueType>& x, std::vector<ValueType> const& b) const {
if(!linEqSolverA) {
linEqSolverA = linearEquationSolverFactory->create(A);
linEqSolverA->setCachingEnabled(true);
}
if (!auxiliaryRowVector.get()) {
@ -263,6 +269,10 @@ namespace storm {
this->scheduler = std::make_unique<storm::storage::TotalScheduler>(std::move(choices));
}
if(!this->isCachingEnabled()) {
clearCache();
}
if(status == Status::Converged || status == Status::TerminatedEarly) {
return true;
} else{
@ -274,6 +284,7 @@ namespace storm {
void StandardMinMaxLinearEquationSolver<ValueType>::repeatedMultiply(OptimizationDirection dir, std::vector<ValueType>& x, std::vector<ValueType>* b, uint_fast64_t n) const {
if(!linEqSolverA) {
linEqSolverA = linearEquationSolverFactory->create(A);
linEqSolverA->setCachingEnabled(true);
}
if (!auxiliaryRowVector.get()) {
@ -288,6 +299,10 @@ namespace storm {
// element of the min/max operator stack.
storm::utility::vector::reduceVectorMinOrMax(dir, multiplyResult, x, this->A.getRowGroupIndices());
}
if(!this->isCachingEnabled()) {
clearCache();
}
}
template<typename ValueType>
@ -324,11 +339,11 @@ namespace storm {
}
template<typename ValueType>
void StandardMinMaxLinearEquationSolver<ValueType>::resetAuxiliaryData() const {
void StandardMinMaxLinearEquationSolver<ValueType>::clearCache() const {
linEqSolverA.reset();
auxiliaryRowVector.reset();
auxiliaryRowGroupVector.reset();
MinMaxLinearEquationSolver<ValueType>::resetAuxiliaryData();
MinMaxLinearEquationSolver<ValueType>::clearCache();
}
template<typename ValueType>

3
src/solver/StandardMinMaxLinearEquationSolver.h

@ -44,7 +44,7 @@ namespace storm {
StandardMinMaxLinearEquationSolverSettings<ValueType> const& getSettings() const;
void setSettings(StandardMinMaxLinearEquationSolverSettings<ValueType> const& newSettings);
virtual void resetAuxiliaryData() const override;
virtual void clearCache() const override;
virtual ValueType getPrecision() const override;
virtual bool getRelative() const override;
@ -58,6 +58,7 @@ namespace storm {
Converged, TerminatedEarly, MaximalIterationsExceeded, InProgress
};
// possibly cached data
mutable std::unique_ptr<storm::solver::LinearEquationSolver<ValueType>> linEqSolverA;
mutable std::unique_ptr<std::vector<ValueType>> auxiliaryRowVector; // A.rowCount() entries
mutable std::unique_ptr<std::vector<ValueType>> auxiliaryRowGroupVector; // A.rowGroupCount() entries

2
src/utility/policyguessing.cpp

@ -68,6 +68,7 @@ namespace storm {
solveLinearEquationSystem(inducedA, x, inducedB, probGreater0States, prob0Value, solver.getPrecision(), solver.getRelative());
solver.setTrackScheduler();
solver.setCachingEnabled(true);
bool resultCorrect = false;
while(!resultCorrect){
solver.solveEquations(goal, x, b);
@ -171,6 +172,7 @@ namespace storm {
storm::utility::vector::selectVectorValues(subB, probGreater0States, b);
std::unique_ptr<storm::solver::GmmxxLinearEquationSolver<ValueType>> linEqSysSolver(static_cast<storm::solver::GmmxxLinearEquationSolver<ValueType>*>(storm::solver::GmmxxLinearEquationSolverFactory<ValueType>().create(eqSysA).release()));
linEqSysSolver->setCachingEnabled(true);
auto eqSettings = linEqSysSolver->getSettings();
eqSettings.setRelativeTerminationCriterion(relative);
eqSettings.setMaximalNumberOfIterations(500);

Loading…
Cancel
Save