pyPDAF.PDAF.assimilate_en3dvar_estkf¶
- pyPDAF.PDAF.assimilate_en3dvar_estkf()¶
It is recommended to use
pyPDAF.PDAF.omi_assimilate_en3dvar_estkf()
orpyPDAF.PDAF.omi_assimilate_en3dvar_estkf_nondiagR()
.PDAF-OMI modules require fewer user-supplied functions and improved efficiency.
3DEnVar for a single DA step. The background error covariance matrix is estimated by an ensemble. The 3DEnVar only calculates the analysis of the ensemble mean. An ESTKF is used along with 3DEnVar to generate ensemble perturbations. This function should be called at each model time step.
The function is a combination of
pyPDAF.PDAF.put_state_en3dvar_estkf()
andpyPDAF.PDAF.get_state()
.- The user-supplied functions are executed in the following sequence:
py__collect_state_pdaf
py__prepoststep_state_pdaf
py__init_dim_obs_pdaf
py__obs_op_pdaf
py__init_obs_pdaf
- the iterative optimisation:
py__cvt_ens_pdaf
py__obs_op_lin_pdaf
py__prodRinvA_pdaf
py__obs_op_adj_pdaf
py__cvt_adj_ens_pdaf
core 3DEnVar algorithm
py__cvt_ens_pdaf
- ESTKF:
py__init_dim_obs_pdaf
py__obs_op_pdaf (for ensemble mean)
py__init_obs_pdaf
py__obs_op_pdaf (for each ensemble member)
py__init_obsvar_pdaf (only relevant for adaptive forgetting factor schemes)
py__prodRinvA_pdaf
core ESTKF algorithm
py__prepoststep_state_pdaf
py__distribute_state_pdaf
py__next_observation_pdaf
Deprecated since version 1.0.0: This function is replaced by
pyPDAF.PDAF.omi_assimilate_en3dvar_estkf()
andpyPDAF.PDAF.omi_assimilate_en3dvar_estkf_nondiagR()
- Parameters:
py__collect_state_pdaf (Callable[dim_p:int, state_p : ndarray[tuple[dim_p], np.float64]]) –
Routine to collect a state vector
- Callback Parameters
- dim_pint
pe-local state dimension
- state_pndarray[tuple[dim_p], np.float64]
local state vector
- Callback Returns
- state_pndarray[tuple[dim_p], np.float64]
local state vector
py__distribute_state_pdaf (Callable[dim_p:int, state_p : ndarray[tuple[dim_p], np.float64]]) –
Routine to distribute a state vector
- Callback Parameters
- dim_pint
pe-local state dimension
- state_pndarray[tuple[dim_p], np.float64]
local state vector
- Callback Returns
- state_pndarray[tuple[dim_p], np.float64]
local state vector
py__init_dim_obs_pdaf (Callable[step:int, dim_obs_p:int]) –
Initialize dimension of observation vector
- Callback Parameters
- stepint
current time step
- dim_obs_pint
dimension of observation vector
- Callback Returns
- dim_obs_pint
dimension of observation vector
py__obs_op_pdaf (Callable[step:int, dim_p:int, dim_obs_p:int, state_p : ndarray[tuple[dim_p], np.float64], m_state_p : ndarray[tuple[dim_obs_p], np.float64]]) –
Observation operator
- Callback Parameters
- stepint
Current time step
- dim_pint
Size of state vector (local part in case of parallel decomposed state)
- dim_obs_pint
Size of observation vector
- state_pndarray[tuple[dim_p], np.float64]
Model state vector
- m_state_pndarray[tuple[dim_obs_p], np.float64]
Observed state vector (i.e. the result after applying the observation operator to state_p)
- Callback Returns
- m_state_pndarray[tuple[dim_obs_p], np.float64]
Observed state vector (i.e. the result after applying the observation operator to state_p)
py__init_obs_pdaf (Callable[step:int, dim_obs_p:int, observation_p : ndarray[tuple[dim_obs_p], np.float64]]) –
Initialize observation vector
- Callback Parameters
- stepint
Current time step
- dim_obs_pint
Size of the observation vector
- observation_pndarray[tuple[dim_obs_p], np.float64]
Vector of observations
- Callback Returns
- observation_pndarray[tuple[dim_obs_p], np.float64]
Vector of observations
py__prodRinvA_pdaf (Callable[step:int, dim_obs_p:int, rank:int, obs_p : ndarray[tuple[dim_obs_p], np.float64], A_p : ndarray[tuple[dim_obs_p, rank], np.float64], C_p : ndarray[tuple[dim_obs_p, rank], np.float64]]) –
Provide product R^-1 A
- Callback Parameters
- stepint
Current time step
- dim_obs_pint
Number of observations at current time step (i.e. the size of the observation vector)
- rankint
Number of the columns in the matrix processes here. This is usually the ensemble size minus one (or the rank of the initial covariance matrix)
- obs_pndarray[tuple[dim_obs_p], np.float64]
Vector of observations
- A_pndarray[tuple[dim_obs_p, rank], np.float64]
Input matrix provided by PDAF
- C_pndarray[tuple[dim_obs_p, rank], np.float64]
Output matrix
- Callback Returns
- C_pndarray[tuple[dim_obs_p, rank], np.float64]
Output matrix
py__cvt_ens_pdaf (Callable[iter:int, dim_p:int, dim_ens:int, dim_cvec_ens:int, ens_p : ndarray[tuple[dim_p, dim_ens], np.float64], v_p : ndarray[tuple[dim_cvec_ens], np.float64], Vv_p : ndarray[tuple[dim_p], np.float64]]) –
Apply control vector transform matrix (ensemble)
- Callback Parameters
- iterint
Iteration of optimization
- dim_pint
PE-local dimension of state
- dim_ensint
Ensemble size
- dim_cvec_ensint
Dimension of control vector
- ens_pndarray[tuple[dim_p, dim_ens], np.float64]
PE-local ensemble
- v_pndarray[tuple[dim_cvec_ens], np.float64]
PE-local control vector
- Vv_pndarray[tuple[dim_p], np.float64]
PE-local state increment
- Callback Returns
- Vv_pndarray[tuple[dim_p], np.float64]
PE-local state increment
py__cvt_adj_ens_pdaf (Callable[iter:int, dim_p:int, dim_ens:int, dim_cv_ens_p:int, ens_p : ndarray[tuple[dim_p, dim_ens], np.float64], Vcv_p : ndarray[tuple[dim_p], np.float64], cv_p : ndarray[tuple[dim_cv_ens_p], np.float64]]) –
Apply adjoint control vector transform matrix (ensemble var)
- Callback Parameters
- iterint
Iteration of optimization
- dim_pint
PE-local observation dimension
- dim_ensint
Ensemble size
- dim_cv_ens_pint
PE-local dimension of control vector
- ens_pndarray[tuple[dim_p, dim_ens], np.float64]
PE-local ensemble
- Vcv_pndarray[tuple[dim_p], np.float64]
PE-local input vector
- cv_pndarray[tuple[dim_cv_ens_p], np.float64]
PE-local result vector
- Callback Returns
- cv_pndarray[tuple[dim_cv_ens_p], np.float64]
PE-local result vector
py__obs_op_lin_pdaf (Callable[step:int, dim_p:int, dim_obs_p:int, state_p : ndarray[tuple[dim_p], np.float64], m_state_p : ndarray[tuple[dim_obs_p], np.float64]]) –
Linearized observation operator
- Callback Parameters
- stepint
Current time step
- dim_pint
PE-local dimension of state
- dim_obs_pint
Dimension of observed state
- state_pndarray[tuple[dim_p], np.float64]
PE-local model state
- m_state_pndarray[tuple[dim_obs_p], np.float64]
PE-local observed state
- Callback Returns
- m_state_pndarray[tuple[dim_obs_p], np.float64]
PE-local observed state
py__obs_op_adj_pdaf (Callable[step:int, dim_p:int, dim_obs_p:int, state_p : ndarray[tuple[dim_p], np.float64], m_state_p : ndarray[tuple[dim_obs_p], np.float64]]) –
Adjoint observation operator
- Callback Parameters
- stepint
Current time step
- dim_pint
PE-local dimension of state
- dim_obs_pint
Dimension of observed state
- state_pndarray[tuple[dim_p], np.float64]
PE-local model state
- m_state_pndarray[tuple[dim_obs_p], np.float64]
PE-local observed state
- Callback Returns
- state_pndarray[tuple[dim_p], np.float64]
PE-local model state
py__init_obsvar_pdaf (Callable[step:int, dim_obs_p:int, obs_p : ndarray[tuple[dim_obs_p], np.float64], meanvar:float]) –
Initialize mean observation error variance
- Callback Parameters
- stepint
Current time step
- dim_obs_pint
Size of observation vector
- obs_pndarray[tuple[dim_obs_p], np.float64]
Vector of observations
- meanvarfloat
Mean observation error variance
- Callback Returns
- meanvarfloat
Mean observation error variance
py__prepoststep_pdaf (Callable[step:int, dim_p:int, dim_ens:int, dim_ens_p:int, dim_obs_p:int, state_p : ndarray[tuple[dim_p], np.float64], uinv : ndarray[tuple[dim_ens-1, dim_ens-1], np.float64], ens_p : ndarray[tuple[dim_p, dim_ens], np.float64], flag:int]) –
User supplied pre/poststep routine
- Callback Parameters
- stepint
current time step (negative for call after forecast)
- dim_pint
pe-local state dimension
- dim_ensint
size of state ensemble
- dim_ens_pint
pe-local size of ensemble
- dim_obs_pint
pe-local dimension of observation vector
- state_pndarray[tuple[dim_p], np.float64]
pe-local forecast/analysis state (the array ‘state_p’ is not generally not initialized in the case of seik. it can be used freely here.)
- uinvndarray[tuple[dim_ens-1, dim_ens-1], np.float64]
inverse of matrix u
- ens_pndarray[tuple[dim_p, dim_ens], np.float64]
pe-local state ensemble
- flagint
pdaf status flag
- Callback Returns
- state_pndarray[tuple[dim_p], np.float64]
pe-local forecast/analysis state (the array ‘state_p’ is not generally not initialized in the case of seik. it can be used freely here.)
- uinvndarray[tuple[dim_ens-1, dim_ens-1], np.float64]
inverse of matrix u
- ens_pndarray[tuple[dim_p, dim_ens], np.float64]
pe-local state ensemble
py__next_observation_pdaf (Callable[stepnow:int, nsteps:int, doexit:int, time:float]) –
Routine to provide time step, time and dimension of next observation
- Callback Parameters
- stepnowint
number of the current time step
- nstepsint
number of time steps until next obs
- doexitint
whether to exit forecasting (1 for exit)
- timefloat
current model (physical) time
- Callback Returns
- nstepsint
number of time steps until next obs
- doexitint
whether to exit forecasting (1 for exit)
- timefloat
current model (physical) time
- Returns:
outflag – Status flag
- Return type:
int