Changes between Initial Version and Version 1 of Implement3DVarAnalysisPDAF3_3DEnVar


Ignore:
Timestamp:
May 26, 2025, 4:52:08 PM (8 weeks ago)
Author:
lnerger
Comment:

--

Legend:

Unmodified
Added
Removed
Modified
  • Implement3DVarAnalysisPDAF3_3DEnVar

    v1 v1  
     1= Implementation of the Analysis Step of 3D Ensemble Var  =
     2
     3{{{
     4#!html
     5<div class="wiki-toc">
     6<h4>Implementation Guide - Analysis Step</h4>
     7<ol><li>Implementing the analysis step</li>
     8<ol>
     9<li><b>Ensemble filters</b></li>
     10<ol>
     11<li> <a href="ImplementFilterAnalysisOverviewPDAF3"> General overview for ensemble filters</a></li>
     12<li><a href="ImplementAnalysisPDAF3Universal">Universal interface </a></li>
     13<li><a href="ImplementAnalysisPDAF3UniversalLocal">Universal interface using g2l/l2g_state</a></li>
     14<li><a href="ImplementanalysisPDAF3Gloval">Interface specific for global filters</a></li>
     15</ol>
     16<li><b>3D-Var methods</b></li>
     17<ol>
     18<li> <a href="Implement3DVarAnalysisOverviewPDAF3"> General overview for 3D-Var methods</a></li>
     19<li><a href="Implement3DVarAnalysisPDAF3Universal">Universal interface for 3D-Var</a></li>
     20<li><a href="Implement3DVarAnalysisPDAF3_3DVar">Implementation for 3D-Var</a></li>
     21<li><a href="Implement3DVarAnalysisPDAF3_3DEnVar">Implementation for 3D Ensemble Var</a></li>
     22<li><a href="Implement3DVarAnalysisPDAF3_Hyb3DVar">Implementation for Hybrid 3D-Var</a></li>
     23</ol>
     24
     25<li><a href="nondiagonal_observation_error_covariance_matrices_PDAF3">Using nondiagonal R-matrices</a></li>
     26<li><a href="PDAF_OMI_Overview">PDAF-OMI Overview</a></li>
     27</ol>
     28</div>
     29}}}
     30
     31
     32[[PageOutline(2-3,Contents of this page)]]
     33
     34== Overview ==
     35
     36This page describes the recommended implementation of the analysis step for the 3D ensemble Var schemes using the PDAF3 interface.
     37
     38|| The interface for 3D ensemble Var are specialized version of the universal interface with a reduced number of arguments. If one implements both 3D ensemble Var and parameterlized 3D-Var we recommend to use the [wiki:Implement3DVarAnalysisPDAF3Universal universal assimilation routines for 3D-Var]. ||
     39
     40For the analysis step of 3D ensemble Var we need different operations related to the observations. These operations are requested by PDAF by call-back routines supplied by the user and provided in the OMI structure. The names of the routines that are provided by the user are specified in the call to the assimilation routines. With regard to the parallelization, all these routines (except `collect_state_pdaf`, `distribute_state_pdaf`, and `next_observation_pdaf`) are executed by the filter processes (`filterpe=.true.`) only.
     41
     42The different 3D-Var methods in PDAF were explained on the [wiki:Implement3DVarAnalysisOverview page providing the verview of the Analysis Step for 3D-Var Methods]. Depending the type of 3D-Var, the background covariance matrix '''B''' is represented either in a parameterized form, by an ensemble, or by a combination of both. The 3D-Var methods that use an ensemble need to transform the ensemble perturbations using an ensemble Kalman filter. PDAF uses for this the error-subspace transform filter ESTKF. There are two variants: The first uses the localized filter LESTKF, while the second uses the global filter ESTKF.
     43
     44For completeness we discuss here all user-supplied routines that are specified in the interface to `PDAF3_assimilate_en3dvar_X`. Thus, some of the user-supplied routines that are explained on the page describing the modification of the model code for the ensemble integration are repeated here.
     45
     46
     47== Analysis Routines ==
     48
     49The general aspects of the filter (or solver) specific routines `PDAF3_assimilate_*` have been described on the page [OnlineModifyModelforEnsembleIntegration_PDAF3 Modification of the model code for the ensemble integration]. Here, we list the full interface of the routine. Subsequently, the user-supplied routines specified in the call is explained.
     50
     51There are two variants that either compute the transformataion of the ensemble transformation using the local LESTKF method, or the global ESTKF.
     52
     53=== `PDAF3_assimilate_en3dvar` ===
     54
     55This routine can be used to run both variants of the 3D ensemble Var, i.e., using the local LESTKF or the global ESTKF for the transformation of the ensemble perturbations.
     56
     57This routine is used both in the ''fully-parallel'' and the ''flexible'' implementation variants of the data assimilation system. (See the page [ModifyModelforEnsembleIntegration Modification of the model code for the ensemble integration] for these variants)
     58
     59The interface is:
     60{{{
     61SUBROUTINE PDAF3_assimilate_en3dvar(collect_state_pdaf, distribute_state_pdaf, &
     62                                 init_dim_obs_pdafomi, obs_op_pdafomi, &
     63                                 cvt_ens_pdaf, cvt_adj_ens_pdaf, &
     64                                 obs_op_lin_pdafomi, obs_op_adj_pdafomi, &
     65                                 init_n_domains_p_pdaf, init_dim_l_pdaf, init_dim_obs_l_pdafomi, &
     66                                 prepoststep_pdaf, next_observation_pdaf, outflag)
     67}}}
     68where all arguments, except the last one, are the names of call-back routines:
     69 * [#collect_state_pdafcollect_state_pdaf.F90 collect_state_pdaf]: 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 `distribute_state` used in `PDAF_init_forecast` as well as here.
     70 * [#distribute_state_pdafdistribute_state_pdaf.F90 distribute_state_pdaf]:  The name of a user supplied routine that initializes the model fields from the array holding the ensemble of model state vectors.
     71 * [#init_dim_obs_pdafomicallback_obs_pdafomi.F90 init_dim_obs_pdafomi]: The name of the user-supplied routine that initializes the observation information and provides the size of observation vector
     72 * [#obs_op_pdafomicallback_obs_pdafomi.F90 obs_op_pdafomi]: The name of the user-supplied routine that acts as the observation operator on some state vector
     73 * [#cvt_ens_pdafcvt_ens_pdaf.F90 cvt_ens_pdaf]: The name of the user-supplied routine that applies the ensemble control-vector transformation (square-root of the B-matrix) on some control vector to obtain a state vector.
     74 * [#cvt_adj_ens_pdafcvt_adj_ens_pdaf.F90 cvt_adj_ens_pdaf]: The name of the user-supplied routine that applies the adjoint ensemble control-vector transformation (with square-root of the B-matrix) on some state vector to obtain the control vector.
     75 * [#obs_op_pdafomicallback_obs_pdafomi.F90 obs_op_lin_pdafomi]: The name of the user-supplied routine that acts as the linearized observation operator on some state vector
     76 * [#obs_op_pdafomicallback_obs_pdafomi.F90 obs_op_lin_pdafomi]: The name of the user-supplied routine that acts as the adjoint observation operator on some state vector
     77 * [#init_n_domains_pdafinit_n_domains_pdaf.F90 init_n_domains_pdaf]: The name of the routine that provides the number of local analysis domains
     78 * [#init_dim_l_pdafinit_dim_l_pdaf.F90 init_dim_l_pdaf]: The name of the routine that provides the state dimension for a local analysis domain
     79 * [#init_dim_obs_l_pdafomicallback_obs_pdafomi.F90 init_dim_obs_l_pdafomi]: The name of the routine that initializes the size of the observation vector for a local analysis domain
     80 * [#prepoststep_pdafprepoststep_ens_pdaf.F90 prepoststep_pdaf]: The name of the pre/poststep routine as in `PDAF_init_forecast`
     81 * [#next_observation_pdafnext_observation.F90 next_observation_pdaf]: The name of a user supplied routine that initializes the variables `nsteps`, `timenow`, and `doexit`. The same routine is also used in `PDAF_init_forecast`.
     82 * `status`: The integer status flag. It is zero, if the routine is exited without errors.
     83
     84=== `PDAF3_assim_offline_3dvar_all` ===
     85
     86For the offline mode of PDAF, the routine `PDAF3_assim_offline_3dvar_all` is used to perform the analysis step.
     87The interface of the routine is identical with that of `PDAF3_assimilate_3dvar_all`, except that the user-supplied routines `U_distribute_state`, `U_collect_state` and `U_next_observation` are missing.
     88
     89The interface when using one of the global filters is the following:
     90{{{
     91  SUBROUTINE PDAF3_assim_offline_en3dvar(&
     92                                 U_init_dim_obs_pdafomi, U_obs_op_pdafomi, &
     93                                 U_cvt_ens, U_cvt_adj_ens, &
     94                                 U_obs_op_lin_pdafomi, U_obs_op_adj_pdafomi, &
     95                                 U_init_n_domains_p, U_init_dim_l, U_init_dim_obs_l_pdafomi, &
     96                                 U_prepoststep, outflag)
     97}}}
     98
     99=== `PDAF3_put_state_en3dvar` ===
     100
     101This routine exists for backward-compatibility. In implementations that were done before the release of PDAF V3.0, a 'put_state' routine was used for the `flexible` parallelization variant and for the offline mode.
     102When the 'flexible' implementation variant is chosen for the assimilation system, the routine. The routine `PDAF3_put_state_3dvar_all` allows to port such implemnetations to the PDAF3 interface with minimal changes.
     103The interface of the routine is identical with that of `PDAF3_assimilate_3dvar_all`, except that the user-supplied routines `U_distribute_state` and `U_next_observation` are missing.
     104
     105The interface when using one of the global filters is the following:
     106{{{
     107  SUBROUTINE PDAF3_put_state_en3dvar(U_collect_state, &
     108                                 U_init_dim_obs_pdafomi, U_obs_op_pdafomi, &
     109                                 U_cvt_ens, U_cvt_adj_ens, &
     110                                 U_obs_op_lin_pdafomi, U_obs_op_adj_pdafomi, &
     111                                 U_init_n_domains_p, U_init_dim_l, U_init_dim_obs_l_pdafomi, &
     112                                 U_prepoststep, outflag)
     113}}}
     114
     115
     116== Analysis Routines specific for using global ESTKF ==
     117
     118=== `PDAF3_assimilate_en3dvar_estkf` ===
     119
     120This routine is particular for the ESTKF. One can use it if one exclusively uses the global filter. In the argument list of this routine, the call-back routine for localization are not present and hence the argument list is shorter than that of `PDAF3_assimilate_en3dvar`.
     121 
     122This routine is used both in the ''fully-parallel'' and the ''flexible'' implementation variants of the data assimilation system. (See the page [ModifyModelforEnsembleIntegration Modification of the model code for the ensemble integration] for these variants)
     123
     124The interface is:
     125{{{
     126SUBROUTINE PDAF3_assimilate_en3dvar_estkf(collect_state_pdaf, distribute_state_pdaf, &
     127                                 init_dim_obs_pdafomi, obs_op_pdafomi, &
     128                                 cvt_ens_pdaf, cvt_adj_ens_pdaf, &
     129                                 obs_op_lin_pdafomi, obs_op_adj_pdafomi, &
     130                                 prepoststep_pdaf, next_observation_pdaf, outflag)
     131}}}
     132where all arguments, except the last one, are the names of call-back routines. See the description of the arguments for `PDAF3_assimilate_en3dvar`.
     133
     134
     135
     136=== `PDAF3_assim_offline_hyb3dvar_estkf` ===
     137
     138This routine is particular for the ESTKF. One can use it if one exclusively uses the global filter. In the argument list of this routine, the call-back routine for localization are not present and hence the argument list is shorter than that of `PDAF3_assim_offline_en3dvar`.
     139
     140This routine is used to perform the analysis step for the offline mode.
     141The interface of the routine is identical with that of `PDAF3_assimilate_en3dvar_estkf`, except that the user-supplied routines `U_distribute_state`, `U_collect_state` and `U_next_observation` are missing.
     142
     143The interface is:
     144{{{
     145SUBROUTINE PDAF3_assimilate_en3dvar_estkf(&
     146                                 init_dim_obs_pdafomi, obs_op_pdafomi, &
     147                                 cvt_ens_pdaf, cvt_adj_ens_pdaf, &
     148                                 obs_op_lin_pdafomi, obs_op_adj_pdafomi, &
     149                                 prepoststep_pdaf, outflag)
     150}}}
     151where all arguments, except the last one, are the names of call-back routines. See the description of the arguments for `PDAF3_assimilate_3dvar_all`.
     152
     153
     154
     155=== `PDAF3_put_state_en3dvar_estkf` ===
     156
     157This routine exists for backward-compatibility. In implementations that were done before the release of PDAF V3.0, a 'put_state' routine was used for the `flexible` parallelization variant and for the offline mode.
     158When the 'flexible' implementation variant is chosen for the assimilation system, the routine. The routine `PDAF3_put_state_en3dvar_estkf` allows to port such implemnetations to the PDAF3 interface with minimal changes.
     159The interface of the routine is identical with that of `PDAF3_assimilate_en3dvar_estkf`, except that the user-supplied routines `U_distribute_state` and `U_next_observation` are missing.
     160
     161This routine is particular for the ESTKF. One can use it if one exclusively uses the global filter. In the argument list of this routine, the call-back routine for localization are not present and hence the argument list is shorter than that of `PDAF3_put_state_en3dvar`.
     162
     163The interface is:
     164{{{
     165SUBROUTINE PDAF3_assimilate_hyb3dvar_estkf(collect_state_pdaf, &
     166                                 init_dim_obs_pdafomi, obs_op_pdafomi, &
     167                                 cvt_ens_pdaf, cvt_adj_ens_pdaf,  &
     168                                 obs_op_lin_pdafomi, obs_op_adj_pdafomi, &
     169                                 prepoststep_pdaf, outflag)
     170}}}
     171where all arguments, except the last one, are the names of call-back routines. See the description of the arguments for `PDAF3_assimilate_en3dvar`.
     172
     173
     174
     175
     176== User-supplied routines ==
     177
     178Here all user-supplied routines are described that are required in the call to the assimilation routines for hybrid 3D-Var. For some of the generic routines, we link to the page on [wiki:OnlineModifyModelforEnsembleIntegration_PDAF3 modifying the model code for the ensemble integration].
     179
     180To indicate user-supplied routines we use the prefix `U_`. In the template directory `templates/` as well as in the tutorial implementations in `tutorial/` these routines exist without the prefix, but with the extension `_pdaf.F90`. The user-routines relating to OMI are collected in the file `callback_obs_pdafomi.F90`. In the section titles below we provide the name of the template file in parentheses.
     181
     182In 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.
     183
     184
     185=== `collect_state_pdaf` (collect_state_pdaf.F90) ===
     186
     187This routine is independent of the filter algorithm used.
     188
     189See the page on [ModifyModelforEnsembleIntegration#collect_state_pdafcollect_state_pdaf.F90 modifying the model code for the ensemble integration] for the description of this routine.
     190
     191=== `distribute_state_pdaf` (distribute_state_pdaf.F90) ===
     192
     193This routine is independent of the filter algorithm used.
     194
     195See the page on [ModifyModelforEnsembleIntegration#distribute_state_pdafdistribute_state_pdaf.F90 modifying the model code for the ensemble integration] for the description of this routine.
     196
     197
     198
     199=== `init_dim_obs_pdafomi` (callback_obs_pdafomi.F90) ===
     200
     201This is a call-back routine initializing the observation information. The routine just calls a routine from the observation module for each observation type.
     202
     203See the [wiki:OMI_Callback_obs_pdafomi documentation on callback_obs_pdafomi.F90] for more information.
     204
     205
     206
     207=== `obs_op_pdafomi` (callback_obs_pdafomi.F90) ===
     208
     209This is a call-back routine applying the observation operator to the state vector. The routine calls a routine from the observation module for each observation type.
     210
     211See the [wiki:OMI_Callback_obs_pdafomi documentation on callback_obs_pdafomi.F90] for more information.
     212
     213
     214
     215=== `cvt_ens_pdaf` (cvt_ens_pdaf.F90) ===
     216
     217The interface for this routine is:
     218{{{
     219SUBROUTINE cvt_ens_pdaf(iter, dim_p, dim_ens, dim_cv_ens_p, ens_p, cv_p, Vcv_p)
     220
     221  INTEGER, INTENT(in) :: iter               ! Iteration of optimization
     222  INTEGER, INTENT(in) :: dim_p              ! PE-local observation dimension
     223  INTEGER, INTENT(in) :: dim_ens            ! Ensemble size
     224  INTEGER, INTENT(in) :: dim_cv_ens_p       ! Dimension of control vector
     225  REAL, INTENT(in) :: ens_p(dim_p, dim_ens) ! PE-local ensemble
     226  REAL, INTENT(in) :: cv_p(dim_cv_ens_p)    ! PE-local control vector
     227  REAL, INTENT(inout) :: Vcv_p(dim_p)       ! PE-local state increment
     228}}}
     229
     230The routine is called during the analysis step during the iterative minimization of the cost function.
     231It has to apply the control vector transformation to the control vector and return the transformed result vector. Usually this transformation is the multiplication with the square-root of the background error covariance matrix '''B'''. For the 3D Ensemble Var, this square root is usually expressed through the ensemble.
     232
     233If the control vector is decomposed in case of parallelization it first needs to the gathered on each processor and afterwards the transformation is computed on the potentially domain-decomposed state vector.
     234
     235
     236=== `cvt_adj_ens_pdaf` (cvt_adj_pdaf.F90) ===
     237
     238The interface for this routine is:
     239{{{
     240SUBROUTINE cvt_adj_ens_pdaf(iter, dim_p, dim_ens, dim_cv_ens_p, ens_p, Vcv_p, cv_p)
     241
     242  INTEGER, INTENT(in) :: iter                ! Iteration of optimization
     243  INTEGER, INTENT(in) :: dim_p               ! PE-local observation dimension
     244  INTEGER, INTENT(in) :: dim_ens             ! Ensemble size
     245  INTEGER, INTENT(in) :: dim_cv_ens_p        ! PE-local dimension of control vector
     246  REAL, INTENT(in) :: ens_p(dim_p, dim_ens)  ! PE-local ensemble
     247  REAL, INTENT(in)    :: Vcv_p(dim_p)        ! PE-local input vector
     248  REAL, INTENT(inout) :: cv_p(dim_cv_ens_p)  ! PE-local result vector
     249}}}
     250
     251The routine is called during the analysis step during the iterative minimization of the cost function.
     252It has to apply the adjoint control vector transformation to a state vector and return the control vector. Usually this transformation is the multiplication with transpose of the square-root of the background error covariance matrix '''B'''. or the 3D Ensemble Var, this square root is usually expressed through the ensemble.
     253
     254If the state vector is decomposed in case of parallelization one needs to take care that the application of the trasformation is complete. This usually requries a comminucation with MPI_Allreduce to obtain a global sun.
     255
     256
     257
     258=== `obs_op_lin_pdafomi` (callback_obs_pdafomi.F90) ===
     259
     260This is a call-back routine for PDAF-OMI applying the linearized observation operator to the state vector. The routine calls a routine from the observation module for each observation type. If the full observation operator is lineaer the same operator can be used here.
     261
     262See the [wiki:OMI_Callback_obs_pdafomi documentation on callback_obs_pdafomi.F90] for more information.
     263
     264
     265=== `obs_op_adj_pdafomi` (callback_obs_pdafomi.F90) ===
     266
     267This is a call-back routine for PDAF-OMI applying the adjoint observation operator to some vector inthe observation space. The routine calls a routine from the observation module for each observation type.
     268
     269See the [wiki:OMI_Callback_obs_pdafomi documentation on callback_obs_pdafomi.F90] for more information.
     270
     271
     272
     273=== `init_n_domains_pdaf` (init_n_domains_pdaf.F90) ===
     274
     275This routine is only used for localization. It is called during the analysis step before the loop over the local analysis domains is entered. It has to provide the number of local analysis domains. In case of a domain-decomposed model, the number of local analysis domain for the model sub-domain of the calling process has to be initialized.
     276
     277The interface for this routine is:
     278{{{
     279SUBROUTINE init_n_domains_pdaf(step, n_domains_p)
     280
     281  INTEGER, INTENT(in)  :: step        ! Current time step
     282  INTEGER, INTENT(out) :: n_domains_p ! Number of analysis domains for local model sub-domain
     283}}}
     284
     285Hints:
     286 * As a simple case, if the localization is only performed horizontally, the local analysis domains can be single vertical columns of the model grid. In this case, `n_domains_p` is simply the number of vertical columns in the process-local model sub-domain.
     287
     288
     289=== `init_dim_l_pdaf` (init_dim_l_pdaf.F90) ===
     290
     291This routine is only used for localization.
     292
     293The interface for this routine is:
     294{{{
     295SUBROUTINE init_dim_l_pdaf(step, domain_p, dim_l)
     296
     297  INTEGER, INTENT(in)  :: step        ! Current time step
     298  INTEGER, INTENT(in)  :: domain_p    ! Current local analysis domain
     299  INTEGER, INTENT(out) :: dim_l       ! Local state dimension
     300}}}
     301
     302The routine is called during the loop over the local analysis domains in the analysis step.
     303
     304It provides in `dim_l` the dimension of the state vector for the local analysis domain with index `domain_p` to PDAF.
     305
     306In the recommended implementation shown in the tutorial and template codes, there are two further initializations:
     3071. The routine has initialize the index array `id_lstate_in_pstate` containing the indices of the elements of the local state vector in the global (or domain-decomposed) state vector. Then it has to provide this array to PDAF by calling `PDAFlocal_set_indices` (see below).
     3082. The routine initializes an array `coords_l` containing the coordinates of the local analysis domain. This is shared with `U_init_dim_obs_l_pdafomi` via the module `mod_assimilation`.
     309
     310Hints:
     311  * The coordinates in `coords_l` have to describe one location in space that is used for localization to compute the distance from observations.
     312  * The coordinates in `coords_l` have the same units as those used for the observations
     313  * For geographic distance computations, the unit of the coordinates needs to be radian, thus (0, 2*pi) or (-pi,pi) for longitude and (-pi/2, pi/2) for latitude.
     314 * Any form of local domain is possible as long as it can be describe as a single location.
     315  * If the local domain is a single grid point, `dim_l` will be the number of model variables at this grid point.
     316  * The local analysis domain can also be a single vertical column of the model grid if observations are only horizontally distributed (a common situation with satellite data in the ocean).
     317   * In this case, `dim_l` will be the number of vertical grid points at this location times the number of model fields that exist in the vertical, plus possible variables at e.g. the surface.
     318   * In this case only the horizontal coordinates are used in `coords_l`.
     319
     320The index array `id_lstate_in_pstate` is an integer array in form of a one-dimensional vector. One initializes this vector by determining the indices of the elements of the local state vector in the global, or domain decomposed, state vector. After initializing `id_lstate_in_pstate`, one has to provided it to PDAF by calling `PDAFlocal_set_indices'. The interface interface is:
     321
     322{{{
     323SUBROUTINE PDAFlocal_set_indices(dim_l, id_lstate_in_pstate)
     324
     325  INTEGER, INTENT(in) :: dim_l                          ! Dimension of local state vector
     326  INTEGER, INTENT(in) :: id_lstate_in_pstate(dim_l)     ! Index array for mapping
     327}}}
     328
     329Hint for `id_lstate_in_pstate`:
     330 * The initialization of the index vector `id_lstate_to_pstate` is analogous to a loop that directly performs the initialization of a local state vector. However, here only the indices are stored.
     331 * See the [wiki:PDAFlocal_overview PDAFlocal overview page] for more information on the functionality of PDAFlocal.
     332
     333
     334=== `init_dim_obs_l_pdafomi` (callback_obs_pdafomi.F90) ===
     335
     336This routine is only used for localization. It is a call-back routine for PDAF-OMI that initializes the local observation vector. The routine calls a routine from the observation module for each observation type.
     337
     338See the [wiki:OMI_Callback_obs_pdafomi documentation on callback_obs_pdafomi.F90] for more information.
     339
     340
     341=== `prepoststep_pdaf` (prepoststep_ens_pdaf.F90) ===
     342
     343The routine has already been described for modifying the model for the ensemble integration and for inserting the analysis step.
     344
     345See the page on [ModifyModelforEnsembleIntegration#distribute_state_pdafdistribute_state_pdaf.F90 modifying the model code for the ensemble integration] for the description of this routine.
     346
     347
     348=== `next_observation_pdaf` (next_observation_pdaf.F90) ===
     349
     350This routine is independent of the filter algorithm used.
     351
     352See the page on [ModifyModelforEnsembleIntegration#distribute_state_pdafdistribute_state_pdaf.F90 modifying the model code for the ensemble integration] for the description of this routine.
     353
     354== Execution order of user-supplied routines ==
     355
     356The 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, `init_dim_obs_pdafomi` prepares an index array that provides the information for executing the observation operator in `obs_op_pdafomi`. How this information is initialized is described in the documentation of OMI.
     357
     358Before the analysis step is called the following routine is executed:
     359 1. [#collect_state_pdafcollect_state_pdaf.F90 collect_state_pdaf]
     360
     361The 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:
     362 1. [#prepoststep_pdafprepoststep_ens_pdaf.F90 prepoststep_pdaf] (Call to act on the forecast ensemble, called with negative value of the time step)
     363 1. [#init_dim_obs_pdafomicallback_obs_pdafomi.F90 init_dim_obs_pdafomi]
     364 1. [#obs_op_pdafomicallback_obs_pdafomi.F90 obs_op_pdafomi] (multiple calls, one for each ensemble member)
     365
     366Inside the analysis step the interative optimization is computed. This involves the repeated call of the routines:
     367 1. [#cvt_ens_pdafcvt_ens_pdaf.F90 cvt_ens_pdaf]
     368 1. [#obs_op_lin_pdafomicallback_obs_pdafomi.F90 obs_op_lin_pdafomi]
     369 1. [#obs_op_adj_pdafomicallback_obs_pdafomi.F90 obs_op_adj_pdafomi]
     370 1. [#cvt_adj_ens_pdafcvt_adj_ens_pdaf.F90 cvt_adj_ens_pdaf]
     371
     372After the iterative optimization the following routines are executes to complte the analysis step:
     373 1. [#cvt_ens_pdafcvt_pdaf.F90 cvt_pdaf] (Call to the parameterized part of the control vector transform to compute the final state vector increment)
     374 1. [#cvt_ens_pdafcvt_ens_pdaf.F90 cvt_ens_pdaf] (Call to the eensemble-part of the control vector transform to compute the final state vector increment)
     375 1. [#prepoststep_pdafprepoststep_ens_pdaf.F90 prepoststep_pdaf] (Call to act on the analysis ensemble, called with (positive) value of the time step)
     376
     377The iterative optimization abovve computes an updated ensemble mean state. Subsequently, the ensemble perturbations are updated using the LESTKF or ESTKF. The execution of the routines for these filters is described on the [wiki:ImplementAnalysisPDAF3Universal page on implementing the local filter analysis step] .
     378
     379In case of the routine `PDAF3_assimilate_en3dvar`, the following routines are executed after the analysis step:
     380 1. [#distribute_state_pdafdistribute_state_pdaf.F90 distribute_state_pdaf]
     381 1. [#next_observation_pdafnext_observation_pdaf.F90 next_observation_pdaf]