#include "GmmxxLinearEquationSolver.h" #include #include #include "src/adapters/GmmxxAdapter.h" #include "src/settings/SettingsManager.h" #include "src/utility/vector.h" #include "src/utility/constants.h" #include "src/exceptions/InvalidStateException.h" #include "src/settings/modules/GmmxxEquationSolverSettings.h" #include "src/solver/NativeLinearEquationSolver.h" #include "src/utility/gmm.h" namespace storm { namespace solver { template GmmxxLinearEquationSolverSettings::GmmxxLinearEquationSolverSettings() { // Get the settings object to customize linear solving. storm::settings::modules::GmmxxEquationSolverSettings const& settings = storm::settings::getModule(); // Get appropriate settings. maximalNumberOfIterations = settings.getMaximalIterationCount(); precision = settings.getPrecision(); relative = settings.getConvergenceCriterion() == storm::settings::modules::GmmxxEquationSolverSettings::ConvergenceCriterion::Relative; restart = settings.getRestartIterationCount(); // Determine the method to be used. storm::settings::modules::GmmxxEquationSolverSettings::LinearEquationMethod methodAsSetting = settings.getLinearEquationSystemMethod(); if (methodAsSetting == storm::settings::modules::GmmxxEquationSolverSettings::LinearEquationMethod::Bicgstab) { method = SolutionMethod::Bicgstab; } else if (methodAsSetting == storm::settings::modules::GmmxxEquationSolverSettings::LinearEquationMethod::Qmr) { method = SolutionMethod::Qmr; } else if (methodAsSetting == storm::settings::modules::GmmxxEquationSolverSettings::LinearEquationMethod::Gmres) { method = SolutionMethod::Gmres; } else if (methodAsSetting == storm::settings::modules::GmmxxEquationSolverSettings::LinearEquationMethod::Jacobi) { method = SolutionMethod::Jacobi; } // Check which preconditioner to use. storm::settings::modules::GmmxxEquationSolverSettings::PreconditioningMethod preconditionAsSetting = settings.getPreconditioningMethod(); if (preconditionAsSetting == storm::settings::modules::GmmxxEquationSolverSettings::PreconditioningMethod::Ilu) { preconditioner = Preconditioner::Ilu; } else if (preconditionAsSetting == storm::settings::modules::GmmxxEquationSolverSettings::PreconditioningMethod::Diagonal) { preconditioner = Preconditioner::Diagonal; } else if (preconditionAsSetting == storm::settings::modules::GmmxxEquationSolverSettings::PreconditioningMethod::None) { preconditioner = Preconditioner::None; } } template void GmmxxLinearEquationSolverSettings::setSolutionMethod(SolutionMethod const& method) { this->method = method; } template void GmmxxLinearEquationSolverSettings::setPreconditioner(Preconditioner const& preconditioner) { this->preconditioner = preconditioner; } template void GmmxxLinearEquationSolverSettings::setPrecision(ValueType precision) { this->precision = precision; } template void GmmxxLinearEquationSolverSettings::setMaximalNumberOfIterations(uint64_t maximalNumberOfIterations) { this->maximalNumberOfIterations = maximalNumberOfIterations; } template void GmmxxLinearEquationSolverSettings::setRelativeTerminationCriterion(bool value) { this->relative = value; } template void GmmxxLinearEquationSolverSettings::setNumberOfIterationsUntilRestart(uint64_t restart) { this->restart = restart; } template typename GmmxxLinearEquationSolverSettings::SolutionMethod GmmxxLinearEquationSolverSettings::getSolutionMethod() const { return method; } template typename GmmxxLinearEquationSolverSettings::Preconditioner GmmxxLinearEquationSolverSettings::getPreconditioner() const { return preconditioner; } template ValueType GmmxxLinearEquationSolverSettings::getPrecision() const { return precision; } template uint64_t GmmxxLinearEquationSolverSettings::getMaximalNumberOfIterations() const { return maximalNumberOfIterations; } template bool GmmxxLinearEquationSolverSettings::getRelativeTerminationCriterion() const { return relative; } template uint64_t GmmxxLinearEquationSolverSettings::getNumberOfIterationsUntilRestart() const { return restart; } template GmmxxLinearEquationSolver::GmmxxLinearEquationSolver(storm::storage::SparseMatrix const& A, GmmxxLinearEquationSolverSettings const& settings) : localA(nullptr), A(A), gmmxxMatrix(storm::adapters::GmmxxAdapter::toGmmxxSparseMatrix(A)), settings(settings) { // Intentionally left empty. } template GmmxxLinearEquationSolver::GmmxxLinearEquationSolver(storm::storage::SparseMatrix&& A, GmmxxLinearEquationSolverSettings const& settings) : localA(std::make_unique>(std::move(A))), A(*localA), gmmxxMatrix(storm::adapters::GmmxxAdapter::toGmmxxSparseMatrix(*localA)), settings(settings) { // Intentionally left empty. } template void GmmxxLinearEquationSolver::solveEquationSystem(std::vector& x, std::vector const& b, std::vector* multiplyResult) const { auto method = this->getSettings().getSolutionMethod(); auto preconditioner = this->getSettings().getPreconditioner(); STORM_LOG_INFO("Using method '" << method << "' with preconditioner '" << preconditioner << "' (max. " << this->getSettings().getMaximalNumberOfIterations() << " iterations)."); if (method == GmmxxLinearEquationSolverSettings::SolutionMethod::Jacobi && preconditioner != GmmxxLinearEquationSolverSettings::Preconditioner::None) { STORM_LOG_WARN("Jacobi method currently does not support preconditioners. The requested preconditioner will be ignored."); } if (method == GmmxxLinearEquationSolverSettings::SolutionMethod::Bicgstab || method == GmmxxLinearEquationSolverSettings::SolutionMethod::Qmr || method == GmmxxLinearEquationSolverSettings::SolutionMethod::Gmres) { // Prepare an iteration object that determines the accuracy and the maximum number of iterations. gmm::iteration iter(this->getSettings().getPrecision(), 0, this->getSettings().getMaximalNumberOfIterations()); if (method == GmmxxLinearEquationSolverSettings::SolutionMethod::Bicgstab) { if (preconditioner == GmmxxLinearEquationSolverSettings::Preconditioner::Ilu) { gmm::bicgstab(*gmmxxMatrix, x, b, gmm::ilu_precond>(*gmmxxMatrix), iter); } else if (preconditioner == GmmxxLinearEquationSolverSettings::Preconditioner::Diagonal) { gmm::bicgstab(*gmmxxMatrix, x, b, gmm::diagonal_precond>(*gmmxxMatrix), iter); } else if (preconditioner == GmmxxLinearEquationSolverSettings::Preconditioner::None) { gmm::bicgstab(*gmmxxMatrix, x, b, gmm::identity_matrix(), iter); } } else if (method == GmmxxLinearEquationSolverSettings::SolutionMethod::Qmr) { if (preconditioner == GmmxxLinearEquationSolverSettings::Preconditioner::Ilu) { gmm::qmr(*gmmxxMatrix, x, b, gmm::ilu_precond>(*gmmxxMatrix), iter); } else if (preconditioner == GmmxxLinearEquationSolverSettings::Preconditioner::Diagonal) { gmm::qmr(*gmmxxMatrix, x, b, gmm::diagonal_precond>(*gmmxxMatrix), iter); } else if (preconditioner == GmmxxLinearEquationSolverSettings::Preconditioner::None) { gmm::qmr(*gmmxxMatrix, x, b, gmm::identity_matrix(), iter); } } else if (method == GmmxxLinearEquationSolverSettings::SolutionMethod::Gmres) { if (preconditioner == GmmxxLinearEquationSolverSettings::Preconditioner::Ilu) { gmm::gmres(*gmmxxMatrix, x, b, gmm::ilu_precond>(*gmmxxMatrix), this->getSettings().getNumberOfIterationsUntilRestart(), iter); } else if (preconditioner == GmmxxLinearEquationSolverSettings::Preconditioner::Diagonal) { gmm::gmres(*gmmxxMatrix, x, b, gmm::diagonal_precond>(*gmmxxMatrix), this->getSettings().getNumberOfIterationsUntilRestart(), iter); } else if (preconditioner == GmmxxLinearEquationSolverSettings::Preconditioner::None) { gmm::gmres(*gmmxxMatrix, x, b, gmm::identity_matrix(), this->getSettings().getNumberOfIterationsUntilRestart(), iter); } } // Check if the solver converged and issue a warning otherwise. if (iter.converged()) { STORM_LOG_INFO("Iterative solver converged after " << iter.get_iteration() << " iterations."); } else { STORM_LOG_WARN("Iterative solver did not converge."); } } else if (method == GmmxxLinearEquationSolverSettings::SolutionMethod::Jacobi) { uint_fast64_t iterations = solveLinearEquationSystemWithJacobi(A, x, b, multiplyResult); // Check if the solver converged and issue a warning otherwise. if (iterations < this->getSettings().getMaximalNumberOfIterations()) { STORM_LOG_INFO("Iterative solver converged after " << iterations << " iterations."); } else { STORM_LOG_WARN("Iterative solver did not converge."); } } } template void GmmxxLinearEquationSolver::performMatrixVectorMultiplication(std::vector& x, std::vector& result, std::vector const* b) const { if (b) { gmm::mult_add(*gmmxxMatrix, x, *b, result); } else { gmm::mult(*gmmxxMatrix, x, result); } } template uint_fast64_t GmmxxLinearEquationSolver::solveLinearEquationSystemWithJacobi(storm::storage::SparseMatrix const& A, std::vector& x, std::vector const& b, std::vector* multiplyResult) const { // Get a Jacobi decomposition of the matrix A. std::pair, std::vector> jacobiDecomposition = A.getJacobiDecomposition(); // Convert the LU matrix to gmm++'s format. std::unique_ptr> gmmLU = storm::adapters::GmmxxAdapter::toGmmxxSparseMatrix(std::move(jacobiDecomposition.first)); // To avoid copying the contents of the vector in the loop, we create a temporary x to swap with. bool multiplyResultProvided = true; std::vector* nextX = multiplyResult; if (nextX == nullptr) { nextX = new std::vector(x.size()); multiplyResultProvided = false; } std::vector const* copyX = nextX; std::vector* currentX = &x; // Target vector for precision calculation. std::vector tmpX(x.size()); // Set up additional environment variables. uint_fast64_t iterationCount = 0; bool converged = false; while (!converged && iterationCount < this->getSettings().getMaximalNumberOfIterations() && !(this->hasCustomTerminationCondition() && this->getTerminationCondition().terminateNow(*currentX))) { // Compute D^-1 * (b - LU * x) and store result in nextX. gmm::mult(*gmmLU, *currentX, tmpX); gmm::add(b, gmm::scaled(tmpX, -storm::utility::one()), tmpX); storm::utility::vector::multiplyVectorsPointwise(jacobiDecomposition.second, tmpX, *nextX); // Now check if the process already converged within our precision. converged = storm::utility::vector::equalModuloPrecision(*currentX, *nextX, this->getSettings().getPrecision(), this->getSettings().getRelativeTerminationCriterion()); // Swap the two pointers as a preparation for the next iteration. std::swap(nextX, currentX); // Increase iteration count so we can abort if convergence is too slow. ++iterationCount; } // 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 == copyX) { std::swap(x, *currentX); } // If the vector for the temporary multiplication result was not provided, we need to delete it. if (!multiplyResultProvided) { delete copyX; } return iterationCount; } template GmmxxLinearEquationSolverSettings& GmmxxLinearEquationSolver::getSettings() { return settings; } template GmmxxLinearEquationSolverSettings const& GmmxxLinearEquationSolver::getSettings() const { return settings; } template std::unique_ptr> GmmxxLinearEquationSolverFactory::create(storm::storage::SparseMatrix const& matrix) const { return std::make_unique>(matrix, settings); } template std::unique_ptr> GmmxxLinearEquationSolverFactory::create(storm::storage::SparseMatrix&& matrix) const { return std::make_unique>(std::move(matrix), settings); } template std::unique_ptr> GmmxxLinearEquationSolverFactory::clone() const { return std::make_unique>(*this); } template GmmxxLinearEquationSolverSettings& GmmxxLinearEquationSolverFactory::getSettings() { return settings; } template GmmxxLinearEquationSolverSettings const& GmmxxLinearEquationSolverFactory::getSettings() const { return settings; } // Explicitly instantiate the solver. template class GmmxxLinearEquationSolverSettings; template class GmmxxLinearEquationSolver; template class GmmxxLinearEquationSolverFactory; } }