pyPDAF.PDAF.init¶
- pyPDAF.PDAF.init()¶
This function initialises the PDAF system.
It is called once at the beginning of the assimilation.
The function specifies the type of DA methods, parameters of the filters, the MPI communicators, and other parallel options. The filter options including filtertype, subtype, param_int, and param_real are introduced in PDAF filter options wiki page. Note that the size of param_int and param_real depends on the filter type and subtype. However, for most filters, they require at least the state vector size and ensemble size for param_int, and the forgetting factor for param_real.
The MPI communicators asked by this function depends on the parallelisation strategy. For the default parallelisation strategy, the user can use the parallelisation module provided under in example directory without modifications. The parallelisation can differ based on online and offline cases. Users can also refer to parallelisation documentation for explanations or modifications.
This function also asks for a user-supplied function
py__init_ens_pdaf()
. This function is designed to provides an initial ensemble to the internal PDAF ensemble array. The internal PDAF ensemble then can be distributed to initialise the model forecast usingpyPDAF.PDAF.get_state()
. This user-supplied function can be empty if the model has already read the ensemble from restart files.- Parameters:
filtertype (int) – type of filter
subtype (int) – sub-type of filter
stepnull (int) – initial time step of assimilation
param_int (ndarray[tuple[dim_pint], np.intc]) –
integer filter parameters in all available filters, the size of state vector, and number of ensemble members must be given; additional required parameters are filter specific See pyPDAF filter parameters
The array dimension dim_pint is dimension of param_int.
param_real (ndarray[tuple[dim_preal], np.float64]) –
real/float filter parameters the forgetting factor must be given; additional required parameters are filter specific See pyPDAF filter parameters
The array dimension dim_preal is dimension of param_real.
COMM_model (int) – model MPI communicator
COMM_filter (int) – filter MPI communicator
COMM_couple (int) – coupling MPI communicator
task_id (int) – index of parallel model task on current process; the index starts from 1
n_modeltasks (int) – number of parallel model tasks See [pyPDAF parallel strategy](https://yumengch.github.io/pyPDAF/parallel.html)
in_filterpe (bool) – True if the current PE is a filter PE else False
py__init_ens_pdaf (Callable[filtertype:int, dim_p:int, dim_ens:int, state_p : ndarray[tuple[dim_p], np.float64], uinv : ndarray[tuple[:, :], np.float64], ens_p : ndarray[tuple[dim_p, dim_ens], np.float64], flag:int]) –
initialise PDAF internal ensemble array ens_p; for SEEK, one also need to fill uinv; this function is called by processes with filterpe = .true. only
- Callback Parameters
- filtertypeint
filter type given in PDAF_init
- dim_pint
PE-local state dimension given by PDAF_init
- dim_ensint
number of ensemble members
- state_pndarray[tuple[dim_p], np.float64]
PE-local model state This array must be filled with the initial state of the model for SEEK, but it is not used for ensemble-based filters.
One can still make use of this array within this function.
- uinvndarray[tuple[:, :], np.float64]
This array is the inverse of matrix formed by right singular vectors of error covariance matrix of ensemble perturbations.
This array has to be filled in SEEK, but it is not used for ensemble-based filters. Nevertheless, one can still make use of this array within this function e.g., for generating an initial ensemble perturbation from a given covariance matrix.
Dimension of this array is determined by the filter type.
(dim_ens, dim_ens) for (L)ETKF, (L)NETF, (L)KNETF, and SEEK
(dim_ens - 1, dim_ens - 1) for (L)SEIK, (L)ESTKF, and 3DVar using ensemble
(1, 1) for (L)EnKF, particle filters and gen_obs
- 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 model state This array must be filled with the initial state of the model for SEEK, but it is not used for ensemble-based filters.
One can still make use of this array within this function.
- uinvndarray[tuple[:, :], np.float64]
This array is the inverse of matrix formed by right singular vectors of error covariance matrix of ensemble perturbations.
This array has to be filled in SEEK, but it is not used for ensemble-based filters. Nevertheless, one can still make use of this array within this function e.g., for generating an initial ensemble perturbation from a given covariance matrix.
Dimension of this array is determined by the filter type.
(dim_ens, dim_ens) for (L)ETKF, (L)NETF, (L)KNETF, and SEEK
(dim_ens - 1, dim_ens - 1) for (L)SEIK, (L)ESTKF, and 3DVar using ensemble
(1, 1) for (L)EnKF, particle filters and gen_obs
- ens_pndarray[tuple[dim_p, dim_ens], np.float64]
PE-local ensemble
- flagint
pdaf status flag
in_screen (int) – Verbosity level of PDAF screen output
- Returns:
param_int (ndarray[tuple[dim_pint], np.intc]) – integer filter parameters in all available filters, the size of state vector, and number of ensemble members must be given; additional required parameters are filter specific See pyPDAF filter parameters
The array dimension dim_pint is dimension of param_int
param_real (ndarray[tuple[dim_preal], np.float64]) – real/float filter parameters the forgetting factor must be given; additional required parameters are filter specific See pyPDAF filter parameters
The array dimension dim_preal is dimension of param_real
flag (int) – Status flag, 0: no error, error codes: