Page MenuHomec4science

linear_solver.py
No OneTemporary

File Metadata

Created
Mon, Jun 3, 02:59

linear_solver.py

# Copyright (C) 2018 by the RROMPy authors
#
# This file is part of RROMPy.
#
# RROMPy is free software: you can redistribute it and/or modify
# it under the terms of the GNU Lesser General Public License as published by
# the Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
#
# RROMPy is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
# GNU Lesser General Public License for more details.
#
# You should have received a copy of the GNU Lesser General Public License
# along with RROMPy. If not, see <http://www.gnu.org/licenses/>.
#
import numpy.linalg as npla
import scipy.sparse.linalg as scspla
from rrompy.utilities.base import purgeDict
from rrompy.utilities.base.types import Tuple, DictAny
from rrompy.utilities.exception_manager import RROMPyException
__all__ = ['RROMPyLinearSolvers', 'setupSolver']
RROMPyLinearSolvers = {
"SOLVE" : (lambda A, b, kwargs: npla.solve(A, b, **kwargs),
[]),
"LSTSQ" : (lambda A, b, kwargs: npla.lstsq(A, b, **kwargs)[0],
["rcond"]),
"SPSOLVE" : (lambda A, b, kwargs: scspla.spsolve(A, b, **kwargs),
["permc_spec", "use_umfpack"]),
"BICG" : (lambda A, b, kwargs: scspla.bicg(A, b, **kwargs)[0],
["x0", "tol", "maxiter", "M", "callback", "atol"]),
"BICGSTAB" : (lambda A, b, kwargs: scspla.bicgstab(A, b, **kwargs)[0],
["x0", "tol", "maxiter", "M", "callback", "atol"]),
"CG" : (lambda A, b, kwargs: scspla.cg(A, b, **kwargs)[0],
["x0", "tol", "maxiter", "M", "callback", "atol"]),
"CGS" : (lambda A, b, kwargs: scspla.cgs(A, b, **kwargs)[0],
["x0", "tol", "maxiter", "M", "callback", "atol"]),
"GMRES" : (lambda A, b, kwargs: scspla.gmres(A, b, **kwargs)[0],
["x0", "tol", "restart", "maxiter", "M", "callback",
"restrt", "atol"]),
"LGMRES" : (lambda A, b, kwargs: scspla.lgmres(A, b, **kwargs)[0],
["x0", "tol", "maxiter", "M", "callback", "inner_m",
"outer_k", "outer_v", "store_outer_Av",
"prepend_outer_v", "atol"]),
"MINRES" : (lambda A, b, kwargs: scspla.minres(A, b, **kwargs)[0],
["x0", "shift", "tol", "maxiter", "M", "callback",
"show", "check"]),
"QMR" : (lambda A, b, kwargs: scspla.qmr(A, b, **kwargs)[0],
["x0", "tol", "maxiter", "M1", "M2", "callback",
"atol"]),
"GCROTMK" : (lambda A, b, kwargs: scspla.gcrotmk(A, b, **kwargs)[0],
["x0", "tol", "maxiter", "M", "callback", "m", "k", "CU",
"discard_C", "truncate", "atol"])
}
def setupSolver(solverType:str, solverArgs : DictAny = {})\
-> Tuple[callable, DictAny]:
solverType = solverType.upper()
if solverType not in RROMPyLinearSolvers.keys():
raise RROMPyException(("Solver type not recognized. Check allowed "
"values in RROMPyLinearSolvers.keys()."))
solver, solverArgsList = RROMPyLinearSolvers[solverType]
solverArgs = purgeDict(solverArgs, solverArgsList,
dictname = "{}.kwargs".format(solverType),
baselevel = 1)
return solver, solverArgs

Event Timeline