ScipyKrylov#

ScipyKrylov is an iterative linear solver that wraps the methods found in scipy.sparse.linalg. The default method is “gmres”, or the Generalized Minimal RESidual method. Support for other scipy.sparse.linalg solvers will be added over time. This linear solver is capable of handling any system topology very effectively. It also solves all subsystems below it in the hierarchy, so assigning different solvers to subsystems will have no effect on the solution at this level.

This is a serial solver, so it should never be used under MPI; use PETScKrylov instead.

Here, we calculate the total derivatives across the Sellar system.

import numpy as np
import openmdao.api as om
from openmdao.test_suite.components.sellar import SellarDis1withDerivatives, SellarDis2withDerivatives

prob = om.Problem()
model = prob.model

model.add_subsystem('d1', SellarDis1withDerivatives(), promotes=['x', 'z', 'y1', 'y2'])
model.add_subsystem('d2', SellarDis2withDerivatives(), promotes=['z', 'y1', 'y2'])

model.add_subsystem('obj_cmp', om.ExecComp('obj = x**2 + z[1] + y1 + exp(-y2)',
                                           z=np.array([0.0, 0.0]), x=0.0),
                    promotes=['obj', 'x', 'z', 'y1', 'y2'])

model.add_subsystem('con_cmp1', om.ExecComp('con1 = 3.16 - y1'), promotes=['con1', 'y1'])
model.add_subsystem('con_cmp2', om.ExecComp('con2 = y2 - 24.0'), promotes=['con2', 'y2'])

model.nonlinear_solver = om.NonlinearBlockGS()

model.linear_solver = om.ScipyKrylov()

prob.setup()

prob.set_val('x', 1.)
prob.set_val('z', np.array([5.0, 2.0]))

prob.run_model()

wrt = ['z']
of = ['obj']

J = prob.compute_totals(of=of, wrt=wrt, return_format='flat_dict')
NL: NLBGS Converged in 8 iterations

ScipyKrylov Options#

OptionDefaultAcceptable ValuesAcceptable TypesDescription
assemble_jacFalse[True, False]['bool']Activates use of assembled jacobian by this solver.
atol1e-12N/AN/Aabsolute error tolerance
err_on_non_convergeFalse[True, False]['bool']When True, AnalysisError will be raised if we don't converge.
iprint1N/A['int']whether to print output
maxiter1000N/A['int']maximum number of iterations
restart20N/A['int']Number of iterations between restarts. Larger values increase iteration cost, but may be necessary for convergence. This option applies only to gmres.
rhs_checkingFalseN/A['bool', 'dict']If True, check RHS vs. cache and/or zero to avoid some solves.Can also be set to a dict of options for the LinearRHSChecker to allow finer control over it. Allowed options are: ('check_zero', 'rtol', 'atol', 'max_cache_entries', 'collect_stats', 'auto', 'verbose')
rtol1e-10N/AN/Arelative error tolerance
solvergmres['gmres']N/Afunction handle for actual solver

ScipyKrylov Constructor#

The call signature for the ScipyKrylov constructor is:

ScipyKrylov.__init__(**kwargs)[source]

Declare the solver option.

ScipyKrylov Option Examples#

maxiter

maxiter lets you specify the maximum number of GMRES iterations to apply. The default maximum is 1000, which is much higher than the other linear solvers because each multiplication by the system Jacobian is considered to be an iteration. You may have to decrease this value if you have a coupled system that is converging very slowly. (Of course, in such a case, it may be better to add a preconditioner.) Alternatively, you may have to raise it if you have an extremely large number of components in your system (a 1000-component ring would need 1000 iterations just to make it around once.)

This example shows what happens if you set maxiter too low. Since iteration is terminated after only 3 iterations, the solver does not converge and the derivatives will not be correct.

from openmdao.test_suite.components.sellar import SellarDis1withDerivatives, SellarDis2withDerivatives

prob = om.Problem()
model = prob.model

model.add_subsystem('d1', SellarDis1withDerivatives(), promotes=['x', 'z', 'y1', 'y2'])
model.add_subsystem('d2', SellarDis2withDerivatives(), promotes=['z', 'y1', 'y2'])

