1

2

3

4

5

6

7

8

9

10

11

12

13

14

15

16

17

18

19

20

21

22

23

24

25

26

27

28

29

30

31

32

33

34

35

36

37

38

39

40

41

42

43

44

45

46

47

48

49

50

51

52

53

54

55

56

57

58

59

60

61

62

63

64

65

66

67

68

69

70

71

72

73

74

75

76

77

78

79

80

81

82

83

84

85

86

87

88

89

90

91

92

93

94

95

96

97

98

99

100

101

102

103

104

105

106

107

108

109

110

111

112

113

114

115

116

117

118

119

120

121

122

123

124

125

126

127

128

129

130

131

132

133

134

135

136

137

138

139

140

141

142

143

144

145

146

147

148

149

150

151

152

153

154

155

156

157

158

159

160

161

162

163

164

165

166

167

168

169

170

171

172

173

174

175

176

177

178

179

180

181

182

183

184

185

186

187

188

189

190

191

192

193

194

195

196

197

198

199

200

201

202

203

204

205

206

207

208

209

210

211

212

213

214

215

216

217

218

219

220

221

222

223

224

225

226

227

228

229

230

231

232

233

234

235

236

237

238

239

240

241

242

243

244

245

246

247

248

249

250

251

252

253

254

255

256

257

258

259

260

261

262

263

264

265

266

267

268

269

270

271

272

273

274

275

276

277

278

279

280

281

282

283

284

285

286

287

288

289

290

291

292

293

294

295

296

297

298

299

300

301

302

303

304

305

306

307

308

309

310

311

312

313

314

315

316

317

318

319

320

321

322

323

324

325

326

327

328

329

330

331

332

333

334

335

336

337

338

339

340

341

342

343

344

345

346

347

348

349

350

351

352

353

354

355

356

357

358

359

360

361

362

363

364

365

366

367

368

369

370

371

372

373

374

375

376

377

378

379

380

381

382

383

384

385

386

387

388

389

390

391

392

393

394

395

396

397

398

399

400

401

402

403

404

405

406

407

408

409

410

411

412

413

414

415

416

417

418

419

420

421

422

423

424

425

426

427

428

429

430

431

432

433

434

435

436

437

438

439

440

441

442

443

444

445

446

447

448

449

450

451

452

453

454

455

456

457

458

459

460

461

462

463

464

465

466

467

468

469

470

471

472

473

474

475

476

477

478

479

480

481

482

483

484

485

486

487

488

489

490

491

492

493

494

495

496

497

498

499

500

501

502

503

504

505

506

507

508

509

510

511

512

513

514

515

516

517

518

519

520

521

522

523

524

525

526

527

528

529

530

531

532

533

534

535

536

537

538

539

540

541

542

543

544

545

546

547

548

549

550

551

552

553

554

555

556

557

558

559

560

561

562

563

564

565

566

567

568

569

570

571

572

573

574

575

576

577

578

579

580

581

582

583

584

585

586

587

588

589

590

591

592

593

594

595

596

597

598

599

600

601

602

603

604

605

606

607

608

609

610

611

612

613

614

615

616

617

618

619

620

621

622

623

624

625

626

627

628

629

630

631

632

633

634

635

636

637

638

639

640

""" 

Unified interfaces to root finding algorithms. 

 

Functions 

--------- 

- root : find a root of a vector function. 

""" 

from __future__ import division, print_function, absolute_import 

 

__all__ = ['root'] 

 

import numpy as np 

 

from scipy._lib.six import callable 

 

from warnings import warn 

 

from .optimize import MemoizeJac, OptimizeResult, _check_unknown_options 

from .minpack import _root_hybr, leastsq 

from ._spectral import _root_df_sane 

from . import nonlin 

 

 

def root(fun, x0, args=(), method='hybr', jac=None, tol=None, callback=None, 

options=None): 

