Changes between Initial Version and Version 1 of OMI_nondiagonal_observation_error_covariance_matrices_PDAF23


Ignore:
Timestamp:
Jun 10, 2025, 12:00:52 PM (7 days ago)
Author:
lnerger
Comment:

--

Legend:

Unmodified
Added
Removed
Modified
  • OMI_nondiagonal_observation_error_covariance_matrices_PDAF23

    v1 v1  
     1= Using non-diagonal R matrices with OMI =
     2
     3{{{
     4#!html
     5<div class="wiki-toc">
     6<h4>PDAF-OMI Guide for PDAF2</h4>
     7<ol><li><a href="PDAF_OMI_Overview_PDAF23">Overview</a></li>
     8<li><a href="OMI_Callback_obs_pdafomi_PDAF23">callback_obs_pdafomi.F90</a></li>
     9<li><a href="OMI_observation_modules_PDAF23">Observation Modules</a></li>
     10<li><a href="OMI_observation_operators_PDAF23">Observation operators</a></li>
     11<li><a href="OMI_error_checking_PDAF23">Checking error status</a></li>
     12<li><a href="OMI_debugging_PDAF23">Debugging functionality</a></li>
     13<li><a href="OMI_ImplementationofAnalysisStep_PDAF23">Implementing the analysis step</a></li>
     14<li>Using nondiagonal R-matrices</li>
     15<li><a href="Porting_to_OMI_PDAF23">Porting an existing implemention to OMI</a></li>
     16<li><a href="PDAFomi_additional_functionality_PDAF23">Additional OMI Functionality</a></li>
     17</ol>
     18</div>
     19}}}
     20
     21[[PageOutline(2-3,Contents of this page)]]
     22
     23This feature was introduced with PDAF V2.3.
     24
     25|| The PDAF3 interface introduced by PDAF3 also provides support for non-diagonal R matrices. For details, see the [wiki:nondiagonal_observation_error_covariance_matrices_PDAF3 page on PDAF3 assimialtion interface for non-diagonal R matrices] ||
     26
     27The default mode of PDAF-OMI is to use a diagonal observation error covariance matrix **R** and specifying the observation error variances, i.e. the diagonalof **R** only. This is in line with the common choice in data assimilation to assume that observation errors are uncorrelated.
     28
     29However, there are also observation types with significant observation error correlations, which should be represented by a non-diagonal observation error covariance matrix. With PDAF V2.3 support for such nondiagonal **R** matrices was added to OMI.
     30
     31To support nondiagonal **R**-matrices consists in giving the user access to the routines that perform operations involving **R**. The operations dedepend on the filter type, e.g. in LESTKF and LETKF a product of some matrix with the inverse of **R** has to be computed, while in the traditional, perturbed observations, EnKF the matrix **R** has to be added to some other matrix. For the particle filter, the NETF and the hybrid filter LKNETF, the computation of the likelihood involves **R**.
     32
     33== Routines to perform the analysis step ==
     34
     35For using a nondiagonal **R** there is a modified variant of the PDAF-OMI routines `PDAFomi_assimilate_*` or `PDAFomi_put_state_*`. Compared to the case of a diagonal **R**, the routines are slightly less generic. Below we provide an overview of the routines including the links to the interface description of each. The last column gives the name of the additional routine(s) provided in the interface to implement the non-diagonal **R**-matrix.
     36
     37||= '''Global Filter'' =||= **diagonal R** =||= **non-diagonal R** =||= **additional routine(s)** =||
     38||= EnKF =|| [wiki:PDAFomi_assimilate_global][[BR]][wiki:PDAFomi_put_state_global] ||  [wiki:PDAFomi_assimilate_enkf_nondiagR][[BR]][wiki:PDAFomi_put_state_enkf_nondiagR]  || add_obs_err_pdafomi [[BR]]init_obscovar_pdafomi ||
     39||= ESTKF =|| [wiki:PDAFomi_assimilate_global][[BR]][wiki:PDAFomi_put_state_global] ||  [wiki:PDAFomi_assimilate_global_nondiagR][[BR]][wiki:PDAFomi_put_state_global_nondiagR]  || prodRinvA_pdafomi ||
     40||= ETKF =|| [wiki:PDAFomi_assimilate_global][[BR]][wiki:PDAFomi_put_state_global] ||  [wiki:PDAFomi_assimilate_global_nondiagR][[BR]][wiki:PDAFomi_put_state_global_nondiagR]  || prodRinvA_pdafomi ||
     41||= PF =|| [wiki:PDAFomi_assimilate_global][[BR]][wiki:PDAFomi_put_state_global] ||  [wiki:PDAFomi_assimilate_nonlin_nondiagR][[BR]][wiki:PDAFomi_put_state_nonlin_nondiagR]  || likelihood_pdafomi ||
     42||= NETF =|| [wiki:PDAFomi_assimilate_global][[BR]][wiki:PDAFomi_put_state_global] ||  [wiki:PDAFomi_assimilate_nonlin_nondiagR][[BR]][wiki:PDAFomi_put_state_nonlin_nondiagR]  || likelihood_pdafomi ||
     43||= SEIK =|| [wiki:PDAFomi_assimilate_global][[BR]][wiki:PDAFomi_put_state_global] ||  [wiki:PDAFomi_assimilate_global_nondiagR][[BR]][wiki:PDAFomi_put_state_global_nondiagR]  || prodRinvA_pdafomi ||
     44[[BR]]
     45||= '''Local Filter'' =||= **diagonal R** =||= **non-diagonal R** =||= **additional routine(s)** =||
     46||= LEnKF =|| [wiki:PDAFomi_assimilate_lenkf][[BR]][wiki:PDAFomi_put_state_lenkf] ||  [wiki:PDAFomi_assimilate_lenkf_nondiagR][[BR]][wiki:PDAFomi_put_state_lenkf_nondiagR]  || add_obs_err_pdafomi [[BR]]init_obscovar_pdafomi ||
     47||= LESTKF =|| [wiki:PDAFomi_assimilate_local][[BR]][wiki:PDAFomi_put_state_local] ||  [wiki:PDAFomi_assimilate_local_nondiagR][[BR]][wiki:PDAFomi_put_state_local_nondiagR]  || prodRinvA_l_pdafomi ||
     48||= LETKF =|| [wiki:PDAFomi_assimilate_local][[BR]][wiki:PDAFomi_put_state_local] ||  [wiki:PDAFomi_assimilate_local_nondiagR][[BR]][wiki:PDAFomi_put_state_local_nondiagR]  || prodRinvA_l_pdafomi ||
     49||= LSEIK =|| [wiki:PDAFomi_assimilate_local][[BR]][wiki:PDAFomi_put_state_local] ||  [wiki:PDAFomi_assimilate_local_nondiagR][[BR]][wiki:PDAFomi_put_state_local_nondiagR]  || prodRinvA_l_pdafomi ||
     50||= LNETF =|| [wiki:PDAFomi_assimilate_local][[BR]][wiki:PDAFomi_put_state_local] ||  [wiki:PDAFomi_assimilate_lnetf_nondiagR][[BR]][wiki:PDAFomi_put_state_lnetf_nondiagR]  || likelihood_l_pdafomi ||
     51||= LKNETF =|| [wiki:PDAFomi_assimilate_local][[BR]][wiki:PDAFomi_put_state_local] ||  [wiki:PDAFomi_assimilate_lknetf_nondiagR][[BR]][wiki:PDAFomi_put_state_lknetf_nondiagR]  ||likelihood_l_pdafomi[[BR]]likelihood_hyb_l_pdafomi[[BR]] prodRinvA_l_pdafomi[[BR]] prodRinvA_hyb_l_pdafomi  ||
     52[[BR]]
     53||= '''3D-Var'' =||= **diagonal R** =||= **non-diagonal R** =||= **additional routine(s)** =||
     54||= 3DVar =|| [wiki:PDAFomi_assimilate_3dvar][[BR]][wiki:PDAFomi_put_state_3dvar] || [wiki:PDAFomi_assimilate_3dvar_nondiagR][[BR]][wiki:PDAFomi_put_state_3dvar_nondiagR ] || prodRinvA_pdafomi ||
     55||= En3DVar ESTKF =|| [wiki:PDAFomi_assimilate_en3dvar_estkf][[BR]][wiki:PDAFomi_put_state_en3dvar_estkf] || [wiki:PDAFomi_assimilate_en3dvar_estkf_nondiagR][[BR]][wiki:PDAFomi_put_state_en3dvar_estkf_nondiagR] || prodRinvA_pdafomi ||
     56||= En3DVar LESTKF =|| [wiki:PDAFomi_assimilate_en3dvar_lestkf][[BR]][wiki:PDAFomi_put_state_en3dvar_lestkf] || [wiki:PDAFomi_assimilate_en3dvar_lestkf_nondiagR][[BR]][wiki:PDAFomi_put_state_en3dvar_lestkf_nondiagR] ||  prodRinvA_l_pdafomi ||
     57||= hyb3DVar ESTKF =|| [wiki:PDAFomi_assimilate_hyb3dvar_estkf][[BR]][wiki:PDAFomi_put_state_hyb3dvar_estkf] || [wiki:PDAFomi_assimilate_hyb3dvar_estkf_nondiagR][[BR]][wiki:PDAFomi_put_state_hyb3dvar_estkf_nondiagR] ||  prodRinvA_pdafomi ||
     58||= hyb3DVar LESTKF =|| [wiki:PDAFomi_assimilate_hyb3dvar_lestkf][[BR]][wiki:PDAFomi_put_state_hyb3dvar_lestkf] || [wiki:PDAFomi_assimilate_hyb3dvar_lestkf_nondiagR][[BR]][wiki:PDAFomi_put_state_hyb3dvar_lestkf_nondiagR] || prodRinvA_l_pdafomi ||
     59
     60== Call-back routines handling observations ==
     61
     62The tables above show the additional routine(s) that need to be implemented the nondiagonal **R** matrix for  the different filters and 3D-Var variants. Here we describe the interface of these routines. Generally the required routine(s) for the chosen DA method should be added to `callback_obs_pdafomi.F90` and to each observation module, where `thisobs` or `thisobs_l` can be included.
     63
     64Note that a mixed implementation is not possible. If one observation type has a non-diagonal **R** matrix, the additional routine(s) have to be implemented for all observations, even if an observation type uses a diagonal **R**. However, for observation wit a diagonal **R**, the corresponding PDAF-OMI routine to be called in the observation module can be called. Below we provide the information where this routine can be found.
     65
     66=== prodRinvA_pdafomi ===
     67
     68The interface for this routine is:
     69{{{
     70SUBROUTINE prodRinvA_pdafomi(step, dim_obs_p, rank, obs_p, A_p, C_p)
     71
     72  INTEGER, INTENT(in) :: step                ! Current time step
     73  INTEGER, INTENT(in) :: dim_obs_p           ! PE-local dimension of obs. vector
     74  INTEGER, INTENT(in) :: rank                ! Rank of initial covariance matrix
     75  REAL, INTENT(in)    :: obs_p(dim_obs_p)    ! PE-local vector of observations
     76  REAL, INTENT(in)    :: A_p(dim_obs_p,rank) ! Input matrix from analysis routine
     77  REAL, INTENT(out)   :: C_p(dim_obs_p,rank) ! Output matrix
     78}}}
     79
     80The routine computes the product of the inverse of the observation error covariance matrix with matrix `A_p`. For the ESTKF this matrix holds the observed part of the ensemble perturbations. The product has to be given as `C_p`.
     81
     82For 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.
     83
     84Hints:
     85 * The matrix `A_p` relates to the observations of all observation types. Thus, one needs to take the offset of an observation in the observation vector of all observation types into account. This offset is given by `thisobs%off_obs_f`.
     86 * 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.
     87 * 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.
     88 * The interface has a difference for ESTKF and ETKF: For ETKF the third argument is the ensemble size (`dim_ens`), while for the ESTKF it is the rank (`rank`) of the covariance matrix (usually ensemble size minus one). In addition, the second dimension of `A_p` and `C_p` has size `dim_ens` for ETKF, while it is `rank` for the ESTKF.  (Practically, one can usually ignore this difference as the fourth argument of the interface can be named arbitrarily in the routine.)
     89 * For diagonal **R** PDAF-OMI uses the routine `PDAFomi_prodRinvA` in `/src/PDAFomi_obs_f.F90` as the routine that is called within the observation module. This routine can serve as a template for an implementation by the user.
     90
     91=== likelihood_pdafomi ===
     92
     93The interface for this routine is:
     94{{{
     95SUBROUTINE likelihood_pdafomi(step, dim_obs_p, obs_p, residual, likely)
     96
     97  INTEGER, INTENT(in) :: step                ! Current time step
     98  INTEGER, INTENT(in) :: dim_obs_p           ! PE-local dimension of obs. vector
     99  REAL, INTENT(in)    :: obs_p(dim_obs_p)    ! PE-local vector of observations
     100  REAL, INTENT(in)    :: residual(dim_obs_p) ! Input vector holding the residual y-Hx
     101  REAL, INTENT(out)   :: likely              ! Output value of the likelihood
     102}}}
     103
     104The routine computes the likelihood of the observations for each ensemble member. The likelihood is computed from the observation-state residual according to the assumed observation error distribution. Commonly, the observation errors are assumed to be Gaussian distributed. In this case, the likelihood is '''exp(-0.5*(y-Hx)^T^*R^-1^*(y-Hx))'''.
     105
     106For a model with domain decomposition, `resid` contains the part of the matrix that resides on the model sub-domain of the calling process. The likelihood has to be computed for the global state vector. Thus some parallel communication might be required to complete the computation.
     107
     108Hints:
     109 * The matrix `residual` relates to the observations of all observation types. Thus, one needs to take the offset of an observation in the observation vector of all observation types into account. This offset is given by `thisobs%off_obs_f`.
     110 * The routine is very similar to the routine [wiki:U_prodRinvA]. The main addition is the computation of the likelihood after computing '''R^-1^*(y-Hx)''', which corresponds to '''R^-1^*A_p''' in [wiki:U_prodRinvA].
     111 * The information about the inverse observation error covariance matrix has to be provided by the user. Possibilities are to read this information from a file, or to use a Fortran module that holds this information, which one could already prepare in init_pdaf.
     112 * The routine does not require that the product is implemented as a real matrix-vector 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 inverse diagonal with the vector `resid` has to be implemented.
     113 * 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.
     114 * For diagonal **R** PDAF-OMI uses the routine `PDAFomi_likelihood` in `/src/PDAFomi_obs_f.F90` as the routine that is called within the observation module. This routine can serve as a template for an implementation by the user.
     115
     116=== add_obs_err_pdafomi ===
     117
     118The interface for this routine is:
     119{{{
     120SUBROUTINE add_obs_err_pdafomi(step, dim_obs, C)
     121
     122  INTEGER, INTENT(in) :: step                ! Current time step
     123  INTEGER, INTENT(in) :: dim_obs             ! Dimension of obs. vector
     124  REAL, INTENT(inout) :: C(dim_obs, dim_obs) ! Matrix to that the observation
     125                                             !    error covariance matrix is added
     126}}}
     127
     128During 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`.
     129
     130The operation is for the global observation space. Thus, it is independent of whether the filter is executed with or without parallelization.
     131
     132Hints:
     133 * Matrix C relates to the observations of all types. Thus, for a single obsevation type one has to take the offset of this observtion type in the observation vector into account. For this, PDAF-OMI initializes an array `obsdims`, which can also be used by a user-implemented routine. See the routine `PDAFomi_add_obs_error` in `/src/PDAFomi_obs_f.F90` for how this can be implemented.
     134 * 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.
     135 * For diagonal **R**, PDAF-OMI uses the routine `PDAFomi_add_obs_error` in `/src/PDAFomi_obs_f.F90` as the routine that is called within the observation module. This routine can serve as a template for an implementation by the user.
     136
     137=== init_obscovar_pdafomi ===
     138
     139The interface for this routine is:
     140{{{
     141SUBROUTINE init_obscovar_pdafomi(step, dim_obs, dim_obs_p, covar, m_state_p, &
     142     isdiag)
     143
     144  INTEGER, INTENT(in) :: step                ! Current time step
     145  INTEGER, INTENT(in) :: dim_obs             ! Dimension of observation vector
     146  INTEGER, INTENT(in) :: dim_obs_p           ! PE-local dimension of observation vector
     147  REAL, INTENT(out) :: covar(dim_obs, dim_obs) ! Observation error covariance matrix
     148  REAL, INTENT(in)  :: m_state_p(dim_obs_p)  ! PE-local observation vector
     149  LOGICAL, INTENT(out) :: isdiag             ! Whether the observation error covar. matrix is diagonal
     150}}}
     151
     152The routine 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.
     153
     154The operation is for the global observation space. Thus, it is independent of whether the filter is executed with or without parallelization.
     155
     156Hints:
     157 * Matrix `covar` relates to the observations of all types. Thus, for a single obsevation type one has to take the offset of this observtion type in the observation vector into account. For this, PDAF-OMI initializes an array `obsdims`, which can also be used by a user-implemented routine. See the routine `PDAFomi_init_obscovar` in `/src/PDAFomi_obs_f.F90` for how this can be implemented.
     158 * 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.
     159 * For diagonal **R**, PDAF-OMI uses the routine `PDAFomi_init_obscovar` in `/src/PDAFomi_obs_f.F90` as the routine that is called within the observation module. This routine can serve as a template for an implementation by the user.
     160
     161=== prodRinvA_l_pdafomi ===
     162
     163The interface for this routine is:
     164{{{
     165SUBROUTINE prodRinvA_l_pdafomi(domain_p, step, dim_obs_l, rank, obs_l, A_l, C_l)
     166
     167  INTEGER, INTENT(in) :: domain_p             ! Current local analysis domain
     168  INTEGER, INTENT(in) :: step                 ! Current time step
     169  INTEGER, INTENT(in) :: dim_obs_l            ! Dimension of local observation vector
     170  INTEGER, INTENT(in) :: rank                 ! Rank of initial covariance matrix
     171  REAL, INTENT(in)    :: obs_l(dim_obs_l)     ! Local vector of observations
     172  REAL, INTENT(inout) :: A_l(dim_obs_l, rank) ! Input matrix from analysis routine
     173  REAL, INTENT(out)   :: C_l(dim_obs_l, rank) ! Output matrix
     174}}}
     175
     176The routine is called during the loop over the local analysis domains. In the algorithm, the product of the inverse of the observation error covariance matrix with some matrix has to be computed. For the SEIK filter this matrix holds the observed part of the ensemble perturbations for the local analysis domain of index `domain_p`. The matrix is provided as `A_l`. The product has to be given as `C_l`.
     177
     178This routine is also the place to perform observation localization. To initialize a vector of weights, the routine `PDAF_local_weight` can be called. The procedure is used in the example implementation and also demonstrated in the template routine.
     179
     180Hints:
     181 * The routine is a local variant of the routine `U_prodRinvA`. Thus if that routine has been implemented before, it can be adapted here for the local filter.
     182 * The matrix `A_p` relates to the observations of all observation types. Thus, one needs to take the offset of an observation in the observation vector of all observation types into account. This offset is given by `thisobs_l%off_obs_l`.
     183 * 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_l` has to be implemented.
     184 * The observation vector `obs_l` is provided through the interface for cases where the observation error variance is relative to the actual value of the observations.
     185 * The interface has a difference for LESTKF and LETKF: For LETKF the third argument is the ensemble size (`dim_ens`), while for LESTKF it is the rank (`rank`) of the covariance matrix (usually ensemble size minus one). In addition, the second dimension of `A_l` and `C_l` has size `dim_ens` for LETKF, while it is `rank` for LESTKF.  (Practically, one can usually ignore this difference as the fourth argument of the interface can be named arbitrarily in the routine.)
     186 * For diagonal **R** PDAF-OMI uses the routine `PDAFomi_prodRinvA_l` in `/src/PDAFomi_obs_k.F90` as the routine that is called within the observation module. This routine can serve as a template for an implementation by the user.
     187
     188=== likelihood_l_pdafomi ===
     189
     190The interface for this routine is:
     191{{{
     192SUBROUTINE likelihood_l_pdafomi(domain_p, step, dim_obs_l, obs_l, resid_l, likely_l)
     193
     194  INTEGER, INTENT(in) :: domain_p             ! Current local analysis domain
     195  INTEGER, INTENT(in) :: step                 ! Current time step
     196  INTEGER, INTENT(in) :: dim_obs_l            ! Dimension of local observation vector
     197  REAL, INTENT(in)    :: obs_l(dim_obs_l)     ! Local vector of observations
     198  REAL, INTENT(inout) :: resid_l(dim_obs_l)   ! Input vector holding the local residual y-Hx
     199  REAL, INTENT(out)   :: likely_l(dim_obs_l)  ! Output value of the likelihood
     200}}}
     201
     202The routine is called during the loop over the local analysis domains. In the NETF, as in other particle filters, the likelihood of the local observations has to be computed for each ensemble member. The likelihood is computed from the observation-state residual according to the assumed observation error distribution. Commonly, the observation errors are assumed to be Gaussian distributed. In this case, the likelihood is '''exp(-0.5*(y-Hx)^T^*R^-1^*(y-Hx))'''.
     203
     204This routine is also the place to perform observation localization. To initialize a vector of weights, the routine `PDAF_local_weight` can be called. The procedure is used in the example implementation and also demonstrated in the template routine.
     205
     206Hints:
     207* The routine is a local variant of the routine `U_likelihood`. Thus if that routine has been implemented before, it can be adapted here for the local filter.
     208 * The matrix `resid_l` relates to the observations of all observation types. Thus, one needs to take the offset of an observation in the observation vector of all observation types into account. This offset is given by `thisobs_l%off_obs_l`.
     209 * The routine is very similar to the routine [wiki:U_prodRinvA_l]. The main addition is the computation of the likelihood after computing '''R^-1^*(y-Hx)''', which corresponds to '''R^-1^*A_p''' in [wiki:U_prodRinvA_l].
     210 * The information about the inverse observation error covariance matrix has to be provided by the user. Possibilities are to read this information from a file, or to use a Fortran module that holds this information, which one could already prepare in init_pdaf.
     211 * The routine does not require that the product is implemented as a real matrix-vector 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 inverse diagonal with the vector `resid` has to be implemented.
     212 * The observation vector `obs_l` is provided through the interface for cases where the observation error variance is relative to the actual value of the observations.
     213 * For diagonal **R** PDAF-OMI uses the routine `PDAFomi_likelihood_l` in `/src/PDAFomi_obs_l.F90` as the routine that is called within the observation module. This routine can serve as a template for an implementation by the user.
     214
     215=== prodRinvA_hyb_l_pdafomi ===
     216
     217The interface for this routine is:
     218{{{
     219SUBROUTINE prodRinvA_hyb_l_pdafomi(domain_p, step, dim_obs_l, dim_ens, obs_l, gamma, A_l, C_l)
     220
     221  INTEGER, INTENT(in) :: domain_p             ! Current local analysis domain
     222  INTEGER, INTENT(in) :: step                 ! Current time step
     223  INTEGER, INTENT(in) :: dim_obs_l            ! Dimension of local observation vector
     224  INTEGER, INTENT(in) :: dim_ens              ! Ensemble size
     225  REAL, INTENT(in)    :: obs_l(dim_obs_l)     ! Local vector of observations
     226  REAL, INTENT(in)    :: gamma                ! Hybrid weight
     227  REAL, INTENT(inout) :: A_l(dim_obs_l, dim_ens) ! Input matrix from analysis routine
     228  REAL, INTENT(out)   :: C_l(dim_obs_l, dim_ens) ! Output matrix
     229}}}
     230
     231The routine is called during the loop over the local analysis domains. In the algorithm, the product of the inverse of the observation error covariance matrix with some matrix has to be computed. This matrix holds the observed part of the ensemble perturbations for the local analysis domain of index `domain_p`. The matrix is provided as `A_l`. The product has to be given as `C_l`.
     232
     233This routine is also the place to perform observation localization. To initialize a vector of weights, the routine `PDAF_local_weight` can be called. The procedure is used in the example implementation and also demonstrated in the template routine.
     234
     235The routine is analogous to `prodRinvA_l_pdafomi, but also has to apply the hybrid weight `gamma` of the hybrid filter LKNETF. This is a simple multiplication with the input value in the loop where `C_l` is initialized.
     236
     237Hints:
     238 * This routine is a simple extension of `prodRinvA_l_pdafomi`. One can implement the hybrid variant by copying this routine and adapting it. `gamma` is computed inside PDAF and provided to the routine.
     239 * Please see the further implementation hints for `prodRinvA_l_pdafomi`.
     240 * For diagonal **R** PDAF-OMI uses the routine `PDAFomi_prodRinvA_hyb_l` in `/src/PDAFomi_obs_l.F90` as the routine that is called within the observation module. This routine can serve as a template for an implementation by the user.
     241
     242=== likelihood_hyb_l_pdafomi ===
     243
     244The interface for this routine is:
     245{{{
     246SUBROUTINE likelihood_hyb_l_pdafomi(domain_p, step, dim_obs_l, obs_l, resid_l, gamma, likely_l)
     247
     248  INTEGER, INTENT(in) :: domain_p             ! Current local analysis domain
     249  INTEGER, INTENT(in) :: step                 ! Current time step
     250  INTEGER, INTENT(in) :: dim_obs_l            ! Dimension of local observation vector
     251  REAL, INTENT(in)    :: obs_l(dim_obs_l)     ! Local vector of observations
     252  REAL, INTENT(inout) :: resid_l(dim_obs_l)   ! Input vector holding the local residual y-Hx
     253  REAL, INTENT(in)    :: gamma                ! Hybrid weight
     254  REAL, INTENT(out)   :: likely_l(dim_obs_l)  ! Output value of the likelihood
     255}}}
     256
     257This routine is a variant for `U_likelihood_l_pdafomi`. See the description of this routine for its functionality. The routine is called during the loop over the local analysis domains. The likelihood of the local observations has to be computed for each ensemble member. The likelihood is computed from the observation-state residual according to the assumed observation error distribution. Commonly, the observation errors are assumed to be Gaussian distributed. In this case, the likelihood is '''exp(-0.5*(y-Hx)^T^*R^-1^*(y-Hx))'''.
     258
     259This routine is also the place to perform observation localization. To initialize a vector of weights, the routine `PDAF_local_weight` can be called. The procedure is used in the example implementation and also demonstrated in the template routine.
     260
     261The routine also has to apply the hybrid weight `gamma`. This is a simple multiplication with `1-gamma` in the loop where `Rinvresid_l` is initialized.
     262
     263Hints:
     264 * This routine is a simple extension of `U_likelihood_l. One can implement the hybrid variant by copying this routine and adapting it. `gamma` is computed inside PDAF and provided to the routine.
     265 * Please see the further implementation hints for `prodRinvA_l_pdafomi`.
     266 * For diagonal **R** PDAF-OMI uses the routine `PDAFomi_likelihood_hyb_l` in `/src/PDAFomi_obs_l.F90` as the routine that is called within the observation module. This routine can serve as a template for an implementation by the user.
     267
     268== Hints ==
     269
     270For now, we do not have an example code that we could provide. However, for the implementation the follonig hints might be useful:
     271 * The matrix **R** could be constructed e.g. from decorrelation length scales or by direct implementation.
     272 * The routine [wiki:PDAF_correlation_function] can be used to obtain values of a correlation function.
     273 * Th eroutine [wiki:PDAFomi_observation_localization_weights] can be used to obtain the values of the localization weights.
     274 * The matrix **R** cannot be stored fully for larger numbers of observations. E.g. for 100 000 observations, the full matric **R** would have 10^10^ entries, which corresponds to 80 GB memory at double precision. Thus, one would need compressed ways of storing the matrix or one can implement the products with the inverse of **R** or the addition of **R** in form of oeprations.
     275 * For the localized operations in `prodRinvA_l_pdafomi` and `likelihood_l_pdafomi` it should be possible to handle the part of **R** that relates only to the local observations used in a local analysis domain. However, the inverse of this local matrix part is not identical to the local part of the global inverse. The different might be mitigated by applying the localization weight to **R**.
     276 * Examples of implementations are also available without PDAF-OMI, e.g. in the tutorial and in /models/lorenz96
     277.