krypy.recycling - Recycling Linear Systems Solvers¶
The recycling module provides functions for the solution of sequences of linear algebraic systems. Once a linear system has been solved, the generated data is examined and a deflation space is determined automatically for the solution of the next linear system. Several selection strategys are available.
- class krypy.recycling.RecyclingCg(*args, **kwargs)¶
Bases: krypy.recycling.linsys._RecyclingSolver
Recycling preconditioned CG method.
See _RecyclingSolver for the documentation of the available parameters.
- class krypy.recycling.RecyclingMinres(*args, **kwargs)¶
Bases: krypy.recycling.linsys._RecyclingSolver
Recycling preconditioned MINRES method.
See _RecyclingSolver for the documentation of the available parameters.
- class krypy.recycling.RecyclingGmres(*args, **kwargs)¶
Bases: krypy.recycling.linsys._RecyclingSolver
Recycling preconditioned GMRES method.
See _RecyclingSolver for the documentation of the available parameters.
- class krypy.recycling.linsys._RecyclingSolver(DeflatedSolver)¶
Base class for recycling solvers.
Initialize recycling solver base.
Parameters: DeflatedSolver – a deflated solver from deflation. After a run of the provided DeflatedSolver via solve(), the resulting instance of the DeflatedSolver is available in the attribute last_solver.
- last_solver = None¶
DeflatedSolver instance from last run of solve().
Instance of DeflatedSolver that resulted from the last call to solve(). Initialized with None before the first run.
- solve(linear_system, vector_factory=None, *args, **kwargs)¶
Solve the given linear system with recycling.
The deflation vectors and the Krylov subspace of the last solve() call are examined and considered for deflation according to the provided deflation_vector_factory.
Parameters: - linear_system – the LinearSystem that is about to be solved.
- vector_factory – (optional) An instance of a subclass of krypy.recycling.factories._DeflationVectorFactory.
All remaining arguments are passed to the DeflatedSolver.
Returns: instance of DeflatedSolver which was used to obtain the approximate solution. The approximate solution is available under the attribute xk.
- timings = None¶
Timings from last run of solve().
Timings of the vector factory runs and the actual solution processes.