""" 

Find a root of a vector function. 

 

Parameters 

---------- 

fun : callable 

A vector function to find a root of. 

x0 : ndarray 

Initial guess. 

args : tuple, optional 

Extra arguments passed to the objective function and its Jacobian. 

method : str, optional 

Type of solver. Should be one of 

 

- 'hybr' :ref:`(see here) <optimize.root-hybr>` 

- 'lm' :ref:`(see here) <optimize.root-lm>` 

- 'broyden1' :ref:`(see here) <optimize.root-broyden1>` 

- 'broyden2' :ref:`(see here) <optimize.root-broyden2>` 

- 'anderson' :ref:`(see here) <optimize.root-anderson>` 

- 'linearmixing' :ref:`(see here) <optimize.root-linearmixing>` 

- 'diagbroyden' :ref:`(see here) <optimize.root-diagbroyden>` 

- 'excitingmixing' :ref:`(see here) <optimize.root-excitingmixing>` 

- 'krylov' :ref:`(see here) <optimize.root-krylov>` 

- 'df-sane' :ref:`(see here) <optimize.root-dfsane>` 

 

jac : bool or callable, optional 

If `jac` is a Boolean and is True, `fun` is assumed to return the 

value of Jacobian along with the objective function. If False, the 

Jacobian will be estimated numerically. 

`jac` can also be a callable returning the Jacobian of `fun`. In 

this case, it must accept the same arguments as `fun`. 

tol : float, optional 

Tolerance for termination. For detailed control, use solver-specific 

options. 

callback : function, optional 

Optional callback function. It is called on every iteration as 

``callback(x, f)`` where `x` is the current solution and `f` 

the corresponding residual. For all methods but 'hybr' and 'lm'. 

options : dict, optional 

A dictionary of solver options. E.g. `xtol` or `maxiter`, see 

:obj:`show_options()` for details. 

 

Returns 

------- 

sol : OptimizeResult 

The solution represented as a ``OptimizeResult`` object. 

Important attributes are: ``x`` the solution array, ``success`` a 

Boolean flag indicating if the algorithm exited successfully and 

``message`` which describes the cause of the termination. See 

`OptimizeResult` for a description of other attributes. 

 

See also 

-------- 

show_options : Additional options accepted by the solvers 

 

Notes 

----- 

This section describes the available solvers that can be selected by the 

'method' parameter. The default method is *hybr*. 

 

Method *hybr* uses a modification of the Powell hybrid method as 

implemented in MINPACK [1]_. 

 

Method *lm* solves the system of nonlinear equations in a least squares 

sense using a modification of the Levenberg-Marquardt algorithm as 

implemented in MINPACK [1]_. 

 

Method *df-sane* is a derivative-free spectral method. [3]_ 

 

Methods *broyden1*, *broyden2*, *anderson*, *linearmixing*, 

*diagbroyden*, *excitingmixing*, *krylov* are inexact Newton methods, 

with backtracking or full line searches [2]_. Each method corresponds 

to a particular Jacobian approximations. See `nonlin` for details. 

 

- Method *broyden1* uses Broyden's first Jacobian approximation, it is 

known as Broyden's good method. 

- Method *broyden2* uses Broyden's second Jacobian approximation, it 

is known as Broyden's bad method. 

- Method *anderson* uses (extended) Anderson mixing. 

- Method *Krylov* uses Krylov approximation for inverse Jacobian. It 

is suitable for large-scale problem. 

- Method *diagbroyden* uses diagonal Broyden Jacobian approximation. 

- Method *linearmixing* uses a scalar Jacobian approximation. 

- Method *excitingmixing* uses a tuned diagonal Jacobian 

approximation. 

 

.. warning:: 

 

The algorithms implemented for methods *diagbroyden*, 

*linearmixing* and *excitingmixing* may be useful for specific 

problems, but whether they will work may depend strongly on the 

problem. 

 

.. versionadded:: 0.11.0 

 

References 

---------- 

.. [1] More, Jorge J., Burton S. Garbow, and Kenneth E. Hillstrom. 

1980. User Guide for MINPACK-1. 

.. [2] C. T. Kelley. 1995. Iterative Methods for Linear and Nonlinear 

Equations. Society for Industrial and Applied Mathematics. 

<http://www.siam.org/books/kelley/fr16/index.php> 

.. [3] W. La Cruz, J.M. Martinez, M. Raydan. Math. Comp. 75, 1429 (2006). 

 

Examples 

-------- 

The following functions define a system of nonlinear equations and its 

jacobian. 

 

>>> def fun(x): 

... return [x[0] + 0.5 * (x[0] - x[1])**3 - 1.0, 

... 0.5 * (x[1] - x[0])**3 + x[1]] 

 

>>> def jac(x): 

... return np.array([[1 + 1.5 * (x[0] - x[1])**2, 

... -1.5 * (x[0] - x[1])**2], 

... [-1.5 * (x[1] - x[0])**2, 

... 1 + 1.5 * (x[1] - x[0])**2]]) 

 

A solution can be obtained as follows. 

 

>>> from scipy import optimize 

>>> sol = optimize.root(fun, [0, 0], jac=jac, method='hybr') 

>>> sol.x 

array([ 0.8411639, 0.1588361]) 

 

""" 

