= Implementation of the Analysis step for the SEEK filter =
{{{
#!html
}}}
[[PageOutline(2-3,Contents of this page)]]
== Overview ==
For the analysis step of the SEEK filter different operations related to the observations are needed. These operations are requested by PDAF by calling user-supplied routines. Intentionally, the operations are split into separate routines in order to keep the operations rather elementary. This procedure should simplify the implementation. The names of the required routines are specified in the call to the routine `PDAF_put_state_seek`. With regard to the parallelization, all these routines are executed by the filter processes (`filterpe=1`) only.
For completeness we discuss here all user-supplied routines that are specified in the interface to PDAF_put_state_seek. Thus, some of the user-supplied routines that are explained on the page explaining the modification of the model code for the ensemble integration for the SEIK filter are repeated here, but specified for the SEEK filter.
The SEEK filter is very similar to the SEIK filter. In fact, the SEIK filter has been introduced as an interpolated (Pham et al., 1998) version of the SEEK filter. Due to the similarity of both filters, the interface to the user-supplied routines is almost identical. Several of the user-suppplied routines can be identical for SEEK and SEIK. Differences are marked in the text below. The implementation of the SEEK filter follows its original description by Pham et al. (1998) as reviewed by Nerger et al. (Tellus, 2005).
There is no localized variant of the SEEK filter in PDAF. In Nerger et al. (Tellus, 2005), the SEIK filter performed much better than the SEEK filter. Due to this, we focused more on the SEIK filter after this comparison study. For real applications, we recommend using SEIK or ETKF, or their local variants LSEIK or LETKF.
== `PDAF_put_state_seek` ==
The general espects of the filter specific routines `PDAF_put_state_*` have been described on the page [ModifyModelforEnsembleIntegration Modification of the model code for the ensemble integration]. Here, we list once more the full interface specifically for the SEEK filter. Subsequently, the full set of user-supplied routines specified in the call to `PDAF_put_state_seek` is explained.
The interface when using the SEEK filter is the following:
{{{
SUBROUTINE PDAF_put_state_seek(U_collect_state, U_init_dim_obs, U_obs_op, &
U_init_obs, U_prepoststep, U_prodRinvA, status)
}}}
with the following arguments:
* [#U_collect_statecollect_state_pdaf.F90 U_collect_state]: The name of the user-supplied routine that initializes a state vector from the array holding the ensemble of model states from the model fields. This is basically the inverse operation to `U_distribute_state` used in `PDAF_get_state`
* [#U_init_dim_obsinit_dim_obs_pdaf.F90 U_init_dim_obs]: The name of the user-supplied routine that provides the size of observation vector
* [#U_obs_opobs_op_pdaf.F90 U_obs_op]: The name of the user-supplied routine that acts as the observation operator on some state vector
* [#U_init_obsinit_obs_pdaf.F90 U_init_obs]: The name of the user-supplied routine that initializes the vector of observations
* [#U_prepoststepprepoststep_seek_pdaf.F90 U_prepoststep]: The name of the pre/poststep routine as in `PDAF_get_state`
* [#U_prodRinvAprodrinva_pdaf.F90 U_prodRinvA]: The name of the user-supplied routine that computes the product of the inverse of the observation error covariance matrix with some matrix provided to the routine by PDAF. This operation occurs during the analysis step of the SEEK filter.
* `status`: The integer status flag. It is zero, if `PDAF_put_state_seek` is exited without errors.
== User-supplied routines ==
Here, all user-supplied routines are described that are required in the call to `PDAF_put_state_seek`. For some of the generic routines, we link to the page on [ModifyModelforEnsembleIntegration modifying the model code for the ensemble integration].
To indicate user-supplied routines we use the prefix `U_`. In the template directory `templates/` these routines are provided in files with the routine's name without this prefix. In the example implementation in `testsuite/src/dummymodel_1D`, the routines exist without the prefix, but with the extension `_pdaf.F90`. In the section titles below we provide the name of the template file in parentheses.
In the subroutine interfaces some variables appear with the suffix `_p`. This suffix indicates that the variable is particular to a model sub-domain, if a domain decomposed model is used. Thus, the value(s) in the variable will be different for different model sub-domains.
=== `U_collect_state` (collect_state_pdaf.F90) ===
This routine is independent of the filter algorithm used.
See the page [ModifyModelforEnsembleIntegration#U_collect_statecollect_state_pdaf.F90 modifying the model code for the ensemble integration] for the description of this routine.
=== `U_init_dim_obs` (init_dim_obs_pdaf.F90) ===
This routine is used by all global filter algorithms (SEEK, SEIK, EnKF, ETKF).
The interface for this routine is:
{{{
SUBROUTINE init_dim_obs(step, dim_obs_p)
INTEGER, INTENT(in) :: step ! Current time step
INTEGER, INTENT(out) :: dim_obs_p ! Dimension of observation vector
}}}
The routine is called at the beginning of each analysis step. It has to initialize the size `dim_obs_p` of the observation vector according to the current time step. Without parallelization `dim_obs_p` will be the size for the full model domain. When a domain-decomposed model is used, `dim_obs_p` will be the size of the observation vector for the sub-domain of the calling process.
Some hints:
* It can be useful to not only determine the size of the observation vector at this point. One can also already gather information about the locations of the observations, which will be used later, e.g. to implement the observation operator. An array for the locations can be defined in a module like `mod_assimilation` of the example implementation.
=== `U_obs_op` (obs_op_pdaf.F90) ===
This routine is used by all global filter algorithms (SEEK, SEIK, EnKF, ETKF).
The interface for this routine is:
{{{
SUBROUTINE obs_op(step, dim_p, dim_obs_p, state_p, m_state_p)
INTEGER, INTENT(in) :: step ! Currrent time step
INTEGER, INTENT(in) :: dim_p ! PE-local dimension of state
INTEGER, INTENT(in) :: dim_obs_p ! Dimension of observed state
REAL, INTENT(in) :: state_p(dim_p) ! PE-local model state
REAL, INTENT(out) :: m_state_p(dim_obs_p) ! PE-local observed state
}}}
The routine is called during the analysis step. It has to perform the operation of the observation operator acting on a state vector that is provided as `state_p`. The observed state has to be returned in `m_state_p`.
For a model using domain decomposition, the operation is on the PE-local sub-domain of the model and has to provide the observed sub-state for the PE-local domain.
Hint:
* If the observation operator involves a global operation, e.g. some global integration, while using domain-decomposition one has to gather the information from the other model domains using MPI communication.
=== `U_init_obs` (init_obs_pdaf.F90) ===
This routine is used by all global filter algorithms (SEEK, SEIK, EnKF, ETKF).
The interface for this routine is:
{{{
SUBROUTINE init_obs(step, dim_obs_p, observation_p)
INTEGER, INTENT(in) :: step ! Current time step
INTEGER, INTENT(in) :: dim_obs_p ! PE-local dimension of obs. vector
REAL, INTENT(out) :: observation_p(dim_obs_p) ! PE-local observation vector
}}}
The routine is called during the analysis step.
It has to provide the vector of observations in `observation_p` for the current time step.
For a model using domain decomposition, the vector of observations that exist on the model sub-domain for the calling process has to be initialized.
=== `U_prepoststep` (prepoststep_seek_pdaf.F90) ===
The routine has already been described on the [ModifyModelforEnsembleIntegration#U_prepoststepprepoststep_seik_pdaf.F90 page on modifying the model code for the ensemble integration] for the SEIK filter. For the SEEK filter there are some differences, because of the fact that the covariance matrix is computed from the modes and eigenvalue matrix rather than from an ensemble of model states.
The interface of the routine is identical for all filters, but sizes can be different. In addition, the particular operations that are performed in the routine can be specific for each filter algorithm.
The interface for this routine is for the SEEK filter
{{{
SUBROUTINE prepoststep(step, dim_p, rank, rank_p, dim_obs_p, &
state_p, Uinv, ens_p, flag)
INTEGER, INTENT(in) :: step ! Current time step
! (When the routine is called before the analysis -step is provided.)
INTEGER, INTENT(in) :: dim_p ! PE-local state dimension
INTEGER, INTENT(in) :: dim_eof ! Rank of covariance matrix (Number of EOF modes)
INTEGER, INTENT(in) :: dim_eof_p ! PE-local rank of covariance matrix/EOF modes
INTEGER, INTENT(in) :: dim_obs_p ! PE-local dimension of observation vector
REAL, INTENT(inout) :: state_p(dim_p) ! PE-local forecast/analysis state
REAL, INTENT(inout) :: Uinv(dim_eof, dim_eof) ! Inverse of matrix U
REAL, INTENT(inout) :: eofV_p(dim_p, dim_ens) ! PE-local mode matrix V
INTEGER, INTENT(in) :: flag ! PDAF status flag
}}}
The routine `U_prepoststep` is called once at the beginning of the assimilation process. In addition, it is called during the assimilation cycles before the analysis step and after the ensemble transformation. The routine is called by all filter processes (that is `filterpe=1`).
The routine provides for the user the full access to the EOF modes of the SEEK filter, as well as the eigenvalue matrix `Uinv` and the state estimate `state_p`. Thus, user-controlled pre- and post-step operations can be performed. For example the forecast and the analysis states and ensemble covariance matrix can be analyzed, e.g. by computing the estimated variances. In addition, the estimates can be written to disk.
Hint:
* If a user considers to perform adjustments to the estimates (e.g. for balances), this routine is the right place for it.
* Only for the SEEK filter the state vector (`state_p`) is initialized. For all other filters, the array is allocated, but it can be used freely during the execution of `U_prepoststep`.
=== `U_prodRinvA` (prodrinva_pdaf.F90) ===
This routine is used by all filter algorithms that use the inverse of the observation error covariance matrix (SEEK, SEIK, and ETKF).
The interface for this routine is:
{{{
SUBROUTINE prodRinvA(step, dim_obs_p, rank, obs_p, A_p, C_p)
INTEGER, INTENT(in) :: step ! Current time step
INTEGER, INTENT(in) :: dim_obs_p ! PE-local dimension of obs. vector
INTEGER, INTENT(in) :: rank ! Rank of initial covariance matrix/Number of EOF modes
REAL, INTENT(in) :: obs_p(dim_obs_p) ! PE-local vector of observations
REAL, INTENT(in) :: A_p(dim_obs_p, rank) ! Input matrix from analysis routine
REAL, INTENT(out) :: C_p(dim_obs_p, rank) ! Output matrix
}}}
The routine is called during the analysis step. In the algorithms the product of the inverse of the observation error covariance matrix with some matrix has to be computed. For the SEEK filter this matrix holds the observed part of the EOF modes. The matrix is provided as `A_p`. The product has to be given as `C_p`.
For a model with domain decomposition, `A_p` contains the part of the matrix that resides on the model sub-domain of the calling process. The product has to be computed for this sub-domain, too.
Hints:
* The routine does not require that the product is implemented as a real matrix-matrix product. Rather, the product can be implemented in its most efficient form. For example, if the observation error covariance matrix is diagonal, only the multiplication of the diagonal with matrix `A_p` has to be implemented.
* The observation vector `obs_p` is provided through the interface for cases where the observation error variance is relative to the actual value of the observations.
== Execution order of user-supplied routines ==
For the SEEK filter, the user-supplied routines are essentially executed in the order they are listed in the interface to `PDAF_put_state_seek`. The order can be important as some routines can perform preparatory work for later routines. For example, `U_init_dim_obs` can prepare an index array that provides the information for executing the observation operator in `U_obs_op`.
Before the analysis step is called the following routine is executed:
1. [#U_collect_statecollect_state_pdaf.F90 U_collect_state]
The analysis step is executed when the ensemble integration of the forecast is completed. During the analysis step the following routines are executed in the given order:
1. [#U_prepoststepprepoststep_seek_pdaf.F90 U_prepoststep] (Call to act on the forecast ensemble, called with negative value of the time step)
1. [#U_init_dim_obsinit_dim_obs_pdaf.F90 U_init_dim_obs]
1. [#U_obs_opobs_op_pdaf.F90 U_obs_op] (A single call to operate on the ensemble mean state)
1. [#U_init_obsinit_obs_pdaf.F90 U_init_obs]
1. [#U_obs_opobs_op_pdaf.F90 U_obs_op] (`dim_eof` calls: one call for each ensemble member)
1. [#U_prodRinvAprodrinva_pdaf.F90 U_prodRinvA]
1. [#U_prepoststepprepoststep_seek_pdaf.F90 U_prepoststep] (Call to act on the analysis ensemble, called with (positive) value of the time step)