Subproblem Solvers

This module provides the machinery to solve 1- and N-dimensional trust-region subproblems.

Functions Summary

dsecular(lam, w, eigvals, eigvecs, delta)

Derivative of the secular equation

dslam(lam, w, eigvals, eigvecs)

Computes the derivative of the solution \(s(\lambda)\) with respect to lambda, where \(s\) is the subproblem solution according to

get_1d_trust_region_boundary_solution(B, g, ...)

quadratic_form(Q, p, x)

Computes the quadratic form \(x^TQx + x^Tp\)

secular(lam, w, eigvals, eigvecs, delta)

Secular equation

slam(lam, w, eigvals, eigvecs)

Computes the solution \(s(\lambda)\) as subproblem solution according to

solve_1d_trust_region_subproblem(B, g, s, ...)

Solves the special case of a one-dimensional subproblem

solve_nd_trust_region_subproblem(B, g, delta)

This function exactly solves the n-dimensional subproblem.

Functions

fides.subproblem.dsecular(lam, w, eigvals, eigvecs, delta)[source]

Derivative of the secular equation

\(\phi(\lambda) = \frac{1}{||s||} - \frac{1}{\Delta}\)

with respect to \(\lambda\)

Parameters:
  • lam (float) – \(\lambda\)

  • w (numpy.ndarray) – precomputed eigenvector coefficients for -g

  • eigvals (numpy.ndarray) – precomputed eigenvalues of B

  • eigvecs (numpy.ndarray) – precomputed eigenvectors of B

  • delta (float) – trust region radius \(\Delta\)

Returns:

\(\frac{\partial \phi(\lambda)}{\partial \lambda}\)

fides.subproblem.dslam(lam, w, eigvals, eigvecs)[source]

Computes the derivative of the solution \(s(\lambda)\) with respect to lambda, where \(s\) is the subproblem solution according to

\(-(B + \lambda I)s = g\)

Parameters:
Returns:

\(\frac{\partial s(\lambda)}{\partial \lambda}\)

fides.subproblem.get_1d_trust_region_boundary_solution(B, g, s, s0, delta)[source]
fides.subproblem.quadratic_form(Q, p, x)[source]

Computes the quadratic form \(x^TQx + x^Tp\)

Parameters:
Return type:

float

Returns:

Value of form

fides.subproblem.secular(lam, w, eigvals, eigvecs, delta)[source]

Secular equation

\(\phi(\lambda) = \frac{1}{||s||} - \frac{1}{\Delta}\)

Subproblem solutions are given by the roots of this equation

Parameters:
  • lam (float) – \(\lambda\)

  • w (numpy.ndarray) – precomputed eigenvector coefficients for -g

  • eigvals (numpy.ndarray) – precomputed eigenvalues of B

  • eigvecs (numpy.ndarray) – precomputed eigenvectors of B

  • delta (float) – trust region radius \(\Delta\)

Returns:

\(\phi(\lambda)\)

fides.subproblem.slam(lam, w, eigvals, eigvecs)[source]

Computes the solution \(s(\lambda)\) as subproblem solution according to

\(-(B + \lambda I)s = g\)

Parameters:
Return type:

numpy.ndarray

Returns:

\(s(\lambda)\)

fides.subproblem.solve_1d_trust_region_subproblem(B, g, s, delta, s0)[source]

Solves the special case of a one-dimensional subproblem

Parameters:
  • B (numpy.ndarray) – Hessian of the quadratic subproblem

  • g (numpy.ndarray) – Gradient of the quadratic subproblem

  • s (numpy.ndarray) – Vector defining the one-dimensional search direction

  • delta (float) – Norm boundary for the solution of the quadratic subproblem

  • s0 (numpy.ndarray) – reference point from where search is started, also counts towards norm of step

Return type:

numpy.ndarray

Returns:

Proposed step-length

fides.subproblem.solve_nd_trust_region_subproblem(B, g, delta, logger=None)[source]

This function exactly solves the n-dimensional subproblem.

\(argmin_s\{s^T B s + s^T g = 0: ||s|| <= \Delta, s \in \mathbb{ R}^n\}\)

The solution is characterized by the equation \(-(B + \lambda I)s = g\). If B is positive definite, the solution can be obtained by \(\lambda = 0`$\) if \(Bs = -g\) satisfies \(||s|| <= \Delta\). If B is indefinite or \(Bs = -g\) satisfies \(||s|| > \Delta\) and an appropriate \(\lambda\) has to be identified via 1D rootfinding of the secular equation

\(\phi(\lambda) = \frac{1}{||s(\lambda)||} - \frac{1}{\Delta} = 0\)

with \(s(\lambda)\) computed according to an eigenvalue decomposition of B. The eigenvalue decomposition, although being more expensive than a cholesky decomposition, has the advantage that eigenvectors are invariant to changes in \(\lambda\) and eigenvalues are linear in \(\lambda\), so factorization only has to be performed once. We perform the linesearch via Newton’s algorithm and Brent-Q as fallback. The hard case is treated separately and serves as general fallback.

Parameters:
Return type:

typing.Tuple[numpy.ndarray, str]

Returns:

s: Selected step, step_type: Type of solution that was obtained