if not isinstance(args, tuple): 

args = (args,) 

 

meth = method.lower() 

if options is None: 

options = {} 

 

if callback is not None and meth in ('hybr', 'lm'): 

warn('Method %s does not accept callback.' % method, 

RuntimeWarning) 

 

# fun also returns the jacobian 

if not callable(jac) and meth in ('hybr', 'lm'): 

if bool(jac): 

fun = MemoizeJac(fun) 

jac = fun.derivative 

else: 

jac = None 

 

# set default tolerances 

if tol is not None: 

options = dict(options) 

if meth in ('hybr', 'lm'): 

options.setdefault('xtol', tol) 

elif meth in ('df-sane',): 

options.setdefault('ftol', tol) 

elif meth in ('broyden1', 'broyden2', 'anderson', 'linearmixing', 

'diagbroyden', 'excitingmixing', 'krylov'): 

options.setdefault('xtol', tol) 

options.setdefault('xatol', np.inf) 

options.setdefault('ftol', np.inf) 

options.setdefault('fatol', np.inf) 

 

if meth == 'hybr': 

sol = _root_hybr(fun, x0, args=args, jac=jac, **options) 

elif meth == 'lm': 

sol = _root_leastsq(fun, x0, args=args, jac=jac, **options) 

elif meth == 'df-sane': 

_warn_jac_unused(jac, method) 

sol = _root_df_sane(fun, x0, args=args, callback=callback, 

**options) 

elif meth in ('broyden1', 'broyden2', 'anderson', 'linearmixing', 

'diagbroyden', 'excitingmixing', 'krylov'): 

_warn_jac_unused(jac, method) 

sol = _root_nonlin_solve(fun, x0, args=args, jac=jac, 

_method=meth, _callback=callback, 

**options) 

else: 

raise ValueError('Unknown solver %s' % method) 

 

return sol 

 

 

def _warn_jac_unused(jac, method): 

if jac is not None: 

warn('Method %s does not use the jacobian (jac).' % (method,), 

RuntimeWarning) 

 

 

def _root_leastsq(func, x0, args=(), jac=None, 

col_deriv=0, xtol=1.49012e-08, ftol=1.49012e-08, 

gtol=0.0, maxiter=0, eps=0.0, factor=100, diag=None, 

**unknown_options): 

""" 

Solve for least squares with Levenberg-Marquardt 

 

Options 

------- 

col_deriv : bool 

non-zero to specify that the Jacobian function computes derivatives 

down the columns (faster, because there is no transpose operation). 

ftol : float 

Relative error desired in the sum of squares. 

xtol : float 

Relative error desired in the approximate solution. 

gtol : float 

Orthogonality desired between the function vector and the columns 

of the Jacobian. 

maxiter : int 

The maximum number of calls to the function. If zero, then 

100*(N+1) is the maximum where N is the number of elements in x0. 

epsfcn : float 

A suitable step length for the forward-difference approximation of 

the Jacobian (for Dfun=None). If epsfcn is less than the machine 

precision, it is assumed that the relative errors in the functions 

are of the order of the machine precision. 

factor : float 

A parameter determining the initial step bound 

(``factor * || diag * x||``). Should be in interval ``(0.1, 100)``. 

diag : sequence 

N positive entries that serve as a scale factors for the variables. 

""" 

 

_check_unknown_options(unknown_options) 

x, cov_x, info, msg, ier = leastsq(func, x0, args=args, Dfun=jac, 

full_output=True, 

col_deriv=col_deriv, xtol=xtol, 

ftol=ftol, gtol=gtol, 

maxfev=maxiter, epsfcn=eps, 

factor=factor, diag=diag) 

sol = OptimizeResult(x=x, message=msg, status=ier, 

success=ier in (1, 2, 3, 4), cov_x=cov_x, 

fun=info.pop('fvec')) 

sol.update(info) 

return sol 

 

 

def _root_nonlin_solve(func, x0, args=(), jac=None, 

_callback=None, _method=None, 

nit=None, disp=False, maxiter=None, 

ftol=None, fatol=None, xtol=None, xatol=None, 

tol_norm=None, line_search='armijo', jac_options=None, 

**unknown_options): 

_check_unknown_options(unknown_options) 

 

