pyPDAF.PDAF.put_state_seik

pyPDAF.PDAF.put_state_seik()

It is recommended to use pyPDAF.PDAF.omi_put_state_global() or pyPDAF.PDAF.omi_put_state_global_nondiagR().

PDAF-OMI modules require fewer user-supplied functions and improved efficiency.

This function will use singular evolutive interpolated Kalman filter [1] for a single DA step.

Compared to pyPDAF.PDAF.assimilate_seik(), this function has no get_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. The pyPDAF.PDAF.get_state() function follows this function call to ensure the sequential DA.

The function should be called at each model step.

The function executes the user-supplied functions 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 (for ensemble mean)

  5. py__init_obs_pdaf

  6. py__obs_op_pdaf (for each ensemble member)

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

  8. py__prodRinvA_pdaf

  9. core DA algorithm

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

References

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 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__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__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__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

Returns:

flag – Status flag

Return type:

int