Changes between Initial Version and Version 1 of ImplementAnalysis_3DEnVar_PDAF23


Ignore:
Timestamp:
Jun 3, 2025, 1:55:53 PM (3 days ago)
Author:
lnerger
Comment:

--

Legend:

Unmodified
Added
Removed
Modified
  • ImplementAnalysis_3DEnVar_PDAF23

    v1 v1  
     1= Implementation of the Analysis Step for 3D Ensemble Var with OMI =
     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="OMI_ImplementationofAnalysisStep">Implementation of the analysis step with OMI</a></li>
     12<ol>
     13<li> <a href="ImplementFilterAnalysisOverview"> General overview for ensemble filters</a></li>
     14<ol>
     15<li><a href="ImplementAnalysisGlobal">Implementation for Global Filters</a></li>
     16<li><a href="ImplementAnalysisLocal">Implementation for Local Filters</a></li>
     17<li><a href="ImplementAnalysislenkfOmi">Implementation for LEnKF</a></li>
     18</ol>
     19<li> <a href="Implement3DVarAnalysisOverview"> General overview for 3D-Var methods</a></li>
     20<ol>
     21<li><a href="ImplementAnalysis_3DVar">Implementation for 3D-Var</a></li>
     22<li>Implementation for 3D Ensemble Var</li>
     23<li><a href="ImplementAnalysis_3DEnVar_untilPDAF221">Implementation for 3D Ensemble Var without PDAFlocal</li>
     24<li><a href="ImplementAnalysis_Hyb3DVar">Implementation for Hybrid 3D-Var</a></li>
     25</ol>
     26<li><a href="OMI_nondiagonal_observation_error_covariance_matrices">Using nondiagonal R-matrices</a></li>
     27<li><a href="PDAF_OMI_Overview">PDAF-OMI Overview</a></li>
     28</ol>
     29<li><a href="AddingMemoryandTimingInformation">Memory and timing information</a></li>
     30<li><a href="EnsembleGeneration">Ensemble Generation</a></li>
     31<li><a href="DataAssimilationDiagnostics">Diagnostics</a></li>
     32</ol>
     33</div>
     34}}}
     35
     36[[PageOutline(2-3,Contents of this page)]]
     37
     38== Overview ==
     39
     40This page describes the recommended implementation of the analysis step of local filters with OMI using the PDAFlocal interface that was introduced with PDAF V2.3. The older approach calling PDAFomi_assimilate_en3dvar_lestkf or PDAFomi_put_state_en3dvar_lestkf is documented on the page on [wiki:ImplementAnalysis_3DEnVar_untilPDAF221 Implementing the Analysis Step for 3D Ensemble Var with OMI without PDAFlocal (until V2.2.1 of PDAF)].
     41
     42There are genenerally three different variants: parameterized 3D-Var, 3D Ensemble Var, and hybrid (parameterized + ensemble) 3D-Var.
     43
     44This page describes the implementation of the analysis step for the 3D Ensemble Var using PDAF-OMI.
     45
     46For the analysis step of 3D-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 routine `PDAFomi_assimilate_en3dvar` in the fully-parallel implementation (or `PDAFomi_put_state_en3dvar` for the 'flexible' implementation) that was discussed before. With regard to the parallelization, all these routines (except `U_collect_state`, `U_distribute_state`, and `U_next_observation`) are executed by the filter processes (`filterpe=.true.`) only.
     47
     48For 3D Ensemble Var the ensemble is used to represent the background covariance matrix '''B'''. This ensemble perturbations need to be transformed by means of 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.
     49
     50For completeness we discuss here all user-supplied routines that are specified in the interface to `PDAFomi_assimilate_en3dvar`. 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.
     51
     52
     53== Analysis Routines ==
     54
     55The general aspects of the filter (or solver) specific routines `PDAF_assimilate_*` have been described on the page [ModifyModelforEnsembleIntegration Modification of the model code for the ensemble integration] and its sub-page on [InsertAnalysisStep inserting the analysis step]. The routine is used in the fully-parallel implementation variant of the data assimilation system. When the 'flexible' implementation variant is used, the routines `PDAF_put_state_*` is used as described further below. Here, we list the full interface of the routine. Subsequently, the user-supplied routines specified in the call is explained.
     56
     57There are two variants that either compute the transformataion of the ensemble transformation using the local LESTKF method, or the global ESTKF.
     58
     59=== `PDAFlocalomi_assimilate_en3dvar_lestkf` ===
     60
     61This routine is called for the case of transforming the ensemble perturbations using the local LESTKF.
     62
     63The interface is:
     64{{{
     65SUBROUTINE PDAFlocalomi_assimilate_en3dvar_lestkf(U_collect_state, U_distribute_state, &
     66                                 U_init_dim_obs_pdafomi, U_obs_op_pdafomi, &
     67                                 U_cvt_ens, U_cvt_adj_ens, U_obs_op_lin_pdafomi, U_obs_op_adj_pdafomi, &
     68                                 U_init_n_domains_p, U_init_dim_l, U_init_dim_obs_l_pdafomi, &
     69                                 U_prepoststep, U_next_observation, outflag)
     70}}}
     71with the following arguments:
     72 * [#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` as well as here.
     73 * [#U_distribute_statedistribute_state_pdaf.F90 U_distribute_state]:  The name of a user supplied routine that initializes the model fields from the array holding the ensemble of model state vectors.
     74 * [#U_init_dim_obs_pdafomicallback_obs_pdafomi.F90 U_init_dim_obs_pdafomi]: The name of the user-supplied routine that initializes the observation information and provides the size of observation vector
     75 * [#U_obs_op_pdafomicallback_obs_pdafomi.F90 U_obs_op_pdafomi]: The name of the user-supplied routine that acts as the observation operator on some state vector
     76 * [#U_cvt_enscvt_ens_pdaf.F90 U_cvt_ens]: 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.
     77 * [#U_cvt_adj_enscvt_adj_ens_pdaf.F90 U_cvt_adj_ens]: 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.
     78 * [#U_obs_op_pdafomicallback_obs_pdafomi.F90 U_obs_op_lin_pdafomi]: The name of the user-supplied routine that acts as the linearized observation operator on some state vector
     79 * [#U_obs_op_pdafomicallback_obs_pdafomi.F90 U_obs_op_lin_pdafomi]: The name of the user-supplied routine that acts as the adjoint observation operator on some state vector
     80 * [#U_init_n_domainsinit_n_domains_pdaf.F90 U_init_n_domains]: The name of the routine that provides the number of local analysis domains
     81 * [#U_init_dim_linit_dim_l_pdaf.F90 U_init_dim_l]: The name of the routine that provides the state dimension for a local analysis domain
     82 * [#U_init_dim_obs_l_pdafomicallback_obs_pdafomi.F90 U_init_dim_obs_l_pdafomi]: The name of the routine that initializes the size of the observation vector for a local analysis domain
     83 * [#U_prepoststepprepoststep_ens_pdaf.F90 U_prepoststep]: The name of the pre/poststep routine as in `PDAF_get_state`
     84 * [#U_next_observationnext_observation.F90 U_next_observation]: The name of a user supplied routine that initializes the variables `nsteps`, `timenow`, and `doexit`. The same routine is also used in `PDAF_get_state`.
     85 * `status`: The integer status flag. It is zero, if the routine is exited without errors.
     86
     87Note:
     88 * If your code shows a call to `PDAFomi_assimilate_en3dvar_lestkf`, it uses the implementation variant without PDAFlocal. This is documented on the page on [wiki:ImplementAnalysis_3DEnVar_untilPDAF221 Implementing the Analysis Step for 3D Ensemble Var with OMI without PDAFlocal (until V2.2.1 of PDAF)].
     89
     90=== `PDAFomi_assimilate_en3dvar_estkf` ===
     91
     92This routine is called for the case of transforming the ensemble perturbations using the global ESTKF. 
     93
     94The interface is:
     95{{{
     96SUBROUTINE PDAFomi_assimilate_en3dvar_estkf(U_collect_state, U_distribute_state, &
     97                                 U_init_dim_obs_pdafomi, U_obs_op_pdafomi, &
     98                                 U_cvt_ens, U_cvt_adj_ens, U_obs_op_lin_pdafomi, U_obs_op_adj_pdafomi, &
     99                                 U_prepoststep, U_next_observation, outflag)
     100}}}
     101with the following arguments:
     102 * [#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` as well as here.
     103 * [#U_distribute_statedistribute_state_pdaf.F90 U_distribute_state]:  The name of a user supplied routine that initializes the model fields from the array holding the ensemble of model state vectors.
     104 * [#U_init_dim_obs_pdafomicallback_obs_pdafomi.F90 U_init_dim_obs_pdafomi]: The name of the user-supplied routine that initializes the observation information and provides the size of observation vector
     105 * [#U_obs_op_pdafomicallback_obs_pdafomi.F90 U_obs_op_pdafomi]: The name of the user-supplied routine that acts as the observation operator on some state vector
     106 * [#U_cvt_enscvt_ens_pdaf.F90 U_cvt_ens]: 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.
     107 * [#U_cvt_adj_enscvt_adj_ens_pdaf.F90 U_cvt_adj_ens]: 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.
     108 * [#U_obs_op_pdafomicallback_obs_pdafomi.F90 U_obs_op_lin_pdafomi]: The name of the user-supplied routine that acts as the linearized observation operator on some state vector
     109 * [#U_obs_op_pdafomicallback_obs_pdafomi.F90 U_obs_op_lin_pdafomi]: The name of the user-supplied routine that acts as the adjoint observation operator on some state vector
     110 * [#U_prepoststepprepoststep_ens_pdaf.F90 U_prepoststep]: The name of the pre/poststep routine as in `PDAF_get_state`
     111 * [#U_next_observationnext_observation.F90 U_next_observation]: The name of a user supplied routine that initializes the variables `nsteps`, `timenow`, and `doexit`. The same routine is also used in `PDAF_get_state`.
     112 * `status`: The integer status flag. It is zero, if the routine is exited without errors.
     113
     114Note that the interface of `PDAFomi_assimilate_en3dvar_estkf` is identical to that of `PDAFomi_assimilate_3dvar` apart from using the routines `U_cvt_ens` and `U_cvt_adj_ens` in case of the ensemble variational method.
     115
     116
     117=== `PDAFlocalomi_put_state_en3dvar_lestkf` ===
     118
     119When the 'flexible' implementation variant is chosen for the assimilation system, the routine `PDAFomi_put_state_*` has to be used instead of `PDAFomi_assimilate_*`. The general aspects of the filter specific routines `PDAF_put_state_*` have been described on the page [ModifyModelforEnsembleIntegration Modification of the model code for the ensemble integration]. The interface of the routine is identical with that of `PDAF_assimilate_*` with the exception the specification of the user-supplied routines `U_distribute_state` and `U_next_observation` are missing.
     120
     121The interface when using one of the global filters is the following:
     122{{{
     123  SUBROUTINE PDAFomi_put_state_en3dvar_lestkf(U_collect_state, &
     124                                 U_init_dim_obs_pdafomi, U_obs_op_pdafomi, &
     125                                 U_cvt_ens, U_cvt_adj_ens, U_obs_op_lin_pdafomi, U_obs_op_adj_pdafomi, &
     126                                 U_init_n_domains_p, U_init_dim_l, U_init_dim_obs_l_pdafomi, &
     127                                 U_prepoststep, outflag)
     128}}}
     129
     130Note:
     131 * If your code shows a call to `PDAFomi_put_state_en3dvar_lestkf`, it uses the implementation variant without PDAFlocal. This is documented on the page on [wiki:ImplementAnalysis_3DEnVar_untilPDAF221 Implementing the Analysis Step for 3D Ensemble Var with OMI without PDAFlocal (until V2.2.1 of PDAF)].
     132
     133=== `PDAFomi_put_state_en3dvar_estkf` ===
     134
     135The interface of this routine is analogous to that of `PDAFomi_assimilate_en3dvar_estkf'. Thus it is identical to this routine with the exception the specification of the user-supplied routines `U_distribute_state` and `U_next_observation` are missing.
     136
     137The interface when using one of the global filters is the following:
     138{{{
     139  SUBROUTINE PDAFomi_put_state_en3dvar_estkf(U_collect_state, &
     140                                 U_init_dim_obs_pdafomi, U_obs_op_pdafomi, &
     141                                 U_cvt_ens, U_cvt_adj_ens, U_obs_op_lin_pdafomi, U_obs_op_adj_pdafomi, &
     142                                 U_prepoststep, outflag)
     143}}}
     144
     145== User-supplied routines ==
     146
     147Here all user-supplied routines are described that are required in the call to `PDAFomi_assimilate_3dvar`. For some of the generic routines, we link to the page on [ModifyModelforEnsembleIntegration modifying the model code for the ensemble integration].
     148
     149To 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.
     150
     151In 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.
     152
     153
     154=== `U_collect_state` (collect_state_pdaf.F90) ===
     155
     156This routine is independent of the filter algorithm used.
     157
     158See the page on [InsertAnalysisStep#U_collect_statecollect_state_pdaf.F90 inserting the analysis step] for the description of this routine.
     159
     160
     161=== `U_distribute_state` (distribute_state_pdaf.F90) ===
     162
     163This routine is independent of the filter algorithm used.
     164
     165See the page on [InsertAnalysisStep#U_distribute_statedistribute_state_pdaf.F90 inserting the analysis step] for the description of this routine.
     166
     167
     168=== `U_init_dim_obs_pdafomi` (callback_obs_pdafomi.F90) ===
     169
     170This is a call-back routine for PDAF-OMI initializing the observation information. The routine just calls a routine from the observation module for each observation type.
     171
     172See the [wiki:OMI_Callback_obs_pdafomi documentation on callback_obs_pdafomi.F90] for more information.
     173
     174
     175
     176=== `U_obs_op_pdafomi` (callback_obs_pdafomi.F90) ===
     177
     178This is a call-back routine for PDAF-OMI applying the observation operator to the state vector. The routine calls a routine from the observation module for each observation type.
     179
     180See the [wiki:OMI_Callback_obs_pdafomi documentation on callback_obs_pdafomi.F90] for more information.
     181
     182
     183
     184
     185=== `U_cvt_ens` (cvt_ens_pdaf.F90) ===
     186
     187The interface for this routine is:
     188{{{
     189SUBROUTINE cvt_ens_pdaf(iter, dim_p, dim_ens, dim_cv_ens_p, ens_p, cv_p, Vcv_p)
     190
     191  INTEGER, INTENT(in) :: iter               ! Iteration of optimization
     192  INTEGER, INTENT(in) :: dim_p              ! PE-local observation dimension
     193  INTEGER, INTENT(in) :: dim_ens            ! Ensemble size
     194  INTEGER, INTENT(in) :: dim_cv_ens_p       ! Dimension of control vector
     195  REAL, INTENT(in) :: ens_p(dim_p, dim_ens) ! PE-local ensemble
     196  REAL, INTENT(in) :: cv_p(dim_cv_ens_p)    ! PE-local control vector
     197  REAL, INTENT(inout) :: Vcv_p(dim_p)       ! PE-local state increment
     198}}}
     199
     200The routine is called during the analysis step during the iterative minimization of the cost function.
     201It 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.
     202
     203If 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.
     204
     205
     206=== `U_cvt_adj` (cvt_adj_pdaf.F90) ===
     207
     208The interface for this routine is:
     209{{{
     210SUBROUTINE cvt_adj_ens_pdaf(iter, dim_p, dim_ens, dim_cv_ens_p, ens_p, Vcv_p, cv_p)
     211
     212  INTEGER, INTENT(in) :: iter                ! Iteration of optimization
     213  INTEGER, INTENT(in) :: dim_p               ! PE-local observation dimension
     214  INTEGER, INTENT(in) :: dim_ens             ! Ensemble size
     215  INTEGER, INTENT(in) :: dim_cv_ens_p        ! PE-local dimension of control vector
     216  REAL, INTENT(in) :: ens_p(dim_p, dim_ens)  ! PE-local ensemble
     217  REAL, INTENT(in)    :: Vcv_p(dim_p)        ! PE-local input vector
     218  REAL, INTENT(inout) :: cv_p(dim_cv_ens_p)  ! PE-local result vector
     219}}}
     220
     221The routine is called during the analysis step during the iterative minimization of the cost function.
     222It 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.
     223
     224If 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.
     225
     226
     227
     228=== `U_obs_op_lin_pdafomi` (callback_obs_pdafomi.F90) ===
     229
     230This 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.
     231
     232See the [wiki:OMI_Callback_obs_pdafomi documentation on callback_obs_pdafomi.F90] for more information.
     233
     234
     235=== `U_obs_op_adj_pdafomi` (callback_obs_pdafomi.F90) ===
     236
     237This 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.
     238
     239See the [wiki:OMI_Callback_obs_pdafomi documentation on callback_obs_pdafomi.F90] for more information.
     240
     241
     242
     243=== `U_init_n_domains` (init_n_domains_pdaf.F90) ===
     244
     245The interface for this routine is:
     246{{{
     247SUBROUTINE init_n_domains(step, n_domains_p)
     248
     249  INTEGER, INTENT(in)  :: step        ! Current time step
     250  INTEGER, INTENT(out) :: n_domains_p ! Number of analysis domains for local model sub-domain
     251}}}
     252
     253The routine is called during the analysis step before the loop over the local analysis domains is entered.
     254It 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.
     255
     256Hints:
     257 * 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 local model sub-domain.
     258
     259
     260
     261
     262=== `U_init_dim_l` (init_dim_l_pdaf.F90) ===
     263
     264The interface for this routine is:
     265{{{
     266SUBROUTINE init_dim_l(step, domain_p, dim_l)
     267
     268  INTEGER, INTENT(in)  :: step        ! Current time step
     269  INTEGER, INTENT(in)  :: domain_p    ! Current local analysis domain
     270  INTEGER, INTENT(out) :: dim_l       ! Local state dimension
     271}}}
     272
     273The routine is called during the loop over the local analysis domains in the analysis step.
     274For PDAF it has to provide in `dim_l` the dimension of the state vector for the local analysis domain with index `domain_p`.
     275
     276In addition, for PDAFlocal the routine has to provide the index array containing the indices of the elements of the local state vector in the global (or domain-decomposed) state vector to PDAFlocal by calling `PDAFlocal_set_indices`. (in the template files, this array is called `id_lstate_in_pstate`)
     277
     278Hints:
     279 * For sharing through the module `mod_assimilation`, we further initialize an array `coords_l` containing the coordinates that describe the local domain.
     280  * These coordinates have to describe one location in space that is used in the OMI observation modules to compute the distance from observations.
     281  * The coordinates in `coords_l` have the same units as those used for the observations
     282  * 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.
     283 * Any form of local domain is possible as long as it can be describe as a single location.
     284  * If the local domain is a single grid point, `dim_l` will be the number of model variables at this grid point.
     285  * 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).
     286   * 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.
     287   * In this case only the horizontal coordinates are used in `coords_l`.
     288
     289The 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 PDAFlocal by calling `PDAFlocal_set_indices'. The interface interface is:
     290
     291{{{
     292SUBROUTINE PDAFlocal_set_indices(dim_l, id_lstate_in_pstate)
     293
     294  INTEGER, INTENT(in) :: dim_l                          ! Dimension of local state vector
     295  INTEGER, INTENT(in) :: id_lstate_in_pstate(dim_l)     ! Index array for mapping
     296}}}
     297
     298Hint for `id_lstate_in_pstate`:
     299 * 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.
     300 * See the [wiki:PDAFlocal_overview PDAFlocal overview page] for more information on the functionality of PDAFlocal.
     301
     302
     303=== `U_init_dim_obs_l_pdafomi` (callback_obs_pdafomi.F90) ===
     304
     305This 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.
     306
     307See the [wiki:OMI_Callback_obs_pdafomi documentation on callback_obs_pdafomi.F90] for more information.
     308
     309
     310=== `U_prepoststep` (prepoststep_ens_pdaf.F90) ===
     311
     312The routine has already been described for modifying the model for the ensemble integration and for inserting the analysis step.
     313
     314See the page on [InsertAnalysisStep#U_prepoststepprepoststep_ens_pdaf.F90 inserting the analysis step] for the description of this routine.
     315
     316
     317=== `U_next_observation` (next_observation_pdaf.F90) ===
     318
     319This routine is independent of the filter algorithm used.
     320
     321See the page on [InsertAnalysisStep#U_next_observationnext_observation_pdaf.F90 inserting the analysis step] for the description of this routine.
     322
     323
     324== Execution order of user-supplied routines ==
     325
     326The user-supplied routines are essentially executed in the order they are listed in the interface to `PDAFomi_assimilate_3dvar`. The order can be important as some routines can perform preparatory work for later routines. For example, `U_init_dim_obs_pdafomi` prepares an index array that provides the information for executing the observation operator in `U_obs_op_pdafomi`. How this information is initialized is described in the documentation of OMI.
     327
     328Before the analysis step is called the following routine is executed:
     329 1. [#U_collect_statecollect_state_pdaf.F90 U_collect_state]
     330
     331The 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:
     332 1. [#U_prepoststepprepoststep_ens_pdaf.F90 U_prepoststep] (Call to act on the forecast ensemble, called with negative value of the time step)
     333 1. [#U_init_dim_obs_pdafomicallback_obs_pdafomi.F90 U_init_dim_obs_pdafomi]
     334 1. [#U_obs_op_pdafomicallback_obs_pdafomi.F90 U_obs_op_pdafomi] (multiple calls, one for each ensemble member)
     335
     336Inside the analysis step the interative optimization is computed. This involves the repeated call of the routines:
     337 1. [#U_cvt_enscvt_ens_pdaf.F90 U_cvt_ens]
     338 1. [#U_obs_op_lin_pdafomicallback_obs_pdafomi.F90 U_obs_op_lin_pdafomi]
     339 1. [#U_obs_op_adj_pdafomicallback_obs_pdafomi.F90 U_obs_op_adj_pdafomi]
     340 1. [#U_cvt_adj_enscvt_adj_ens_pdaf.F90 U_cvt_adj_ens]
     341
     342After the iterative optimization the following routines are executes to complte the analysis step:
     343 1. [#U_cvt_enscvt_ens_pdaf.F90 U_cvt_ens] (Call to the control vector transform to compute the final state vector increment
     344 1. [#U_prepoststepprepoststep_ens_pdaf.F90 U_prepoststep] (Call to act on the analysis ensemble, called with (positive) value of the time step)
     345
     346The 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 for the LESTKF on the [wiki:ImplementAnalysisLocal page on implementing the local filter analysis step] and for the ESTKF on the [wiki:ImplementAnalysisGlobal page on implementing the global filter analysis step].
     347
     348In case of the routine `PDAFomi_assimilate_*`, the following routines are executed after the analysis step:
     349 1. [#U_distribute_statedistribute_state_pdaf.F90 U_distribute_state]
     350 1. [#U_next_observationnext_observation_pdaf.F90 U_next_observation]