continuate.newton

Basic Newton methods using Krylov subspace method.

class continuate.newton.Hessian(func, x0, jacobi_alpha, **opt)[source]

Bases: object

Calculate the deviation from linear approximation

Methods

__call__(v)
deviation(v)
trusted_region(v, eps[, r0, p]) Estimate the trusted region in which the deviation is smaller than eps.
deviation(v)[source]
logger = <continuate.misc.Logger object>
trusted_region(v, eps, r0=None, p=2)[source]

Estimate the trusted region in which the deviation is smaller than eps.

Parameters:

eps : float

Destination value of deviation

p : float, optional (default=2)

Iteration will end if the deviation is in [eps/p, eps*p].

Returns:

r : float

radius of the trusted region

continuate.newton.Jacobi(func, x0, jacobi_alpha, fx=None, **opt)[source]

Jacobi oprator \(DF(x0)\), where

\[DF(x0)dx = (F(x0 + \alpha dx / |dx|) - F(x0))/(\alpha/|dx|)\]
Parameters:

func : numpy.array -> numpy.array

A function to be differentiated

x0 : numpy.array

A position where the Jacobian is evaluated

fx : numpy.array, optional

func(x0)

Examples

>>> import numpy as np
>>> f = lambda x: np.array([x[1]**2, x[0]**2])
>>> x0 = np.array([1, 2])
>>> J = Jacobi(f, x0, 1e-7)
continuate.newton.default_options = {'newton_krylov_tol_ratio': 0.1, 'newton_maxiter': 100, 'jacobi_alpha': 1e-07, 'newton_tol': 1e-07, 'hook_maxiter': 100, 'hook_tol': 0.1, 'trusted_region': 0.1}

default values of options

You can get these values through continuate.get_default_options()

Parameters:

newton_tol : float

Tolerrance of Newton step

newton_krylov_tol_ratio : float

relative tolerrance of Krylov method (GMRES)

newton_maxiter : int

Max iteration number of Newton method

jacobi_alpha : float

Infinitesimal for calculating Jacobian matrix

trusted_region : float

Radius in model trusted region approach

hook_maxiter : int

Max iteration number of hook step

hook_tol : int

Relative tolerance of hook step iteration

continuate.newton.hook_step(A, b, trusted_region, hook_maxiter, hook_tol, nu=0, **opt)[source]

optimal hook step based on trusted region approach

Return \(x\) which minimizes \(|Ax-b|\) under the constraint \(|x| < r\), where \(r\) denotes the radius of trusted region.

Parameters:

A : numpy.array (2d)

\(A\)

b : numpy.array (1d)

\(b\)

nu : float, optional (default=0.0)

initial value of Lagrange multiplier

Returns:

x : numpy.array

argmin of x

mu : float

Lagrange multiplier

References

continuate.newton.newton_krylov(func, x, *args, **kwds)[source]

solve multi-dimensional equation \(F(x) = 0\) using Newton-Krylov method.

Parameters:

func : numpy.array -> numpy.array

\(F\) in the problem \(F(x) = 0\)

x0 : numpy.array

initial guess of the solution

Returns:

x : numpy.array

The solution \(x\) satisfies \(F(x)=0\)

Examples

>>> import continuate
>>> opt = continuate.get_default_options()
>>> opt["newton_tol"] = 1e-7
>>> f = lambda x : (x-1)**2
>>> x0 = np.array([1.1, 0.9])
>>> x = newton_krylov(f, x0, **opt)
>>> np.allclose(f(x), np.zeros_like(x), atol=1e-7)
True
continuate.newton.newton_krylov_hook(func, x, *args, **kwds)[source]

Solve multi-dimensional nonlinear equation \(F(x)=0\)

Parameters:

func : numpy.array -> numpy.array

The function \(F\) of the problem.

x0 : numpy.array

Initial guess of the solution

opt : dict

Options for calculation. Please generate by continuate.get_default_options()

Examples

>>> import continuate
>>> opt = continuate.get_default_options()
>>> opt["newton_tol"] = 1e-7
>>> f = lambda x : (x-1)**2
>>> x0 = np.array([1.1, 0.9])
>>> x = newton_krylov_hook(f, x0, **opt)
>>> np.allclose(f(x), np.zeros_like(x), atol=1e-7)
True
continuate.newton.newton_krylov_hook_gen(func, x, *args, **kwds)[source]

Generator of Newton-Krylov-hook iteration

Yields:

x : numpy.array

\(x_n\)

residual : float

\(|F(x_n)|\)

fx : numpy.array

\(F(x_n)\)