f_tol = fatol 

f_rtol = ftol 

x_tol = xatol 

x_rtol = xtol 

verbose = disp 

if jac_options is None: 

jac_options = dict() 

 

jacobian = {'broyden1': nonlin.BroydenFirst, 

'broyden2': nonlin.BroydenSecond, 

'anderson': nonlin.Anderson, 

'linearmixing': nonlin.LinearMixing, 

'diagbroyden': nonlin.DiagBroyden, 

'excitingmixing': nonlin.ExcitingMixing, 

'krylov': nonlin.KrylovJacobian 

}[_method] 

 

if args: 

if jac: 

def f(x): 

return func(x, *args)[0] 

else: 

def f(x): 

return func(x, *args) 

else: 

f = func 

 

x, info = nonlin.nonlin_solve(f, x0, jacobian=jacobian(**jac_options), 

iter=nit, verbose=verbose, 

maxiter=maxiter, f_tol=f_tol, 

f_rtol=f_rtol, x_tol=x_tol, 

x_rtol=x_rtol, tol_norm=tol_norm, 

line_search=line_search, 

callback=_callback, full_output=True, 

raise_exception=False) 

sol = OptimizeResult(x=x) 

sol.update(info) 

return sol 

 

def _root_broyden1_doc(): 

""" 

Options 

------- 

nit : int, optional 

Number of iterations to make. If omitted (default), make as many 

as required to meet tolerances. 

disp : bool, optional 

Print status to stdout on every iteration. 

maxiter : int, optional 

Maximum number of iterations to make. If more are needed to 

meet convergence, `NoConvergence` is raised. 

ftol : float, optional 

Relative tolerance for the residual. If omitted, not used. 

fatol : float, optional 

Absolute tolerance (in max-norm) for the residual. 

If omitted, default is 6e-6. 

xtol : float, optional 

Relative minimum step size. If omitted, not used. 

xatol : float, optional 

Absolute minimum step size, as determined from the Jacobian 

approximation. If the step size is smaller than this, optimization 

is terminated as successful. If omitted, not used. 

tol_norm : function(vector) -> scalar, optional 

Norm to use in convergence check. Default is the maximum norm. 

line_search : {None, 'armijo' (default), 'wolfe'}, optional 

Which type of a line search to use to determine the step size in 

the direction given by the Jacobian approximation. Defaults to 

'armijo'. 

jac_options : dict, optional 

Options for the respective Jacobian approximation. 

alpha : float, optional 

Initial guess for the Jacobian is (-1/alpha). 

reduction_method : str or tuple, optional 

Method used in ensuring that the rank of the Broyden 

matrix stays low. Can either be a string giving the 

name of the method, or a tuple of the form ``(method, 

param1, param2, ...)`` that gives the name of the 

method and values for additional parameters. 

 

Methods available: 

- ``restart``: drop all matrix columns. Has no 

extra parameters. 

- ``simple``: drop oldest matrix column. Has no 

extra parameters. 

- ``svd``: keep only the most significant SVD 

components. 

Extra parameters: 

- ``to_retain``: number of SVD components to 

retain when rank reduction is done. 

Default is ``max_rank - 2``. 

max_rank : int, optional 

Maximum rank for the Broyden matrix. 

Default is infinity (ie., no rank reduction). 

""" 

pass 

 

def _root_broyden2_doc(): 

""" 

Options 

------- 

nit : int, optional 

Number of iterations to make. If omitted (default), make as many 

as required to meet tolerances. 

disp : bool, optional 

Print status to stdout on every iteration. 

maxiter : int, optional 

Maximum number of iterations to make. If more are needed to 

meet convergence, `NoConvergence` is raised. 

ftol : float, optional 

Relative tolerance for the residual. If omitted, not used. 

fatol : float, optional 

Absolute tolerance (in max-norm) for the residual. 

If omitted, default is 6e-6. 

xtol : float, optional 

Relative minimum step size. If omitted, not used. 

xatol : float, optional 

Absolute minimum step size, as determined from the Jacobian 

approximation. If the step size is smaller than this, optimization 

is terminated as successful. If omitted, not used. 

tol_norm : function(vector) -> scalar, optional 

Norm to use in convergence check. Default is the maximum norm. 

line_search : {None, 'armijo' (default), 'wolfe'}, optional 

Which type of a line search to use to determine the step size in 

the direction given by the Jacobian approximation. Defaults to 

'armijo'. 

jac_options : dict, optional 

Options for the respective Jacobian approximation. 

 

alpha : float, optional 

Initial guess for the Jacobian is (-1/alpha). 

reduction_method : str or tuple, optional 

Method used in ensuring that the rank of the Broyden 

matrix stays low. Can either be a string giving the 

name of the method, or a tuple of the form ``(method, 

param1, param2, ...)`` that gives the name of the 

method and values for additional parameters. 

 

Methods available: 

- ``restart``: drop all matrix columns. Has no 

extra parameters. 

- ``simple``: drop oldest matrix column. Has no 

extra parameters. 

- ``svd``: keep only the most significant SVD 

components. 

Extra parameters: 

- ``to_retain``: number of SVD components to 

retain when rank reduction is done. 

Default is ``max_rank - 2``. 

max_rank : int, optional 

Maximum rank for the Broyden matrix. 

Default is infinity (ie., no rank reduction). 

""" 