model.add_subsystem('obj_cmp', om.ExecComp('obj = x**2 + z[1] + y1 + exp(-y2)',
                                           z=np.array([0.0, 0.0]), x=0.0),
                    promotes=['obj', 'x', 'z', 'y1', 'y2'])

model.add_subsystem('con_cmp1', om.ExecComp('con1 = 3.16 - y1'), promotes=['con1', 'y1'])
model.add_subsystem('con_cmp2', om.ExecComp('con2 = y2 - 24.0'), promotes=['con2', 'y2'])

model.nonlinear_solver = om.NonlinearBlockGS()

model.linear_solver = om.ScipyKrylov()
model.linear_solver.options['maxiter'] = 3

prob.setup()

prob.set_val('x', 1.)
prob.set_val('z', np.array([5.0, 2.0]))

prob.run_model()

wrt = ['z']
of = ['obj']

J = prob.compute_totals(of=of, wrt=wrt, return_format='flat_dict')
NL: NLBGS Converged in 8 iterations
LN: SCIPYSolver 'LN: SCIPY' on system '' failed to converge in 3 iterations.
print(J['obj', 'z'][0][0])
print(J['obj', 'z'][0][1])
9.265405444309767
1.8724662355885655

Depending on your situation, you may want to set the err_on_non_converge option to True, so that this will raise an AnalysisError.

from openmdao.test_suite.components.sellar import SellarDis1withDerivatives, SellarDis2withDerivatives
from openmdao.api import AnalysisError

prob = om.Problem()
model = prob.model

model.add_subsystem('d1', SellarDis1withDerivatives(), promotes=['x', 'z', 'y1', 'y2'])
model.add_subsystem('d2', SellarDis2withDerivatives(), promotes=['z', 'y1', 'y2'])

model.add_subsystem('obj_cmp', om.ExecComp('obj = x**2 + z[1] + y1 + exp(-y2)',
                                           z=np.array([0.0, 0.0]), x=0.0),
                    promotes=['obj', 'x', 'z', 'y1', 'y2'])

model.add_subsystem('con_cmp1', om.ExecComp('con1 = 3.16 - y1'), promotes=['con1', 'y1'])
model.add_subsystem('con_cmp2', om.ExecComp('con2 = y2 - 24.0'), promotes=['con2', 'y2'])

model.nonlinear_solver = om.NonlinearBlockGS()

model.linear_solver = om.ScipyKrylov()
model.linear_solver.options['maxiter'] = 3
model.linear_solver.options['err_on_non_converge'] = True

prob.setup()

prob.set_val('x', 1.)
prob.set_val('z', np.array([5.0, 2.0]))

prob.run_model()

wrt = ['z']
of = ['obj']

try:
    J = prob.compute_totals(of=of, wrt=wrt, return_format='flat_dict')
except AnalysisError:
    # Special handling for a failure to converge when computing derivatives
    print("===\nLinear Solver failed to converge, derivatives are not valid.\n===")
NL: NLBGS Converged in 8 iterations
LN: SCIPYSolver 'LN: SCIPY' on system '' failed to converge in 3 iterations.
===
Linear Solver failed to converge, derivatives are not valid.
===

atol

Here, we set the absolute tolerance to a much tighter value (default is 1.0e-12) to show what happens. In practice, the tolerance serves a dual role in GMRES. In addition to being a termination criteria, the tolerance also defines what GMRES considers to be tiny. Tiny numbers are replaced by zero when the argument vector is normalized at the start of each new matrix-vector product. The end result here is that we iterate longer to get a marginally better answer.

You may need to adjust this setting if you have abnormally large or small values in your global Jacobian.

from openmdao.test_suite.components.sellar import SellarDis1withDerivatives, SellarDis2withDerivatives

prob = om.Problem()
model = prob.model

model.add_subsystem('d1', SellarDis1withDerivatives(), promotes=['x', 'z', 'y1', 'y2'])
model.add_subsystem('d2', SellarDis2withDerivatives(), promotes=['z', 'y1', 'y2'])

model.add_subsystem('obj_cmp', om.ExecComp('obj = x**2 + z[1] + y1 + exp(-y2)',
                                           z=np.array([0.0, 0.0]), x=0.0),
                    promotes=['obj', 'x', 'z', 'y1', 'y2'])

