Changes between Initial Version and Version 1 of ImplementAnalysisenkf


Ignore:
Timestamp:
May 17, 2011, 5:27:25 PM (13 years ago)
Author:
lnerger
Comment:

Added page on EnKF based on SEIK

Legend:

Unmodified
Added
Removed
Modified
  • ImplementAnalysisenkf

    v1 v1  
     1= Implementation of the Analysis step for the EnKF (Ensemble Kalman Filter) =
     2
     3{{{
     4#!html
     5<div class="wiki-toc">
     6<h4>Implementation Guide</h4>
     7<ol><li><a href="ImplementationGuide">Main page</a></li>
     8<li><a href="AdaptParallelization">Adaptation of the parallelization</a></li>
     9<li><a href="InitPdaf">Initialization of PDAF</a></li>
     10<li><a href="ModifyModelforEnsembleIntegration">Modifications for ensemble integration</a></li>
     11<li><a href="ImplementationofAnalysisStep">Implementation of the analysis step</a></li>
     12<ol>
     13<li><a href="ImplementAnalysisseik">Implementation for SEIK</a></li>
     14<li><a href="ImplementAnalysislseik">Implementation for LSEIK</a></li>
     15<li><a href="ImplementAnalysisetkf">Implementation for ETKF</a></li>
     16<li><a href="ImplementAnalysisletkf">Implementation for LETKF</a></li>
     17<li><a href="ImplementAnalysisseek">Implementation for SEEK</a></li>
     18<li>Implementation for EnKF</li>
     19</ol>
     20<li><a href="AddingMemoryandTimingInformation">Memory and timing information</a></li>
     21</ol>
     22</div>
     23}}}
     24
     25[[PageOutline(2-3,Contents of this page)]]
     26
     27== Overview ==
     28
     29For the analysis step of the EnKF 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_enkf`. With regard to the parallelization, all these routines are executed by the filter processes (`filterpe=1`) only.
     30
     31For completeness we discuss here all user-supplied routines that are specified in the interface to PDAF_put_state_enkf. Thus, some of the user-supplied routines that are explained on the page explaining the modification of the model code for the ensemble integration are repeated here.
     32
     33The EnKF implemented in PDAF follows the original EnKF by Evensen (1994) including the correction for perturbed observations (Burgers et al. 1998). The EnKF implemented in PDAF is reviewed by Nerger et al (2005) and described in more detail by Nerger (2004). (See the [PublicationsandPresentations page on publications and presentations] for publications and presenations involving and about PDAF)
     34
     35There is no localized variant of the EnKF in PDAF. In our studies (Nerger et al. 2005, Neger et al. 2007), the EnKF showed deficiencies compared to the SEIK filter performed. Due to this, we focused more on the SEIK filter and the ETKF after these comparison studies. For real applications, we generally recommend using SEIK or ETKF, or their local variants LSEIK or LETKF.
     36
     37== `PDAF_put_state_enkf` ==
     38
     39The 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 the full interface specifically for the EnKF. Subsequently, the full set of user-supplied routines specified in the call to `PDAF_put_state_enkf` is explained.
     40
     41The interface when using the EnKF is the following:
     42{{{
     43  SUBROUTINE PDAF_put_state_enkf(U_collect_state, U_init_dim_obs, U_obs_op, &
     44                                 U_init_obs, U_prepoststep, U_add_obs_err, U_init_obscovar, status)
     45}}}
     46with the following arguments:
     47 * [#U_collect_statecollect_state.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`
     48 * [#U_init_dim_obsinit_dim_obs.F90 U_init_dim_obs]: The name of the user-supplied routine that provides the size of observation vector
     49 * [#U_obs_opobs_op.F90 U_obs_op]: The name of the user-supplied routine that acts as the observation operator on some state vector
     50 * [#U_init_obsinit_obs.F90 U_init_obs]: The name of the user-supplied routine that initializes the vector of observations
     51 * [#U_prepoststepprepoststep_seik.F90 U_prepoststep]: The name of the pre/poststep routine as in `PDAF_get_state`
     52 * [#U_add_obs_erradd_obs_err.F90 U_add_obs_err]: The name of the user-supplied routine that adds the observation error covariance matrix to the ensemble covariance matrix projected onto the observation space.
     53 * [#U_init_obscovarinit_obscovar.F90 U_init_obscovar]: The name of the user-supplied routine that initializes the observation error covariance matrix.
     54 * `status`: The integer status flag. It is zero, if `PDAF_put_state_enkf` is exited without errors.
     55
     56Note:
     57 * The order of the routine names does not show the order in which these routines are executed. See the [#Executionorderofuser-suppliedroutines section on the order of the execution] at the bottom of this page.
     58
     59
     60== User-supplied routines ==
     61
     62Here all user-supplied routines are described that are required in the call to `PDAF_put_state_enkf`. For some of the generic routines, we link to the page on [ModifyModelforEnsembleIntegration modifying the model code for the ensemble integration].
     63
     64To 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 `_dummy_D.F90`. In the section titles below we provide the name of the template file in parentheses.
     65
     66In 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.
     67
     68
     69=== `U_collect_state` (collect_state.F90) ===
     70
     71This routine is independent of the filter algorithm used.
     72See the page [ModifyModelforEnsembleIntegration#U_collect_statecollect_state.F90 modifying the model code for the ensemble integration] for the description of this routine.
     73
     74
     75=== `U_init_dim_obs` (init_dim_obs.F90) ===
     76
     77This routine is used by all global filter algorithms (SEEK, SEIK, EnKF, ETKF).
     78
     79The interface for this routine is:
     80{{{
     81SUBROUTINE init_dim_obs(step, dim_obs_p)
     82
     83  INTEGER, INTENT(in)  :: step       ! Current time step
     84  INTEGER, INTENT(out) :: dim_obs_p  ! Dimension of observation vector
     85}}}
     86
     87The 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.
     88
     89Some hints:
     90 * 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.
     91
     92
     93=== `U_obs_op` (obs_op.F90) ===
     94
     95This routine is used by all global filter algorithms (SEEK, SEIK, EnKF, ETKF).
     96
     97The interface for this routine is:
     98{{{
     99SUBROUTINE obs_op(step, dim_p, dim_obs_p, state_p, m_state_p)
     100
     101  INTEGER, INTENT(in) :: step               ! Currrent time step
     102  INTEGER, INTENT(in) :: dim_p              ! PE-local dimension of state
     103  INTEGER, INTENT(in) :: dim_obs_p          ! Dimension of observed state
     104  REAL, INTENT(in)    :: state_p(dim_p)     ! PE-local model state
     105  REAL, INTENT(out) :: m_state_p(dim_obs_p) ! PE-local observed state
     106}}}
     107
     108The 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`.
     109
     110For 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.
     111
     112Hint:
     113 * 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.
     114
     115
     116=== `U_init_obs` (init_obs.F90) ===
     117
     118This routine is used by all global filter algorithms (SEEK, SEIK, EnKF, ETKF).
     119
     120The interface for this routine is:
     121{{{
     122SUBROUTINE init_obs(step, dim_obs_p, observation_p)
     123
     124  INTEGER, INTENT(in) :: step             ! Current time step
     125  INTEGER, INTENT(in) :: dim_obs_p        ! PE-local dimension of obs. vector
     126  REAL, INTENT(out)   :: observation_p(dim_obs_p) ! PE-local observation vector
     127}}}
     128
     129The routine is called during the analysis step.
     130It has to provide the vector of observations in `observation_p` for the current time step.
     131
     132For a model using domain decomposition, the vector of observations that exist on the model sub-domain for the calling process has to be initialized.
     133
     134
     135=== `U_prepoststep` (prepoststep_enkf.F90) ===
     136
     137The general aspects of this routines have already been described on the [ModifyModelforEnsembleIntegration#U_prepoststepprepoststep_seik.F90 page on modifying the model code for the ensemble integration] for the SEIK filter. For completeness, the description is repeated specifically for the EnKF:
     138
     139The interface of the routine is identical for all filters, but sizes can vary. Also, the particular operations that are performed in the routine can be specific for each filter algorithm.
     140
     141The interface for this routine is for the EnKF
     142{{{
     143SUBROUTINE prepoststep(step, dim_p, dim_ens, dim_ens_p, dim_obs_p, &
     144                       state_p, Uinv, ens_p, flag)
     145
     146  INTEGER, INTENT(in) :: step        ! Current time step
     147                         ! (When the routine is called before the analysis -step is provided.)
     148  INTEGER, INTENT(in) :: dim_p       ! PE-local state dimension
     149  INTEGER, INTENT(in) :: dim_ens     ! Size of state ensemble
     150  INTEGER, INTENT(in) :: dim_ens_p   ! PE-local size of ensemble
     151  INTEGER, INTENT(in) :: dim_obs_p   ! PE-local dimension of observation vector
     152  REAL, INTENT(inout) :: state_p(dim_p) ! PE-local forecast/analysis state
     153                                     ! The array 'state_p' is not generally not initialized in the case of SEIK/EnKF/ETKF.
     154                                     ! It can be used freely in this routine.
     155  REAL, INTENT(inout) :: Uinv(1, 1)  ! Not used not EnKF
     156  REAL, INTENT(inout) :: ens_p(dim_p, dim_ens)      ! PE-local state ensemble
     157  INTEGER, INTENT(in) :: flag        ! PDAF status flag
     158}}}
     159
     160The 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`).
     161
     162The routine provides for the user the full access to the ensemble of model states. 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.
     163
     164Hint:
     165 * If a user considers to perform adjustments to the estimates (e.g. for balances), this routine is the right place for it.
     166 * 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`.
     167 * The array `Uinv` is not used in the EnKF. Internally to PDAF, it is allocated to be of size (1,1).
     168 * Apart from the size of the array `Uinv`, the interface is identical for all ensemble filters (SEIK/ETKF/EnKF/LSEIK/LETKF). In general it should be possible to use an identical pre/poststep routine for all these filters.
     169
     170
     171=== `U_add_obs_err` (add_obs_err.F90) ===
     172
     173This routine is only used for the EnKF.
     174
     175The interface for this routine is:
     176{{{
     177SUBROUTINE add_obs_err(step, dim_obs, C)
     178
     179  INTEGER, INTENT(in) :: step                ! Current time step
     180  INTEGER, INTENT(in) :: dim_obs             ! Dimension of obs. vector
     181  REAL, INTENT(inout) :: C(dim_obs, dim_obs) ! Matrix to that the observation
     182                                             !    error covariance matrix is added
     183}}}
     184
     185The routine is called during the analysis step. During the analysis step of the EnKF, the projection of the ensemble covariance onto the observation space is computed. This matrix is provided to the routine as `C_p`. The routine has to add the observation error covariance matrix to `C_p`.
     186
     187The operation is for the global observation space. Thus, it is independent of whether the filter is executed with or without parallelization.
     188
     189Hints:
     190 * The routine does not require that the observation error covariance matrix is added as a full matrix. If the matrix is diagonal, only the diagonal elements have to be added.
     191
     192
     193
     194=== `U_init_obscovar` (init_obscovar.F90) ===
     195
     196This routine is only used for the EnKF.
     197
     198The interface for this routine is:
     199{{{
     200SUBROUTINE init_obscovar(step, dim_obs, dim_obs_p, covar, m_state_p, &
     201     isdiag)
     202
     203  INTEGER, INTENT(in) :: step                ! Current time step
     204  INTEGER, INTENT(in) :: dim_obs             ! Dimension of observation vector
     205  INTEGER, INTENT(in) :: dim_obs_p           ! PE-local dimension of observation vector
     206  REAL, INTENT(out) :: covar(dim_obs, dim_obs) ! Observation error covariance matrix
     207  REAL, INTENT(in)  :: m_state_p(dim_obs_p)  ! PE-local observation vector
     208  LOGICAL, INTENT(out) :: isdiag             ! Whether the observation error covar. matrix is diagonal
     209}}}
     210
     211The routine is called during the analysis step and is required for the generation of an ensemble of observations. It has to initialize the global observation error covariance matrix `covar`. In addition, the flag `isdiag` has to be initialized to provide the information, whether the observation error covariance matrix is diagonal.
     212
     213The operation is for the global observation space. Thus, it is independent of whether the filter is executed with or without parallelization.
     214
     215Hints:
     216 * The local observation vector `m_state_p` is provided to the routine for the case that the observation errors are relative to the value of the observation.
     217
     218
     219== Execution order of user-supplied routines ==
     220
     221The user-supplied routines are executed in the order listed below. 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`.
     222
     223Before the analysis step is called the following routine is executed:
     224 1. [#U_collect_statecollect_state.F90 U_collect_state]
     225
     226The 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:
     227 1. [#U_prepoststepprepoststep_seik.F90 U_prepoststep] (Call to act on the forecast ensemble, called with negative value of the time step)
     228 1. [#U_init_dim_obsinit_dim_obs.F90 U_init_dim_obs]
     229 1. [#U_obs_opobs_op.F90 U_obs_op] (`dim_ens` calls: one call for each ensemble member)
     230 1. [#U_add_obs_erradd_obs_err.F90 U_add_obs_err]
     231 1. [#U_init_obsinit_obs.F90 U_init_obs]
     232 1. [#U_init_obscovarinit_obscovar.F90 U_init_obscovar]
     233 1. [#U_obs_opobs_op.F90 U_obs_op] (`dim_ens` calls: one call for each ensemble member, repeated to reduce storage)
     234 1. [#U_prepoststepprepoststep_seik.F90 U_prepoststep] (Call to act on the analysis ensemble, called with (positive) value of the time step)
     235