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:
H – A
Hessian
defining \(H\).M – A
firedrake.function.Function
orfiredrake.cofunction.Cofunction
, or aSequence
offiredrake.function.Function
orfiredrake.cofunction.Cofunction
objects, defining the control.nullspace – A
Nullspace
or aSequence
ofNullspace
objects defining the nullspace and left nullspace of the Hessian matrix. None indicates aNoneNullspace
.comm – A communicator.
- solve(u, b, **kwargs)
Solve a linear system involving a Hessian matrix,
\[H u = b.\]- Parameters:
u – A
firedrake.function.Function
orfiredrake.cofunction.Cofunction
, or aSequence
offiredrake.function.Function
orfiredrake.cofunction.Cofunction
objects, defining the solution \(u\).b – A
firedrake.function.Function
orfiredrake.cofunction.Cofunction
, or aSequence
offiredrake.function.Function
orfiredrake.cofunction.Cofunction
objects, defining the conjugate of the right-hand-side \(b\).
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:
H – A
Hessian
.m – A
firedrake.function.Function
orfiredrake.cofunction.Cofunction
defining the control.B_inv_action – A callable accepting a
firedrake.function.Function
orfiredrake.cofunction.Cofunction
defining \(v\) and computing the conjugate of the action of \(B^{-1}\) on \(v\), returning the result as afiredrake.function.Function
orfiredrake.cofunction.Cofunction
.B_action – A callable accepting a
firedrake.function.Function
orfiredrake.cofunction.Cofunction
defining \(v\) and computing the action of \(B\) on the conjugate of \(v\), returning the result as afiredrake.function.Function
orfiredrake.cofunction.Cofunction
.nullspace – A
Nullspace
defining the nullspace and left nullspace of \(H\) and \(B^{-1}\).problem_type – The eigenproblem type – see slepc4py.SLEPc.EPS.ProblemType. Defaults to slepc4py.SLEPc.EPS.ProblemType.GHEP in the real case and slepc4py.SLEPc.EPS.ProblemType.GNHEP in the complex case.
pre_callback – A callable accepting a single slepc4py.SLEPc.EPS argument. Used for detailed manual configuration. Called after all other configuration options are set, but before the slepc4py.SLEPc.EPS.setUp method is called.
correct_eigenvectors – Whether to apply a nullspace correction to the eigenvectors.
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:
B_inv_action – A callable accepting a
firedrake.function.Function
orfiredrake.cofunction.Cofunction
defining \(v\) and computing the action of \(B^{-1}\) on \(v\), returning the result as afiredrake.function.Function
orfiredrake.cofunction.Cofunction
.V – A
Sequence
offiredrake.function.Function
orfiredrake.cofunction.Cofunction
objects to test for \(B^{-1}\)-orthonormality.
- Returns:
A
tuple
(max_diagonal_error_norm, max_off_diagonal_error_norm) withmax_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:
B_action – A callable accepting a
firedrake.function.Function
orfiredrake.cofunction.Cofunction
defining \(v\) and computing the action of \(B\) on \(v\), returning the result as afiredrake.function.Function
orfiredrake.cofunction.Cofunction
.Lam – A
Sequence
defining the diagonal of \(\Lambda\).V – A
Sequence
offiredrake.function.Function
orfiredrake.cofunction.Cofunction
objects defining the columns of \(V\).
- Returns:
A callable suitable for use as the pc_fn argument to
HessianSystem.solve()
.