pyPDAF.PDAF.assimilate_en3dvar_estkf

pyPDAF.PDAF.assimilate_en3dvar_estkf()

It is recommended to use pyPDAF.PDAF.omi_assimilate_en3dvar_estkf() or pyPDAF.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() and pyPDAF.PDAF.get_state().

User-supplied functions are executed in the following sequence:
  1. py__collect_state_pdaf

  2. py__prepoststep_state_pdaf

  3. py__init_dim_obs_pdaf

  4. py__obs_op_pdaf

  5. py__init_obs_pdaf

  6. the iterative optimisation:
    1. py__cvt_ens_pdaf

    2. py__obs_op_lin_pdaf

    3. py__prodRinvA_pdaf

    4. py__obs_op_adj_pdaf

    5. py__cvt_adj_ens_pdaf

    6. core 3DEnVar algorithm

  7. py__cvt_ens_pdaf

  8. ESTKF:
    1. py__init_dim_obs_pdaf

    2. py__obs_op_pdaf (for ensemble mean)

    3. py__init_obs_pdaf

    4. py__obs_op_pdaf (for each ensemble member)

    5. py__init_obsvar_pdaf (only relevant for adaptive forgetting factor schemes)

    6. py__prodRinvA_pdaf

    7. core ESTKF algorithm

  9. py__prepoststep_state_pdaf

  10. py__distribute_state_pdaf

  11. py__next_observation_pdaf

Deprecated since version 1.0.0: This function is replaced by pyPDAF.PDAF.omi_assimilate_en3dvar_estkf() and pyPDAF.PDAF.omi_assimilate_en3dvar_estkf_nondiagR()

Parameters:
  • py__collect_state_pdaf (Callable[dim_p:int, state_p : ndarray[tuple[dim_p], np.float64]]) –

    Collect state vector from model/any arrays to pdaf arrays

    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]]) –

    distribute a state vector from pdaf to the model/any arrays

    Callback Parameters
    • dim_pint
      • PE-local state dimension

    • state_pndarray[tuple[dim_p], np.float64]
      • PE-local state vector

    Callback Returns
    • state_pndarray[tuple[dim_p], np.float64]
      • PE-local state vector

  • py__init_dim_obs_pdaf (Callable[step:int, dim_obs_p:int]) –

    The primary purpose of this function is to obtain the dimension of the observation vector. In OMI, in this function, one also sets the properties of obs_f, read the observation vector from files, setting the observation error variance when diagonal observation error covariance matrix is used. The pyPDAF.PDAF.omi_gather_obs function is also called here.

    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 PE-local 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_l: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]) –

    Preprocesse the ensemble before analysis and postprocess the ensemble before distributing to the model for next forecast

    Callback Parameters
    • stepint
      • current time step (negative for call before analysis/preprocessing)

    • dim_pint
      • PE-local state vector dimension

    • dim_ensint
      • number of ensemble members

    • dim_ens_lint
      • number of ensemble members run serially on each model task

    • 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 generally not initialised in the case of ESTKF/ETKF/EnKF/SEIK, so it can be used freely here.)

    • uinvndarray[tuple[dim_ens-1, dim_ens-1], np.float64]
      • Inverse of the transformation matrix in ETKF and ESKTF; inverse of matrix formed by right singular vectors of error covariance matrix of ensemble perturbations in SEIK/SEEK. not used in EnKF.

    • ens_pndarray[tuple[dim_p, dim_ens], np.float64]
      • PE-local ensemble

    • flagint
      • pdaf status flag

    Callback Returns
    • state_pndarray[tuple[dim_p], np.float64]
      • pe-local forecast/analysis state (the array ‘state_p’ is generally not initialised in the case of ESTKF/ETKF/EnKF/SEIK, so it can be used freely here.)

    • uinvndarray[tuple[dim_ens-1, dim_ens-1], np.float64]
      • Inverse of the transformation matrix in ETKF and ESKTF; inverse of matrix formed by right singular vectors of error covariance matrix of ensemble perturbations in SEIK/SEEK. not used in EnKF.

    • ens_pndarray[tuple[dim_p, dim_ens], np.float64]
      • PE-local ensemble

  • py__next_observation_pdaf (Callable[stepnow:int, nsteps:int, doexit:int, time:float]) –

    Routine to provide number of forecast time steps until next assimilations, model physical time and end of assimilation cycles

    Callback Parameters
    • stepnowint
      • the current time step given by PDAF

    • nstepsint
      • number of forecast time steps until next assimilation; this can also be interpreted as number of assimilation function calls to perform a new assimilation

    • doexitint
      • whether to exit forecasting (1 for exit)

    • timefloat
      • current model (physical) time

    Callback Returns
    • nstepsint
      • number of forecast time steps until next assimilation; this can also be interpreted as number of assimilation function calls to perform a new assimilation

    • doexitint
      • whether to exit forecasting (1 for exit)

    • timefloat
      • current model (physical) time

Returns:

outflag – Status flag

Return type:

int