model.add_subsystem('con_cmp1', om.ExecComp('con1 = 3.16 - y1'), promotes=['con1', 'y1'])
model.add_subsystem('con_cmp2', om.ExecComp('con2 = y2 - 24.0'), promotes=['con2', 'y2'])

model.nonlinear_solver = om.NonlinearBlockGS()

model.linear_solver = om.ScipyKrylov()
model.linear_solver.options['atol'] = 1.0e-20

prob.setup()

prob.set_val('x', 1.)
prob.set_val('z', np.array([5.0, 2.0]))

prob.run_model()

wrt = ['z']
of = ['obj']

J = prob.compute_totals(of=of, wrt=wrt, return_format='flat_dict')
NL: NLBGS Converged in 8 iterations
print(J['obj', 'z'][0][0])
print(J['obj', 'z'][0][1])
9.610010556989941
1.7844853356313644

rtol

The ‘rtol’ setting is not supported by SciPy GMRES for versions earlier than SciPy 1.12.0.

Specifying a Preconditioner#

You can specify a preconditioner to improve the convergence of the iterative linear solution by setting the precon attribute. The motivation for using a preconditioner is the observation that iterative methods have better convergence properties if the linear system has a smaller condition number, so the goal of the preconditioner is to improve the condition number in part or all of the Jacobian.

Here, we add a Gauss-Seidel preconditioner to a problem that contains two subgroups, each with an implicit component that implements a quadratic equation. These are solved together by a Newton solver at the top. The goal of the preconditioner here is to solve the smaller linear systems for the quadratic components independently, and use that solution to precondition the full system solution. This is accomplished by setting up the linear solver hierarchy so that the preconditioner is LinearBlockGS and the subsytems sub1 and sub2 contain a DirectSolver.

Note that the number of GMRES iterations is lower when using the preconditioner.

from openmdao.test_suite.components.quad_implicit import QuadraticComp

prob = om.Problem()
model = prob.model

sub1 = model.add_subsystem('sub1', om.Group())
sub1.add_subsystem('q1', QuadraticComp())
sub1.add_subsystem('z1', om.ExecComp('y = -6.0 + .01 * x'))
sub2 = model.add_subsystem('sub2', om.Group())
sub2.add_subsystem('q2', QuadraticComp())
sub2.add_subsystem('z2', om.ExecComp('y = -6.0 + .01 * x'))

model.connect('sub1.q1.x', 'sub1.z1.x')
model.connect('sub1.z1.y', 'sub2.q2.c')
model.connect('sub2.q2.x', 'sub2.z2.x')
model.connect('sub2.z2.y', 'sub1.q1.c')

model.nonlinear_solver = om.NewtonSolver(solve_subsystems=False)
model.linear_solver = om.ScipyKrylov()

prob.setup()

model.sub1.linear_solver = om.DirectSolver()
model.sub2.linear_solver = om.DirectSolver()

model.linear_solver.precon = om.LinearBlockGS()

