fdvar.preconditioners package

Submodules

fdvar.preconditioners.allatonce module

class fdvar.preconditioners.allatonce.AllAtOnceRFGaussSeidelPC

Bases: PCBase

Python preconditioner to approximate the inverse of the tangent linear or adjoint model of an AllAtOnceReducedFunctional.

The tangent linear \(L\) (or adjoint \(L^{T}\)) is block lower (upper) triangular so can be solved exactly with forward (backward) substitution, i.e. block Gauss-Seidel.

\[\begin{split}L = \begin{pmatrix} I & 0 & 0 & 0 \\ -M & I & 0 & 0 \\ 0 & -M & I & 0 \\ 0 & 0 & -M & I \\ \end{pmatrix}\end{split}\]

The default is to use \(M\) in the Gauss-Seidel iteration, in which case this PC is an exact solve. Two approximations of \(M\) commonly used in weak-constraint 4DVar are also available:

  1. \(M=I\). Approximating the propagator with identity is equivalent to computing an inclusive sum of the right hand side.

  2. \(M=0\). This “approximation” is equivalent to -pc_type none.

Use the -pc_aaogs_type option to select which approximation to use.

PETSc Options

  • -pc_aaogs_type (model|identity|zero) - Whether to use \(M\) or to approximate it with \(I\) or \(0\).

Notes

The Mat for this PC must be a ReducedFunctionalMat for an AllAtOnceReducedFunctional.

See also

AllAtOnceReducedFunctional, WC4DVarSchurPC

prefix = 'aaogs_'

The options prefix of this PC.

needs_python_pmat = True

Set this to False if the P matrix needs to be Python (matfree).

initialize(pc: PC)

Initialize any state in this preconditioner.

This method is only called on the first time that the setUp method is called.

apply(pc: PC, x: Vec, y: Vec)

Apply the preconditioner to x, putting the result in y.

Both x and y are PETSc Vecs, y is not guaranteed to be zero on entry.

applyTranspose(pc: PC, x: Vec, y: Vec)

Apply the preconditioner transpose to x, putting the result in y.

Both x and y are PETSc Vecs, y is not guaranteed to be zero on entry.

apply_tlm(pc: PC, x: Vec, y: Vec)
apply_adjoint(pc: PC, x: Vec, y: Vec)
update(pc: PC)

Update any state in this preconditioner.

This method is called the on second and later times that the setUp method is called.

This method is not needed for all preconditioners and can often be a no-op.

view(pc: PC, viewer: Viewer | None = None)

Write a basic description of this PC.

fdvar.preconditioners.wcsaddle module

fdvar.preconditioners.wcsaddle.getSubWC4DVarSaddleMat(mat: Mat, sub: str | None = None)

Return a sub matrix of the saddle point MatNest. Options are ‘D’, ‘R’, ‘L’, ‘LT’, ‘H’, ‘HT’, or None to return all sub matrices.

Parameters:
  • mat – The MatNest for the saddle point system returned by WC4DVarSaddleMat.

  • sub – Which sub matrix to return.

Returns:

The sub Mat requested or a tuple of all sub Mats.

Return type:

tuple[Mat] | Mat

fdvar.preconditioners.wcsaddle.WC4DVarSaddleMat(Jhat: WC4DVarReducedFunctional)

Mat of type MatNest for the saddle point formulation of Weak Constraint 4DVar.

\[\begin{split}A = \begin{pmatrix} D & 0 & L \\ 0 & R & H \\ L^{T} & H^{T} & 0 \end{pmatrix}\end{split}\]

where \(L\) is the all-at-once system, and H, D, and R are the block-diagonal matrices with the observation operators, observation error covariances, and model error covariances at each observation time respectively.

Parameters:

Jhat – The reduced functional to construct the saddle point matrix from.

Returns:

The 3x3 MatNest for the saddle point system.

Return type:

Mat

Raises:

TypeError : – If Jhat is not a WC4DVarReducedFunctional.

fdvar.preconditioners.wcsaddle.WC4DVarSaddleKSP(Jhat: WC4DVarReducedFunctional, Jphat: WC4DVarReducedFunctional | None = None, *, solver_parameters: dict | None = None, options_prefix: str | None = None)

A KSP for the saddle point formulation of Weak Constraint 4DVar.

Parameters:
  • Jhat – The reduced functional to construct the WC4DVarSaddleMat.

  • Jphat – The reduced functional to construct a WC4DVarSaddleMat for the Pmat operator of the KSP. If not provided then Jhat is used for both Amat and Pmat.

  • solver_parameters – PETSc options for the KSP.

  • options_prefix – Options prefix for the KSP.

Returns:

The KSP for the saddle point system.

Return type:

KSP

Raises:

TypeError : – If Jhat is not a WC4DVarReducedFunctional.

class fdvar.preconditioners.wcsaddle.WC4DVarSaddlePC

Bases: PCBase

