fdvar package¶
- class fdvar.SC4DVarReducedFunctional(control: Control, background: OverloadedType, background_covariance: CovarianceOperatorBase, observation_covariance: CovarianceOperatorBase, observation_error: Callable[[OverloadedType], OverloadedType], tape: Tape | None = None)¶
Bases:
AbstractReducedFunctionalReducedFunctional for strong constraint 4DVar data assimilation.
- Parameters:
control – The
Functionfor the control x_{0} at the initial condition.background_covariance – The inner product to calculate the background error functional from the background error \(x_{0} - x_{b}\). Can include the error covariance matrix.
observation_covariance – The inner product to calculate the observation error functional from the observation error \(y_{0} - \mathcal{H}_{0}(x)\). Can include the error covariance matrix. Must be provided if observation_error is provided.
observation_error – Given a state \(x\), returns the observations error \(y_{0} - \mathcal{H}_{0}(x)\) where \(y_{0}\) are the observations at the initial time and \(\mathcal{H}_{0}\) is the observation operator for the initial time.
background – The background (prior) data for the initial condition \(x_{b}\). If not provided, the value of the control will be used.
See also
- property controls¶
Return the list of controls on which this functional depends.
- property functional¶
- background()¶
- __call__(values: OverloadedType)¶
Compute the reduced functional with supplied control value.
- Parameters:
values ([pyadjoint.OverloadedType]) – If you have multiple controls this should be a list of new values for each control in the order you listed the controls to the constructor. If you have a single control it can either be a list or a single object. Each new value should have the same type as the corresponding control.
- Returns:
- The computed value. Often of type
- Return type:
- derivative(adj_input: float = 1.0, apply_riesz: bool = False)¶
Return the derivative of the functional w.r.t. the control.
Using the adjoint method, the derivative of the functional with respect to the control, around the last supplied value of the control, is computed and returned.
- Parameters:
adj_input – The adjoint value to the functional result. Required if the functional is not scalar-valued, or if the functional is an intermediate result in the computation of an outer functional.
apply_riesz – If True, apply the Riesz map of each control in order to return a primal gradient rather than a derivative in the dual space.
- Returns:
- The derivative with respect to the control.
If apply_riesz is False, should be an instance of the type dual to that of the control. If apply_riesz is True should have the same type as the control.
- Return type:
- tlm(m_dot: OverloadedType)¶
Return the action of the tangent linear model of the functional.
The tangent linear model is evaluated w.r.t. the control on a vector m_dot, around the last supplied value of the control.
- Parameters:
m_dot ([pyadjoint.OverloadedType]) – The direction in which to compute the action of the tangent linear model.
- Returns:
- The action of the tangent linear model in the direction m_dot.
Should be an instance of the same type as the functional.
- Return type:
- hessian(m_dot: OverloadedType, hessian_input: OverloadedType | None = None, evaluate_tlm: bool = True, apply_riesz: bool = False)¶
Return the action of the Hessian of the functional.
The Hessian is evaluated w.r.t. the control on a vector m_dot.
Using the second-order adjoint method, the action of the Hessian of the functional with respect to the control, around the last supplied value of the control, is computed and returned.
- Parameters:
m_dot ([pyadjoint.OverloadedType]) – The direction in which to compute the action of the Hessian.
hessian_input (pyadjoint.OverloadedType) – The Hessian value for the functional result. Required if the functional is not scalar-valued, or if the functional is an intermediate result in the computation of an outer functional.
evaluate_tlm (bool) – If True, will evaluate the tangent linear model before evaluating the Hessian. If False, assumes that the tape is already populated with the tangent linear values and does not evaluate them.
apply_riesz (bool) – If True, apply the Riesz map of each control in order to return the (primal) Riesz representer of the Hessian action.
- Returns:
- The action of the Hessian in the direction m_dot.
If apply_riesz is False, should be an instance of the type dual to that of the control. If apply_riesz is True should have the same type as the control.
- Return type:
- recording_stages(nstages=None, **stage_kwargs)¶
- class fdvar.WC4DVarReducedFunctional(control: Control, background_covariance: CovarianceOperatorBase, observation_covariance: CovarianceOperatorBase, observation_error: Callable[[OverloadedType], OverloadedType] = None, background: OverloadedType | None = None, gauss_newton: bool = False)¶
Bases:
AbstractReducedFunctionalReducedFunctional for weak constraint 4DVar data assimilation.
- Parameters:
control – The
EnsembleFunctionfor the control x_{i} at the initial condition and at the end of each observation stage.background_covariance – The inner product to calculate the background error functional from the background error \(x_{0} - x_{b}\). Can include the error covariance matrix. Only used on ensemble rank 0.
observation_covariance – The inner product to calculate the observation error functional from the observation error \(y_{0} - \mathcal{H}_{0}(x)\). Can include the error covariance matrix. Must be provided if observation_error is provided. Only used on ensemble rank 0
observation_error – Given a state \(x\), returns the observations error \(y_{0} - \mathcal{H}_{0}(x)\) where \(y_{0}\) are the observations at the initial time and \(\mathcal{H}_{0}\) is the observation operator for the initial time. Only used on ensemble rank 0. Optional.
background – The background (prior) data for the initial condition \(x_{b}\). If not provided, the value of the first subfunction on the first ensemble member of the control
EnsembleFunctionwill be used.
See also
- property controls¶
Return the list of controls on which this functional depends.
- __call__(values: OverloadedType)¶
Computes the reduced functional with supplied control value.
- Parameters:
values – If you have multiple controls this should be a list of new values for each control in the order you listed the controls to the constructor. If you have a single control it can either be a list or a single object. Each new value should have the same type as the corresponding control.
- Returns:
The computed value. Typically of instance of
pyadjoint.AdjFloat.- Return type:
- derivative(adj_input: float = 1.0, apply_riesz: bool = False)¶
Returns the derivative of the functional w.r.t. the control. Using the adjoint method, the derivative of the functional with respect to the control, around the last supplied value of the control, is computed and returned.
- Parameters:
adj_input – The adjoint input.
options – Additional options for the derivative computation.
- Returns:
The derivative with respect to the control. Should be an instance of the same type as the control.
- Return type:
- tlm(m_dot: OverloadedType)¶
Return the action of the tangent linear model of the functional.
The tangent linear model is evaluated w.r.t. the control on a vector m_dot, around the last supplied value of the control.
- Parameters:
m_dot ([pyadjoint.OverloadedType]) – The direction in which to compute the action of the tangent linear model.
- Returns:
- The action of the tangent linear model in the direction m_dot.
Should be an instance of the same type as the functional.
- Return type:
- hessian(m_dot: OverloadedType, hessian_input: OverloadedType | None = None, evaluate_tlm: bool = True, apply_riesz: bool = False)¶
Returns the action of the Hessian of the functional w.r.t. the control on a vector m_dot.
Using the second-order adjoint method, the action of the Hessian of the functional with respect to the control, around the last supplied value of the control, is computed and returned.
- Parameters:
m_dot – The direction in which to compute the action of the Hessian.
- Returns:
The action of the Hessian in the direction m_dot. Should be an instance of the same type as the control.
- Return type:
- recording_stages(sequential=True, nstages=None, **stage_kwargs)¶
- class fdvar.AllAtOnceReducedFunctional(functional, control, propagator_rfs, background=None)¶
Bases:
AbstractReducedFunctional- property controls¶
Return the list of controls on which this functional depends.
- __call__(values: OverloadedType)¶
Compute the reduced functional with supplied control value.
- Parameters:
values ([pyadjoint.OverloadedType]) – If you have multiple controls this should be a list of new values for each control in the order you listed the controls to the constructor. If you have a single control it can either be a list or a single object. Each new value should have the same type as the corresponding control.
- Returns:
- The computed value. Often of type
- Return type:
- derivative(adj_input: float = 1.0, apply_riesz: bool = False)¶
Return the derivative of the functional w.r.t. the control.
Using the adjoint method, the derivative of the functional with respect to the control, around the last supplied value of the control, is computed and returned.
- Parameters:
adj_input – The adjoint value to the functional result. Required if the functional is not scalar-valued, or if the functional is an intermediate result in the computation of an outer functional.
apply_riesz – If True, apply the Riesz map of each control in order to return a primal gradient rather than a derivative in the dual space.
- Returns:
- The derivative with respect to the control.
If apply_riesz is False, should be an instance of the type dual to that of the control. If apply_riesz is True should have the same type as the control.
- Return type:
- tlm(m_dot: OverloadedType)¶
Return the action of the tangent linear model of the functional.
The tangent linear model is evaluated w.r.t. the control on a vector m_dot, around the last supplied value of the control.
- Parameters:
m_dot ([pyadjoint.OverloadedType]) – The direction in which to compute the action of the tangent linear model.
- Returns:
- The action of the tangent linear model in the direction m_dot.
Should be an instance of the same type as the functional.
- Return type:
- hessian(m_dot: OverloadedType, hessian_input: OverloadedType = None, evaluate_tlm: bool = True, apply_riesz: bool = False)¶
Return the action of the Hessian of the functional.
The Hessian is evaluated w.r.t. the control on a vector m_dot.
Using the second-order adjoint method, the action of the Hessian of the functional with respect to the control, around the last supplied value of the control, is computed and returned.
- Parameters:
m_dot ([pyadjoint.OverloadedType]) – The direction in which to compute the action of the Hessian.
hessian_input (pyadjoint.OverloadedType) – The Hessian value for the functional result. Required if the functional is not scalar-valued, or if the functional is an intermediate result in the computation of an outer functional.
evaluate_tlm (bool) – If True, will evaluate the tangent linear model before evaluating the Hessian. If False, assumes that the tape is already populated with the tangent linear values and does not evaluate them.
apply_riesz (bool) – If True, apply the Riesz map of each control in order to return the (primal) Riesz representer of the Hessian action.
- Returns:
- The action of the Hessian in the direction m_dot.
If apply_riesz is False, should be an instance of the type dual to that of the control. If apply_riesz is True should have the same type as the control.
- Return type:
- class fdvar.AllAtOnceRFGaussSeidelPC¶
Bases:
PCBasePython 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:
\(M=I\). Approximating the propagator with identity is equivalent to computing an inclusive sum of the right hand side.
\(M=0\). This “approximation” is equivalent to
-pc_type none.
Use the
-pc_aaogs_typeoption 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
Matfor this PC must be aReducedFunctionalMatfor anAllAtOnceReducedFunctional.See also
- 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
setUpmethod 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.
- class fdvar.WC4DVarSchurPC¶
Bases:
PCBasePreconditioner 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
ReducedFunctionalMatfor theAllAtOnceReducedFunctional, and for \(\tilde{D}\) using anEnsembleBlockDiagonalMatwhere each block is aCovarianceMat.PETSc Options
-wcschur_l- Options for theKSPfor \(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 theKSPfor \(D^{-1}\)
Notes
Identical solver options should be used for \(\tilde{L}\) and \(\tilde{L}^{T}\) to ensure symmetry of \(\tilde{S}^{-1}\).
The
Pmatoperator for this preconditioner must be aReducedFunctionalMatfor aWC4DVarReducedFunctional.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
setUpmethod 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.
- class fdvar.WC4DVarSaddlePC¶
Bases:
PCBasePreconditioner 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
ReducedFunctionalMatfor aWC4DVarReducedFunctional. It solves the larger saddle point system and returns just the \(\delta x\) part of the solution.PETSc Options
-wcsaddle- Options for theKSPfor 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
See also
- 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
setUpmethod 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.
Subpackages¶
- fdvar.adjoint package
EnsembleAdjVecEnsembleTransformEnsembleReduceEnsembleBcastEnsembleShiftReducedFunctionalPipeline- Submodules
- fdvar.adjoint.ensemble_adjvec module
- fdvar.adjoint.ensemble_operations module
reductionbroadcastEnsembleCommunicationBaseEnsembleCommunicationBase.controlsEnsembleCommunicationBase.functionalEnsembleCommunicationBase.srcEnsembleCommunicationBase.dstEnsembleCommunicationBase.ensembleEnsembleCommunicationBase.rootEnsembleCommunicationBase.__call__EnsembleCommunicationBase.tlmEnsembleCommunicationBase.derivativeEnsembleCommunicationBase.hessian
EnsembleReduceEnsembleBcastEnsembleTransformEnsembleShift
- fdvar.adjoint.reduced_functional_pipeline module
- fdvar.preconditioners package
- Submodules
- fdvar.preconditioners.allatonce module
AllAtOnceRFGaussSeidelPCAllAtOnceRFGaussSeidelPC.prefixAllAtOnceRFGaussSeidelPC.needs_python_pmatAllAtOnceRFGaussSeidelPC.initializeAllAtOnceRFGaussSeidelPC.applyAllAtOnceRFGaussSeidelPC.applyTransposeAllAtOnceRFGaussSeidelPC.apply_tlmAllAtOnceRFGaussSeidelPC.apply_adjointAllAtOnceRFGaussSeidelPC.updateAllAtOnceRFGaussSeidelPC.view
- fdvar.preconditioners.wcsaddle module
- fdvar.preconditioners.wcschur module
Submodules¶
fdvar.allatonce_reduced_functional module¶
- class fdvar.allatonce_reduced_functional.AllAtOnceReducedFunctional(functional, control, propagator_rfs, background=None)¶
Bases:
AbstractReducedFunctional- property controls¶
Return the list of controls on which this functional depends.
- __call__(values: OverloadedType)¶
Compute the reduced functional with supplied control value.
- Parameters:
values ([pyadjoint.OverloadedType]) – If you have multiple controls this should be a list of new values for each control in the order you listed the controls to the constructor. If you have a single control it can either be a list or a single object. Each new value should have the same type as the corresponding control.
- Returns:
- The computed value. Often of type
- Return type:
- derivative(adj_input: float = 1.0, apply_riesz: bool = False)¶
Return the derivative of the functional w.r.t. the control.
Using the adjoint method, the derivative of the functional with respect to the control, around the last supplied value of the control, is computed and returned.
- Parameters:
adj_input – The adjoint value to the functional result. Required if the functional is not scalar-valued, or if the functional is an intermediate result in the computation of an outer functional.
apply_riesz – If True, apply the Riesz map of each control in order to return a primal gradient rather than a derivative in the dual space.
- Returns:
- The derivative with respect to the control.
If apply_riesz is False, should be an instance of the type dual to that of the control. If apply_riesz is True should have the same type as the control.
- Return type:
- tlm(m_dot: OverloadedType)¶
Return the action of the tangent linear model of the functional.
The tangent linear model is evaluated w.r.t. the control on a vector m_dot, around the last supplied value of the control.
- Parameters:
m_dot ([pyadjoint.OverloadedType]) – The direction in which to compute the action of the tangent linear model.
- Returns:
- The action of the tangent linear model in the direction m_dot.
Should be an instance of the same type as the functional.
- Return type:
- hessian(m_dot: OverloadedType, hessian_input: OverloadedType = None, evaluate_tlm: bool = True, apply_riesz: bool = False)¶
Return the action of the Hessian of the functional.
The Hessian is evaluated w.r.t. the control on a vector m_dot.
Using the second-order adjoint method, the action of the Hessian of the functional with respect to the control, around the last supplied value of the control, is computed and returned.
- Parameters:
m_dot ([pyadjoint.OverloadedType]) – The direction in which to compute the action of the Hessian.
hessian_input (pyadjoint.OverloadedType) – The Hessian value for the functional result. Required if the functional is not scalar-valued, or if the functional is an intermediate result in the computation of an outer functional.
evaluate_tlm (bool) – If True, will evaluate the tangent linear model before evaluating the Hessian. If False, assumes that the tape is already populated with the tangent linear values and does not evaluate them.
apply_riesz (bool) – If True, apply the Riesz map of each control in order to return the (primal) Riesz representer of the Hessian action.
- Returns:
- The action of the Hessian in the direction m_dot.
If apply_riesz is False, should be an instance of the type dual to that of the control. If apply_riesz is True should have the same type as the control.
- Return type:
fdvar.sc4dvar_reduced_functional module¶
- class fdvar.sc4dvar_reduced_functional.SC4DVarReducedFunctional(control: Control, background: OverloadedType, background_covariance: CovarianceOperatorBase, observation_covariance: CovarianceOperatorBase, observation_error: Callable[[OverloadedType], OverloadedType], tape: Tape | None = None)¶
Bases:
AbstractReducedFunctionalReducedFunctional for strong constraint 4DVar data assimilation.
- Parameters:
control – The
Functionfor the control x_{0} at the initial condition.background_covariance – The inner product to calculate the background error functional from the background error \(x_{0} - x_{b}\). Can include the error covariance matrix.
observation_covariance – The inner product to calculate the observation error functional from the observation error \(y_{0} - \mathcal{H}_{0}(x)\). Can include the error covariance matrix. Must be provided if observation_error is provided.
observation_error – Given a state \(x\), returns the observations error \(y_{0} - \mathcal{H}_{0}(x)\) where \(y_{0}\) are the observations at the initial time and \(\mathcal{H}_{0}\) is the observation operator for the initial time.
background – The background (prior) data for the initial condition \(x_{b}\). If not provided, the value of the control will be used.
See also
- property controls¶
Return the list of controls on which this functional depends.
- property functional¶
- background()¶
- __call__(values: OverloadedType)¶
Compute the reduced functional with supplied control value.
- Parameters:
values ([pyadjoint.OverloadedType]) – If you have multiple controls this should be a list of new values for each control in the order you listed the controls to the constructor. If you have a single control it can either be a list or a single object. Each new value should have the same type as the corresponding control.
- Returns:
- The computed value. Often of type
- Return type:
- derivative(adj_input: float = 1.0, apply_riesz: bool = False)¶
Return the derivative of the functional w.r.t. the control.
Using the adjoint method, the derivative of the functional with respect to the control, around the last supplied value of the control, is computed and returned.
- Parameters:
adj_input – The adjoint value to the functional result. Required if the functional is not scalar-valued, or if the functional is an intermediate result in the computation of an outer functional.
apply_riesz – If True, apply the Riesz map of each control in order to return a primal gradient rather than a derivative in the dual space.
- Returns:
- The derivative with respect to the control.
If apply_riesz is False, should be an instance of the type dual to that of the control. If apply_riesz is True should have the same type as the control.
- Return type:
- tlm(m_dot: OverloadedType)¶
Return the action of the tangent linear model of the functional.
The tangent linear model is evaluated w.r.t. the control on a vector m_dot, around the last supplied value of the control.
- Parameters:
m_dot ([pyadjoint.OverloadedType]) – The direction in which to compute the action of the tangent linear model.
- Returns:
- The action of the tangent linear model in the direction m_dot.
Should be an instance of the same type as the functional.
- Return type:
- hessian(m_dot: OverloadedType, hessian_input: OverloadedType | None = None, evaluate_tlm: bool = True, apply_riesz: bool = False)¶
Return the action of the Hessian of the functional.
The Hessian is evaluated w.r.t. the control on a vector m_dot.
Using the second-order adjoint method, the action of the Hessian of the functional with respect to the control, around the last supplied value of the control, is computed and returned.
- Parameters:
m_dot ([pyadjoint.OverloadedType]) – The direction in which to compute the action of the Hessian.
hessian_input (pyadjoint.OverloadedType) – The Hessian value for the functional result. Required if the functional is not scalar-valued, or if the functional is an intermediate result in the computation of an outer functional.
evaluate_tlm (bool) – If True, will evaluate the tangent linear model before evaluating the Hessian. If False, assumes that the tape is already populated with the tangent linear values and does not evaluate them.
apply_riesz (bool) – If True, apply the Riesz map of each control in order to return the (primal) Riesz representer of the Hessian action.
- Returns:
- The action of the Hessian in the direction m_dot.
If apply_riesz is False, should be an instance of the type dual to that of the control. If apply_riesz is True should have the same type as the control.
- Return type:
- recording_stages(nstages=None, **stage_kwargs)¶
- class fdvar.sc4dvar_reduced_functional.SC4DVarObservationStageSequence(control: OverloadedType, rf: SC4DVarReducedFunctional, stage_kwargs: dict, nstages: int)¶
Bases:
object
- class fdvar.sc4dvar_reduced_functional.SC4DVarObservationStage(control: OverloadedType, rf: SC4DVarReducedFunctional, index: int)¶
Bases:
objectRecord an observation for strong constraint 4DVar at the time of state.
- Parameters:
control – The state at the beginning of this stage.
rf – The SC4DVarReducedFunctional.
index – The index of this stage, numbered from 0.
- set_observation(state: OverloadedType, observation_error: Callable[[OverloadedType], OverloadedType], observation_covariance: CovarianceOperatorBase)¶
Record an observation at the time of state.
- Parameters:
state – The state at the current observation time.
observation_error – Given a state \(x\), returns the observations error \(y_{i} - \mathcal{H}_{i}(x)\) where \(y_{i}\) are the observations at the current observation time and \(\mathcal{H}_{i}\) is the observation operator for the current observation time.
observation_covariance – The inner product to calculate the observation error functional from the observation error \(y_{i} - \mathcal{H}_{i}(x)\). Can include the error covariance matrix.
fdvar.wc4dvar_reduced_functional module¶
- class fdvar.wc4dvar_reduced_functional.WC4DVarReducedFunctional(control: Control, background_covariance: CovarianceOperatorBase, observation_covariance: CovarianceOperatorBase, observation_error: Callable[[OverloadedType], OverloadedType] = None, background: OverloadedType | None = None, gauss_newton: bool = False)¶
Bases:
AbstractReducedFunctionalReducedFunctional for weak constraint 4DVar data assimilation.
- Parameters:
control – The
EnsembleFunctionfor the control x_{i} at the initial condition and at the end of each observation stage.background_covariance – The inner product to calculate the background error functional from the background error \(x_{0} - x_{b}\). Can include the error covariance matrix. Only used on ensemble rank 0.
observation_covariance – The inner product to calculate the observation error functional from the observation error \(y_{0} - \mathcal{H}_{0}(x)\). Can include the error covariance matrix. Must be provided if observation_error is provided. Only used on ensemble rank 0
observation_error – Given a state \(x\), returns the observations error \(y_{0} - \mathcal{H}_{0}(x)\) where \(y_{0}\) are the observations at the initial time and \(\mathcal{H}_{0}\) is the observation operator for the initial time. Only used on ensemble rank 0. Optional.
background – The background (prior) data for the initial condition \(x_{b}\). If not provided, the value of the first subfunction on the first ensemble member of the control
EnsembleFunctionwill be used.
See also
- property controls¶
Return the list of controls on which this functional depends.
- __call__(values: OverloadedType)¶
Computes the reduced functional with supplied control value.
- Parameters:
values – If you have multiple controls this should be a list of new values for each control in the order you listed the controls to the constructor. If you have a single control it can either be a list or a single object. Each new value should have the same type as the corresponding control.
- Returns:
The computed value. Typically of instance of
pyadjoint.AdjFloat.- Return type:
- derivative(adj_input: float = 1.0, apply_riesz: bool = False)¶
Returns the derivative of the functional w.r.t. the control. Using the adjoint method, the derivative of the functional with respect to the control, around the last supplied value of the control, is computed and returned.
- Parameters:
adj_input – The adjoint input.
options – Additional options for the derivative computation.
- Returns:
The derivative with respect to the control. Should be an instance of the same type as the control.
- Return type:
- tlm(m_dot: OverloadedType)¶
Return the action of the tangent linear model of the functional.
The tangent linear model is evaluated w.r.t. the control on a vector m_dot, around the last supplied value of the control.
- Parameters:
m_dot ([pyadjoint.OverloadedType]) – The direction in which to compute the action of the tangent linear model.
- Returns:
- The action of the tangent linear model in the direction m_dot.
Should be an instance of the same type as the functional.
- Return type:
- hessian(m_dot: OverloadedType, hessian_input: OverloadedType | None = None, evaluate_tlm: bool = True, apply_riesz: bool = False)¶
Returns the action of the Hessian of the functional w.r.t. the control on a vector m_dot.
Using the second-order adjoint method, the action of the Hessian of the functional with respect to the control, around the last supplied value of the control, is computed and returned.
- Parameters:
m_dot – The direction in which to compute the action of the Hessian.
- Returns:
The action of the Hessian in the direction m_dot. Should be an instance of the same type as the control.
- Return type:
- recording_stages(sequential=True, nstages=None, **stage_kwargs)¶
- class fdvar.wc4dvar_reduced_functional.WC4DVarObservationStageSequence(controls: Control, global_index: int, observation_index: int, stage_kwargs: dict = None)¶
Bases:
object
- class fdvar.wc4dvar_reduced_functional.WC4DVarObservationStage(control: Control, local_index: int | None = None, global_index: int | None = None, observation_index: int | None = None)¶
Bases:
objectA single stage for weak constraint 4DVar at the time of state.
Records the time propagator from the control at the beginning of the stage, and the model and observation errors at the end of the stage.
- Parameters:
control – The control x_{i-1} at the beginning of the stage
local_index – The index of this stage in the timeseries on the local ensemble member.
global_index – The index of this stage in the global timeseries.
observation_index – The index of the observation at the end of this stage in the global timeseries. May be different from global_index if an observation is taken at the initial time.
- set_observation(state: OverloadedType, observation_error: Callable[[OverloadedType], OverloadedType], observation_covariance: CovarianceOperatorBase, forward_model_covariance: CovarianceOperatorBase)¶
Record an observation at the time of state.
- Parameters:
state – The state at the current observation time.
observation_error – Given a state \(x\), returns the observations error \(y_{i} - \mathcal{H}_{i}(x)\) where \(y_{i}\) are the observations at the current observation time and \(\mathcal{H}_{i}\) is the observation operator for the current observation time.
observation_covariance – The inner product to calculate the observation error functional from the observation error \(y_{i} - \mathcal{H}_{i}(x)\). Can include the error covariance matrix.
forward_model_covariance – The inner product to calculate the model error functional from the model error \(x_{i} - \mathcal{M}_{i}(x_{i-1})\). Can include the error covariance matrix.