pass 

 

def _root_anderson_doc(): 

""" 

Options 

------- 

nit : int, optional 

Number of iterations to make. If omitted (default), make as many 

as required to meet tolerances. 

disp : bool, optional 

Print status to stdout on every iteration. 

maxiter : int, optional 

Maximum number of iterations to make. If more are needed to 

meet convergence, `NoConvergence` is raised. 

ftol : float, optional 

Relative tolerance for the residual. If omitted, not used. 

fatol : float, optional 

Absolute tolerance (in max-norm) for the residual. 

If omitted, default is 6e-6. 

xtol : float, optional 

Relative minimum step size. If omitted, not used. 

xatol : float, optional 

Absolute minimum step size, as determined from the Jacobian 

approximation. If the step size is smaller than this, optimization 

is terminated as successful. If omitted, not used. 

tol_norm : function(vector) -> scalar, optional 

Norm to use in convergence check. Default is the maximum norm. 

line_search : {None, 'armijo' (default), 'wolfe'}, optional 

Which type of a line search to use to determine the step size in 

the direction given by the Jacobian approximation. Defaults to 

'armijo'. 

jac_options : dict, optional 

Options for the respective Jacobian approximation. 

 

alpha : float, optional 

Initial guess for the Jacobian is (-1/alpha). 

M : float, optional 

Number of previous vectors to retain. Defaults to 5. 

w0 : float, optional 

Regularization parameter for numerical stability. 

Compared to unity, good values of the order of 0.01. 

""" 

pass 

 

def _root_linearmixing_doc(): 

""" 

Options 

------- 

nit : int, optional 

Number of iterations to make. If omitted (default), make as many 

as required to meet tolerances. 

disp : bool, optional 

Print status to stdout on every iteration. 

maxiter : int, optional 

Maximum number of iterations to make. If more are needed to 

meet convergence, ``NoConvergence`` is raised. 

ftol : float, optional 

Relative tolerance for the residual. If omitted, not used. 

fatol : float, optional 

Absolute tolerance (in max-norm) for the residual. 

If omitted, default is 6e-6. 

xtol : float, optional 

Relative minimum step size. If omitted, not used. 

xatol : float, optional 

Absolute minimum step size, as determined from the Jacobian 

approximation. If the step size is smaller than this, optimization 

is terminated as successful. If omitted, not used. 

tol_norm : function(vector) -> scalar, optional 

Norm to use in convergence check. Default is the maximum norm. 

line_search : {None, 'armijo' (default), 'wolfe'}, optional 

Which type of a line search to use to determine the step size in 

the direction given by the Jacobian approximation. Defaults to 

'armijo'. 

jac_options : dict, optional 

Options for the respective Jacobian approximation. 

 

alpha : float, optional 

initial guess for the jacobian is (-1/alpha). 

""" 

pass 

 

def _root_diagbroyden_doc(): 

""" 

Options 

------- 

nit : int, optional 

Number of iterations to make. If omitted (default), make as many 

as required to meet tolerances. 

disp : bool, optional 

Print status to stdout on every iteration. 

maxiter : int, optional 

Maximum number of iterations to make. If more are needed to 

meet convergence, `NoConvergence` is raised. 

ftol : float, optional 

Relative tolerance for the residual. If omitted, not used. 

fatol : float, optional 

Absolute tolerance (in max-norm) for the residual. 

If omitted, default is 6e-6. 

xtol : float, optional 

Relative minimum step size. If omitted, not used. 

xatol : float, optional 

Absolute minimum step size, as determined from the Jacobian 

approximation. If the step size is smaller than this, optimization 

is terminated as successful. If omitted, not used. 

tol_norm : function(vector) -> scalar, optional 

Norm to use in convergence check. Default is the maximum norm. 

line_search : {None, 'armijo' (default), 'wolfe'}, optional 

Which type of a line search to use to determine the step size in 

the direction given by the Jacobian approximation. Defaults to 

'armijo'. 

jac_options : dict, optional 

Options for the respective Jacobian approximation. 

 

alpha : float, optional 

initial guess for the jacobian is (-1/alpha). 

""" 

