Skip to content

Levenberg-Marquardt

Given a function \(F : \mathbb{R}^n \rightarrow \mathbb{R}^m\) of class \(\mathcal{C}^1\), the Levenberg-Marquardt optimizer iteratively solves:

\[ \min_{x \in \mathbb{R}^n} \frac{1}{2} ||F(x)||^2\]

Principle

Let \(x_k\) be the variable vector at iteration \(k\). Define \(m_\lambda^k : \mathbb{R}^n \rightarrow \mathbb{R}\) by

\[\begin{align} m_\lambda^k(x) &= \frac{1}{2}|| F_k + J_k (x - x_k) ||^2 + \frac{\lambda}{2} || x - x_k||^2 \\ &= \frac{1}{2}|| F_k + (J_k + \lambda I) (x - x_k) ||^2 \end{align} \]

where \(F_k = F(x_k)\) is the value of the function at point \(x_k\) and \(J_k = J(x_k)\) is the value of its Jacobian matrix at point \(x_k\).

Set \(x_{k+1} = \text{argmin}_x m_\lambda^k(x)\) for some \(\lambda = \lambda_k \in \mathbb{R}^+\)

And iterate until convergence.

Example

import mouette as M
import numpy as np
import scipy.sparse as sp

def function_to_optimize(x, a, b):
    # Rosenbrock function and gradient
    f = np.array([(a-x[0])**2 + b*(x[1] - x[0]**2)**2])
    J = sp.csr_matrix([2*(x[0]-a) - 4*b*x[0]*(x[1]-x[0]*x[0]),  2*b*(x[1]-x[0]*x[0])]).reshape((1,2))
    return f,J

optimizer = M.optimize.LevenbergMarquardt()

# adjust hyperparameters from the HP attribute
optimizer.HP.ENERGY_MIN = 1E-10
optimizer.HP.MIN_GRAD_NORM = 0. 

# Register the function to optimize
optimizer.register_function(lambda x : function_to_optimize(x, 1., 10.), name="Rosenbrock")
optimizer.run(x_init=np.random.random((2,))) # run the optimization starting from a random position

print("Local minimum found:" , optimizer.X)
print("Expected minimum :", [A,A])

Implementation

LMParameters(N_ITER_MAX=1000, ENERGY_MIN=1e-07, MIN_STEP_NORM=1e-05, MIN_GRAD_NORM=1e-05, MIN_DELTA_E=1e-05, MU_MAX=100000000.0, MU_MIN=1e-08, alpha=0.5, beta=2.0) dataclass

Hyper parameters for the Levenberg-Marquardt algorithm

LevenbergMarquardt(HP=LMParameters(), verbose=LMVerboseOptions(), **kwargs)

Bases: Logger

Levenberg-Marquardt algorithm for non-linear least-square minimization.

References

[1] https://en.wikipedia.org/wiki/Levenberg%E2%80%93Marquardt_algorithm

[2] Constrained Levenberg-Marquardt Method with global complexity bounds, Marumo et al.

Example

https://github.com/GCoiffier/mouette/blob/main/examples/optimization/levenberg_marquardt_example.py

optimize(x_init)

Alias of the 'run' method.

Parameters:

Name Type Description Default
x_init ndarray

initial values for the variables.

required

Returns:

Name Type Description
float float

final value of the energy

register_constraints(A, l=None, u=None)

Adds linear constraints to the optimization problem :

min_X ||F(X)||² s.t. l <= AX <= u

Parameters:

Name Type Description Default
cstMat spmatrix

constraint matrix, either dense or sparse

required
cstL ndarray

vector l. Defaults to None.

required
cstR ndarray

vector u. Defaults to None.

required

register_function(fun, fun_noJ=None, weight=1.0, name=None)

Adds a function to minimize.

Parameters:

Name Type Description Default
fun python callable

A function taking a single argument X (np.array) returning a vector of values F and a (sparse or dense) jacobian matrix J

required
fun_noJ python callable

The same function that avoids the computation of the jacobian for speed purposes. If not provided, will use the function provided above and ignore the jacobian.

None
weight float

real weight to be applied to the function. Defaults to 1..

1.0
name str

Name of the function in the logs. If not specified, the function will be given a default name.

None

run(x_init)

Runs the optimizer

Parameters:

Name Type Description Default
x_init ndarray

initial values for the variables.

required

Returns:

Name Type Description
float

final value of the energy