prob.set_solver_print(level=2)
prob.run_model()
NL: Newton 0 ; 10 1
|  | precon: LNBGS 0 ; 10 1
|  | precon: LNBGS 1 ; 6.9399 0.69399
|  | precon: LNBGS 2 ; 0.00069399 6.9399e-05
|  | precon: LNBGS 3 ; 6.93990003e-08 6.93990003e-09
|  | precon: LNBGS 4 ; 6.93933799e-12 6.93933799e-13
|  | precon:LN: LNBGS Converged
|  | precon: LNBGS 0 ; 6.93933799e-12 1
|  | precon:LN: LNBGS Converged
|  | precon: LNBGS 0 ; 0 0
|  | precon:LN: LNBGS Converged
|  LN: SCIPY 0 ; 2.86880433e-16 1
NL: Newton 1 ; 49.9085269 4.99085269
|  | precon: LNBGS 0 ; 148.641386 1
|  | precon: LNBGS 1 ; 6.91321409 0.046509349
|  | precon: LNBGS 2 ; 4.16646906e-06 2.80303432e-08
|  | precon: LNBGS 3 ; 2.50822592e-12 1.68743443e-14
|  | precon:LN: LNBGS Converged
|  | precon: LNBGS 0 ; 247.507374 1
|  | precon: LNBGS 1 ; 13.853804 0.0559732979
|  | precon: LNBGS 2 ; 8.349437e-06 3.37340939e-08
|  | precon: LNBGS 3 ; 5.03064257e-12 2.0325223e-14
|  | precon:LN: LNBGS Converged
|  | precon: LNBGS 0 ; 0 0
|  | precon:LN: LNBGS Converged
|  LN: SCIPY 0 ; 0 0
NL: Newton 2 ; 10.5986027 1.05986027
|  | precon: LNBGS 0 ; 18.1127188 1
|  | precon: LNBGS 1 ; 0.0172701646 0.000953482735
|  | precon: LNBGS 2 ; 3.14865982e-08 1.73836952e-09
|  | precon: LNBGS 3 ; 5.7738429e-14 3.18772846e-15
|  | precon:LN: LNBGS Converged
|  | precon: LNBGS 0 ; 46.8240402 1
|  | precon: LNBGS 1 ; 0.044645914 0.000953482735
|  | precon: LNBGS 2 ; 8.13974843e-08 1.73836952e-09
|  | precon: LNBGS 3 ; 1.49213975e-13 3.18669585e-15
|  | precon:LN: LNBGS Converged
|  | precon: LNBGS 0 ; 0 0
|  | precon:LN: LNBGS Converged
|  LN: SCIPY 0 ; 0 0
NL: Newton 3 ; 1.44424034 0.144424034
|  | precon: LNBGS 0 ; 6.26588166 1
|  | precon: LNBGS 1 ; 0.00821261904 0.0013106885
|  | precon: LNBGS 2 ; 2.83220809e-08 4.52004721e-09
|  | precon: LNBGS 3 ; 9.76996262e-14 1.55923191e-14
|  | precon:LN: LNBGS Converged
|  | precon: LNBGS 0 ; 13.9760037 1
|  | precon: LNBGS 1 ; 0.0183181873 0.0013106885
|  | precon: LNBGS 2 ; 6.31721964e-08 4.52004721e-09
|  | precon: LNBGS 3 ; 2.17603713e-13 1.55698094e-14
|  | precon:LN: LNBGS Converged
|  | precon: LNBGS 0 ; 0 0
|  | precon:LN: LNBGS Converged
|  LN: SCIPY 0 ; 0 0
NL: Newton 4 ; 0.0506752521 0.00506752521
|  | precon: LNBGS 0 ; 1.29221459 1
|  | precon: LNBGS 1 ; 0.00182151595 0.00140960795
|  | precon: LNBGS 2 ; 7.26768965e-09 5.62421266e-09
|  | precon: LNBGS 3 ; 2.89976376e-14 2.24402649e-14
|  | precon:LN: LNBGS Converged
|  | precon: LNBGS 0 ; 2.63510443 1
|  | precon: LNBGS 1 ; 0.00371446414 0.00140960795
|  | precon: LNBGS 2 ; 1.48203876e-08 5.62421262e-09
|  | precon: LNBGS 3 ; 5.90639693e-14 2.24142803e-14
|  | precon:LN: LNBGS Converged
|  | precon: LNBGS 0 ; 0 0
|  | precon:LN: LNBGS Converged
|  LN: SCIPY 0 ; 2.39466931e-15 1
NL: Newton 5 ; 7.21617369e-05 7.21617369e-06
|  | precon: LNBGS 0 ; 0.0504587669 1
|  | precon: LNBGS 1 ; 7.13302261e-05 0.00141363395
|  | precon: LNBGS 2 ; 2.86232687e-10 5.67260567e-09
|  | precon: LNBGS 3 ; 1.14859701e-15 2.27630811e-14
|  | precon:LN: LNBGS Converged
|  | precon: LNBGS 0 ; 0.100989696 1
|  | precon: LNBGS 1 ; 0.000142762462 0.00141363395
|  | precon: LNBGS 2 ; 5.7287472e-10 5.67260567e-09
|  | precon: LNBGS 3 ; 2.29677388e-15 2.27426558e-14
|  | precon:LN: LNBGS Converged
|  | precon: LNBGS 0 ; 0 0
|  | precon:LN: LNBGS Converged
|  LN: SCIPY 0 ; 0 0
NL: Newton 6 ; 1.47167278e-10 1.47167278e-11
NL: Newton Converged
print(prob.get_val('sub1.q1.x'))
print(prob.get_val('sub2.q2.x'))
[1.9960048]
[1.9960048]