Preconditioner for Weak Constraint 4DVar using the saddle point formulation.

\[\begin{split}\begin{pmatrix} D & 0 & L \\ 0 & R & H \\ L^{T} & H^{T} & 0 \end{pmatrix} \begin{pmatrix} \eta \\ \lambda \\ \delta x \end{pmatrix} = \begin{pmatrix} b \\ d \\ 0 \end{pmatrix}\end{split}\]

This PC acts on a ReducedFunctionalMat for a WC4DVarReducedFunctional. It solves the larger saddle point system and returns just the \(\delta x\) part of the solution.

PETSc Options

  • -wcsaddle - Options for the KSP for the saddle point system.

Notes

The \((3, 3)\) Schur complement of the saddle point system is the WC4DVar Hessian \(L^{T}D^{-1}L+H^{T}R^{-1}H\) and is often approximated with the WC4DVarSchurPC. PCFieldsplit can be used to eliminate down to the \(\delta x\) component with the following options:

-pc_type fieldsplit
-pc_fieldsplit_type schur
-pc_fieldsplit_0_fields 0,1
-pc_fieldsplit_1_fields 2

References

Fisher M. and Gurol S., 2017: “Parallelization in the time dimension of four-dimensional variational data assimilation”. Q.J.R. Meteorol. Soc. 142: 1136–1147, DOI:10.1002/qj.2997

needs_python_amat = True

Set this to True if the A matrix needs to be Python (matfree).

needs_python_pmat = True

Set this to False if the P matrix needs to be Python (matfree).

prefix = 'wcsaddle_'

The options prefix of this PC.

initialize(pc: PC)

Initialize any state in this preconditioner.

This method is only called on the first time that the setUp method is called.

apply(pc: PC, x: Vec, y: Vec)

Apply the preconditioner to x, putting the result in y.

Both x and y are PETSc Vecs, y is not guaranteed to be zero on entry.

update(pc: PC)

Update any state in this preconditioner.

This method is called the on second and later times that the setUp method is called.

This method is not needed for all preconditioners and can often be a no-op.

view(pc: PC, viewer: Viewer | None = None)

Write a basic description of this PC.

fdvar.preconditioners.wcschur module

class fdvar.preconditioners.wcschur.WC4DVarSchurPC

Bases: PCBase

Preconditioner to approximate the inverse of the Schur complement of the saddle point formulation of the weak constraint 4DVar, which is equivalent to the Gauss-Newton Hessian of the primal WC4DVar formulation.

The exact Schur complement \(S\) and the approximation \(\tilde{S}\) that this PC applies are:

\[ \begin{align}\begin{aligned}S & = L^{T}D^{-1}L + H^{T}R^{-1}H\\\tilde{S}^{-1} & = \tilde{L}^{-T}\tilde{D}\tilde{L}^{-1}\end{aligned}\end{align} \]

where \(L\) is the all-at-once system, and H, D, and R are the block-diagonal matrices with the observation operators, observation error covariances, and model error covariances at each observation time respectively.

KSPs are created for \(\tilde{L}\) and \(\tilde{L}^{-T}\) using a ReducedFunctionalMat for the AllAtOnceReducedFunctional, and for \(\tilde{D}\) using an EnsembleBlockDiagonalMat where each block is a CovarianceMat.

PETSc Options

  • -wcschur_l - Options for the KSP for \(L\) and \(L^{T}\).

  • -wcschur_ltlm - Options solely for \(L\), e.g. monitors.

  • -wcschur_ladj - Options solely for \(L^{T}\), e.g. monitors.

  • -wcschur_d - Options for the KSP for \(D^{-1}\)

Notes

Identical solver options should be used for \(\tilde{L}\) and \(\tilde{L}^{T}\) to ensure symmetry of \(\tilde{S}^{-1}\).

The Pmat operator for this preconditioner must be a ReducedFunctionalMat for a WC4DVarReducedFunctional.

References

Fisher M. and Gurol S., 2017: “Parallelization in the time dimension of four-dimensional variational data assimilation”. Q.J.R. Meteorol. Soc. 142: 1136–1147, DOI:10.1002/qj.2997

See also

WC4DVarReducedFunctional, AllAtOnceReducedFunctional, EnsembleBlockDiagonalMat, CovarianceMat

prefix = 'wcschur_'

The options prefix of this PC.

initialize(pc: PC)

Initialize any state in this preconditioner.

This method is only called on the first time that the setUp method is called.

apply(pc: PC, x: Vec, y: Vec)

Apply the preconditioner to x, putting the result in y.

Both x and y are PETSc Vecs, y is not guaranteed to be zero on entry.

update(pc: PC)

Update any state in this preconditioner.

This method is called the on second and later times that the setUp method is called.

This method is not needed for all preconditioners and can often be a no-op.

view(pc: PC, viewer: Viewer | None = None)

Write a basic description of this PC.