pyPDAF.PDAF.omi_put_state_hyb3dvar_lestkf_nondiagR¶
- pyPDAF.PDAF.omi_put_state_hyb3dvar_lestkf_nondiagR()¶
It is recommended to use
pyPDAF.PDAF.localomi_put_state_hyb3dvar_lestkf_nondiagR()
orpyPDAF.PDAF.localomi_put_state_hyb3dvar_lestkf()
.PDAFlocal-OMI modules require fewer user-supplied functions and improved efficiency.
Hybrid 3DEnVar for a single DA step using non-diagnoal observation error covariance matrix without post-processing, distributing analysis, and setting next observation step, where the background error covariance is hybridised by a static background error covariance, and a flow-dependent background error covariance estimated from ensemble.
Compared to
pyPDAF.PDAF.omi_assimilate_hyb3dvar_lestkf_nondiagR()
, this function has noget_state()
call. This means that the analysis is not post-processed, and distributed to the model forecast by user-supplied functions. The next DA step will not be assigned by user-supplied functions as well. This function is typically used when there are not enough CPUs to run the ensemble in parallel, and some ensemble members have to be run serially. ThepyPDAF.PDAF.get_state()
function follows this function call to ensure the sequential DA.The 3DVar generates an ensemble mean and the ensemble perturbation is generated by LESTKF in this implementation. This function should be called at each model time step.
- 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
- The iterative optimisation:
py__cvt_pdaf
py__cvt_ens_pdaf
py__obs_op_lin_pdaf
py__prodRinvA_pdaf
py__obs_op_adj_pdaf
py__cvt_adj_pdaf
py__cvt_adj_ens_pdaf
core DA algorithm
py__cvt_pdaf
py__cvt_ens_pdaf
- Perform LESTKF:
py__init_n_domains_p_pdaf
py__init_dim_obs_pdaf
py__obs_op_pdaf (for each ensemble member)
- loop over each local domain:
py__init_dim_l_pdaf
py__init_dim_obs_l_pdaf
py__g2l_state_pdaf
py__prodRinvA_l_pdaf
core DA algorithm
py__l2g_state_pdaf
Deprecated since version 1.0.0: This function is replaced by
pyPDAF.PDAF.localomi_put_state_hyb3dvar_lestkf()
andpyPDAF.PDAF.localomi_put_state_hyb3dvar_lestkf_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__init_dim_obs_pdaf (Callable[step:int, dim_obs_p:int]) –
Initialize dimension of full 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]]) –
Full 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__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 to control vector
- 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
- 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__cvt_pdaf (Callable[iter:int, dim_p:int, dim_cvec:int, cv_p : ndarray[tuple[dim_cvec], np.float64], Vv_p : ndarray[tuple[dim_p], np.float64]]) –
Apply control vector transform matrix to control vector
- Callback Parameters
- iterint
Iteration of optimization
- dim_pint
PE-local observation dimension
- dim_cvecint
Dimension of control vector
- cv_pndarray[tuple[dim_cvec], np.float64]
PE-local control vector
- Vv_pndarray[tuple[dim_p], np.float64]
PE-local result vector (state vector increment)
- Callback Returns
- Vv_pndarray[tuple[dim_p], np.float64]
PE-local result vector (state vector increment)
py__cvt_adj_pdaf (Callable[iter:int, dim_p:int, dim_cvec:int, Vcv_p : ndarray[tuple[dim_p], np.float64], cv_p : ndarray[tuple[dim_cvec], np.float64]]) –
Apply adjoint control vector transform matrix
- Callback Parameters
- iterint
Iteration of optimization
- dim_pint
PE-local observation dimension
- dim_cvecint
Dimension of control vector
- Vcv_pndarray[tuple[dim_p], np.float64]
PE-local result vector (state vector increment)
- cv_pndarray[tuple[dim_cvec], np.float64]
PE-local control vector
- Callback Returns
- cv_pndarray[tuple[dim_cvec], np.float64]
PE-local control 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__prodRinvA_l_pdaf (Callable[domain_p:int, step:int, dim_obs_l:int, rank:int, obs_l : ndarray[tuple[dim_obs_l], np.float64], A_l : ndarray[tuple[dim_obs_l, rank], np.float64], C_l : ndarray[tuple[dim_obs_l, rank], np.float64]]) –
Provide product R^-1 A for analysis domain
- Callback Parameters
- domain_pint
Index of current local analysis domain
- stepint
Current time step
- dim_obs_lint
Number of local observations at current time step (i.e. the size of the local 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_lndarray[tuple[dim_obs_l], np.float64]
Local vector of observations
- A_lndarray[tuple[dim_obs_l, rank], np.float64]
Input matrix provided by PDAF
- C_lndarray[tuple[dim_obs_l, rank], np.float64]
Output matrix
- Callback Returns
- C_lndarray[tuple[dim_obs_l, rank], np.float64]
Output matrix
py__init_n_domains_p_pdaf (Callable[step:int, n_domains_p:int]) –
Provide number of local analysis domains
- Callback Parameters
- stepint
current time step
- n_domains_pint
pe-local number of analysis domains
- Callback Returns
- n_domains_pint
pe-local number of analysis domains
py__init_dim_l_pdaf (Callable[step:int, domain_p:int, dim_l:int]) –
Init state dimension for local ana. domain
- Callback Parameters
- stepint
current time step
- domain_pint
current local analysis domain
- dim_lint
local state dimension
- Callback Returns
- dim_lint
local state dimension
py__init_dim_obs_l_pdaf (Callable[domain_p:int, step:int, dim_obs_f:int, dim_obs_l:int]) –
Initialize local dimimension of obs. vector
- Callback Parameters
- domain_pint
index of current local analysis domain
- stepint
current time step
- dim_obs_fint
full dimension of observation vector
- dim_obs_lint
local dimension of observation vector
- Callback Returns
- dim_obs_lint
local dimension of observation vector
py__g2l_state_pdaf (Callable[step:int, domain_p:int, dim_p:int, state_p : ndarray[tuple[dim_p], np.float64], dim_l:int, state_l : ndarray[tuple[dim_l], np.float64]]) –
Get state on local ana. domain from full state
- Callback Parameters
- stepint
current time step
- domain_pint
current local analysis domain
- dim_pint
pe-local full state dimension
- state_pndarray[tuple[dim_p], np.float64]
pe-local full state vector
- dim_lint
local state dimension
- state_lndarray[tuple[dim_l], np.float64]
state vector on local analysis domain
- Callback Returns
- state_lndarray[tuple[dim_l], np.float64]
state vector on local analysis domain
py__l2g_state_pdaf (Callable[step:int, domain_p:int, dim_l:int, state_l : ndarray[tuple[dim_l], np.float64], dim_p:int, state_p : ndarray[tuple[dim_p], np.float64]]) –
Init full state from local state
- Callback Parameters
- stepint
current time step
- domain_pint
current local analysis domain
- dim_lint
local state dimension
- state_lndarray[tuple[dim_l], np.float64]
state vector on local analysis domain
- dim_pint
pe-local full state dimension
- state_pndarray[tuple[dim_p], np.float64]
pe-local full state vector
- Callback Returns
- state_pndarray[tuple[dim_p], np.float64]
pe-local full state vector
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
outflag (int) – Status flag
- Returns:
outflag – Status flag
- Return type:
int