pass 

 

def _root_excitingmixing_doc(): 

""" 

Options 

------- 

nit : int, optional 

Number of iterations to make. If omitted (default), make as many 

as required to meet tolerances. 

disp : bool, optional 

Print status to stdout on every iteration. 

maxiter : int, optional 

Maximum number of iterations to make. If more are needed to 

meet convergence, `NoConvergence` is raised. 

ftol : float, optional 

Relative tolerance for the residual. If omitted, not used. 

fatol : float, optional 

Absolute tolerance (in max-norm) for the residual. 

If omitted, default is 6e-6. 

xtol : float, optional 

Relative minimum step size. If omitted, not used. 

xatol : float, optional 

Absolute minimum step size, as determined from the Jacobian 

approximation. If the step size is smaller than this, optimization 

is terminated as successful. If omitted, not used. 

tol_norm : function(vector) -> scalar, optional 

Norm to use in convergence check. Default is the maximum norm. 

line_search : {None, 'armijo' (default), 'wolfe'}, optional 

Which type of a line search to use to determine the step size in 

the direction given by the Jacobian approximation. Defaults to 

'armijo'. 

jac_options : dict, optional 

Options for the respective Jacobian approximation. 

 

alpha : float, optional 

Initial Jacobian approximation is (-1/alpha). 

alphamax : float, optional 

The entries of the diagonal Jacobian are kept in the range 

``[alpha, alphamax]``. 

""" 

pass 

 

def _root_krylov_doc(): 

""" 

Options 

------- 

nit : int, optional 

Number of iterations to make. If omitted (default), make as many 

as required to meet tolerances. 

disp : bool, optional 

Print status to stdout on every iteration. 

maxiter : int, optional 

Maximum number of iterations to make. If more are needed to 

meet convergence, `NoConvergence` is raised. 

ftol : float, optional 

Relative tolerance for the residual. If omitted, not used. 

fatol : float, optional 

Absolute tolerance (in max-norm) for the residual. 

If omitted, default is 6e-6. 

xtol : float, optional 

Relative minimum step size. If omitted, not used. 

xatol : float, optional 

Absolute minimum step size, as determined from the Jacobian 

approximation. If the step size is smaller than this, optimization 

is terminated as successful. If omitted, not used. 

tol_norm : function(vector) -> scalar, optional 

Norm to use in convergence check. Default is the maximum norm. 

line_search : {None, 'armijo' (default), 'wolfe'}, optional 

Which type of a line search to use to determine the step size in 

the direction given by the Jacobian approximation. Defaults to 

'armijo'. 

jac_options : dict, optional 

Options for the respective Jacobian approximation. 

 

rdiff : float, optional 

Relative step size to use in numerical differentiation. 

method : {'lgmres', 'gmres', 'bicgstab', 'cgs', 'minres'} or function 

Krylov method to use to approximate the Jacobian. 

Can be a string, or a function implementing the same 

interface as the iterative solvers in 

`scipy.sparse.linalg`. 

 

The default is `scipy.sparse.linalg.lgmres`. 

inner_M : LinearOperator or InverseJacobian 

Preconditioner for the inner Krylov iteration. 

Note that you can use also inverse Jacobians as (adaptive) 

preconditioners. For example, 

 

>>> jac = BroydenFirst() 

>>> kjac = KrylovJacobian(inner_M=jac.inverse). 

 

If the preconditioner has a method named 'update', it will 

be called as ``update(x, f)`` after each nonlinear step, 

with ``x`` giving the current point, and ``f`` the current 

function value. 

inner_tol, inner_maxiter, ... 

Parameters to pass on to the "inner" Krylov solver. 

See `scipy.sparse.linalg.gmres` for details. 

outer_k : int, optional 

Size of the subspace kept across LGMRES nonlinear 

iterations. 

 

See `scipy.sparse.linalg.lgmres` for details. 

""" 

pass