tlm_adjoint.firedrake.hessian_system

Module Contents

class tlm_adjoint.firedrake.hessian_system.HessianSystem(H, M, *, nullspace=None, comm=None)

Defines a linear system involving a Hessian matrix,

\[H u = b.\]
Parameters:
solve(u, b, **kwargs)

Solve a linear system involving a Hessian matrix,

\[H u = b.\]
Parameters:

Remaining arguments are handed to the base class System.solve() method.

tlm_adjoint.firedrake.hessian_system.hessian_eigendecompose(H, m, B_inv_action, B_action, *, nullspace=None, problem_type=None, pre_callback=None, correct_eigenvectors=True, **kwargs)

Interface with SLEPc via slepc4py, for the matrix free solution of generalized eigenproblems

\[H v = \lambda B^{-1} v,\]

where \(H\) is a Hessian matrix.

Despite the notation \(B^{-1}\) may be singular, defining an inverse operator only on an appropriate subspace.

Parameters:

Remaining keyword arguments are passed to eigendecompose().

tlm_adjoint.firedrake.hessian_system.B_inv_orthonormality_test(V, B_inv_action)

Check for \(B^{-1}\)-orthonormality.

Requires real spaces.

Parameters:
Returns:

A tuple (max_diagonal_error_norm, max_off_diagonal_error_norm) with

  • max_diagonal_error_norm: The maximum \(B^{-1}\) normalization error magnitude.

  • max_diagonal_error_norm: The maximum \(B^{-1}\) orthogonality error magnitude.

tlm_adjoint.firedrake.hessian_system.hessian_eigendecomposition_pc(B_action, Lam, V)

Construct a Hessian matrix preconditioner using a partial spectrum generalized eigendecomposition. Assumes that the Hessian matrix consists of two terms

\[H = R^{-1} + B^{-1},\]

where \(R\) and \(B\) are symmetric.

Assumes real spaces. Despite the notation \(R^{-1}\) and \(B^{-1}\) (and later \(H^{-1}\)) may be singular, defining inverse operators only on an appropriate subspace. \(B\) is assumed to define a symmetric positive definite operator on that subspace.

The approximation is defined via

\[H^{-1} \approx B + V \Lambda \left( I + \Lambda \right)^{-1} V^T\]

where

\[R^{-1} V = B^{-1} V \Lambda,\]

and where \(\Lambda\) is a diagonal matrix and \(V\) has \(B^{-1}\)-orthonormal columns, \(V^T B^{-1} V = I\).

This low rank update approximation for the Hessian matrix inverse is described in

  • Tobin Isaac, Noemi Petra, Georg Stadler, and Omar Ghattas, ‘Scalable and efficient algorithms for the propagation of uncertainty from data through inference to prediction for large-scale problems, with application to flow of the Antarctic ice sheet’, Journal of Computational Physics, 296, pp. 348–368, 2015, doi: 10.1016/j.jcp.2015.04.047

See in particular their equation (20).

Parameters:
Returns:

A callable suitable for use as the pc_fn argument to HessianSystem.solve().