atmos_flux_inversion.variational module¶
Functions implementing 3D-Var.
- Signatures follow the functions in
Note
Forces in-memory computations. BFGS method requires this, and there are odd shape-mismatch errors if I just change to dask arrays. Conjugate gradient solvers may work better for dask arrays if we drop the covariance matrix from the return values.
-
atmos_flux_inversion.variational.
incr_chol
(background, background_covariance, observations, observation_covariance, observation_operator, reduced_background_covariance=None, reduced_observation_operator=None)[source]¶ Feed everything to scipy’s minimizer.
Use the change from the background to try to avoid precision loss. Also use Cholesky factorization of the covariances to speed solution of matrix equations.
Assumes everything follows a multivariate normal distribution with the specified covariance matrices. Under this assumption analysis_covariance is exact, and analysis is the Maximum Likelihood Estimator and the Best Linear Unbiased Estimator for the underlying state in the frequentist framework, and specify the posterior distribution for the state in the Bayesian framework. If these are not satisfied, these still form the Generalized Least Squares estimates for the state and an estimated uncertainty.
- Parameters
background (array_like[N]) – The background state estimate.
background_covariance (array_like[N, N]) – Covariance of background state estimate across realizations/ensemble members. “Ensemble” is here interpreted in the sense used in statistical mechanics or frequentist statistics, and may not be derived from a sample as in meteorological ensemble Kalman filters
observations (array_like[M]) – The observations constraining the background estimate.
observation_covariance (array_like[M, M]) – Covariance of observations across realizations/ensemble members. “Ensemble” again has the statistical meaning.
observation_operator (array_like[M, N]) – The relationship between the state and the observations.
reduced_background_covariance (array_like[Nred, Nred], optional) – The covariance for a smaller state space, usually obtained by reducing resolution in space and time. Note that reduced_observation_operator must also be provided
reduced_observation_operator (array_like[M, Nred], optional) – The relationship between the reduced state space and the observations. Note that reduced_background_covariance must also be provided.
- Returns
analysis (array_like[N]) – Analysis state estimate
analysis_covariance (array_like[Nred, Nred] or array_like[N, N]) – Estimated uncertainty of analysis across realizations/ensemble members. Calculated using reduced_background_covariance and reduced_observation_operator if possible
- Raises
ConvergenceError – If iterative solver does not converge
Notes
minimizes
\[(dx)^T P_B^{-1} (dx) + (y - h(x_0) - H dx)^T R^{-1} (y - h(x_0) - H dx)\]which has gradient
\[P_B^{-1} (dx) - H^T R^{-1} (y - h(x) - H dx)\]where \(x = x_0 + dx\)
-
atmos_flux_inversion.variational.
incremental
(background, background_covariance, observations, observation_covariance, observation_operator, reduced_background_covariance=None, reduced_observation_operator=None)[source]¶ Feed everything to scipy’s minimizer.
Use the change from the background to try to avoid precision loss.
Assumes everything follows a multivariate normal distribution with the specified covariance matrices. Under this assumption analysis_covariance is exact, and analysis is the Maximum Likelihood Estimator and the Best Linear Unbiased Estimator for the underlying state in the frequentist framework, and specify the posterior distribution for the state in the Bayesian framework. If these are not satisfied, these still form the Generalized Least Squares estimates for the state and an estimated uncertainty.
- Parameters
background (array_like[N]) – The background state estimate.
background_covariance (array_like[N, N]) – Covariance of background state estimate across realizations/ensemble members. “Ensemble” is here interpreted in the sense used in statistical mechanics or frequentist statistics, and may not be derived from a sample as in meteorological ensemble Kalman filters
observations (array_like[M]) – The observations constraining the background estimate.
observation_covariance (array_like[M, M]) – Covariance of observations across realizations/ensemble members. “Ensemble” again has the statistical meaning.
observation_operator (array_like[M, N]) – The relationship between the state and the observations.
reduced_background_covariance (array_like[Nred, Nred], optional) – The covariance for a smaller state space, usually obtained by reducing resolution in space and time. Note that reduced_observation_operator must also be provided
reduced_observation_operator (array_like[M, Nred], optional) – The relationship between the reduced state space and the observations. Note that reduced_background_covariance must also be provided.
- Returns
analysis (array_like[N]) – Analysis state estimate
analysis_covariance (array_like[Nred, Nred] or array_like[N, N]) – Estimated uncertainty of analysis across realizations/ensemble members. Calculated using reduced_background_covariance and reduced_observation_operator if possible
- Raises
ConvergenceError – If iterative solver does not converge
Notes
minimizes
\[(dx)^T P_B^{-1} (dx) + (y - h(x_0) - H dx)^T R^{-1} (y - h(x_0) - H dx)\]which has gradient
\[P_B^{-1} (dx) - H^T R^{-1} (y - h(x) - H dx)\]where \(x = x_0 + dx\)
-
atmos_flux_inversion.variational.
simple
(background, background_covariance, observations, observation_covariance, observation_operator, reduced_background_covariance=None, reduced_observation_operator=None)[source]¶ Feed everything to scipy’s minimizer.
Assumes everything follows a multivariate normal distribution with the specified covariance matrices. Under this assumption analysis_covariance is exact, and analysis is the Maximum Likelihood Estimator and the Best Linear Unbiased Estimator for the underlying state in the frequentist framework, and specify the posterior distribution for the state in the Bayesian framework. If these are not satisfied, these still form the Generalized Least Squares estimates for the state and an estimated uncertainty.
- Parameters
background (array_like[N]) – The background state estimate.
background_covariance (array_like[N, N]) – Covariance of background state estimate across realizations/ensemble members. “Ensemble” is here interpreted in the sense used in statistical mechanics or frequentist statistics, and may not be derived from a sample as in meteorological ensemble Kalman filters
observations (array_like[M]) – The observations constraining the background estimate.
observation_covariance (array_like[M, M]) – Covariance of observations across realizations/ensemble members. “Ensemble” again has the statistical meaning.
observation_operator (array_like[M, N]) – The relationship between the state and the observations.
reduced_background_covariance (array_like[Nred, Nred], optional) – The covariance for a smaller state space, usually obtained by reducing resolution in space and time. Note that reduced_observation_operator must also be provided
reduced_observation_operator (array_like[M, Nred], optional) – The relationship between the reduced state space and the observations. Note that reduced_background_covariance must also be provided.
- Returns
analysis (array_like[N]) – Analysis state estimate
analysis_covariance (array_like[Nred, Nred] or array_like[N, N]) – Estimated uncertainty of analysis across realizations/ensemble members. Calculated using reduced_background_covariance and reduced_observation_operator if possible
- Raises
ConvergenceError – If iterative solver does not converge
Notes
minimizes
\[(x - x_0)^T P_B^{-1} (x - x_0) + (y - h(x))^T R^{-1} (y - h(x))\]which has gradient
\[P_B^{-1} (x - x_0) + H^T R^{-1} (y - h(x))\]