= Implementation of the Analysis Step for Hybrid 3D-Var without using OMI = {{{ #!html

Implementation Guide

  1. Main page
  2. Adaptation of the parallelization
  3. Initialization of PDAF
  4. Modifications for ensemble integration
  5. Implementation of the analysis step
    1. Implementation for ESTKF
    2. Implementation for LESTKF
    3. Implementation for ETKF
    4. Implementation for LETKF
    5. Implementation for SEIK
    6. Implementation for LSEIK
    7. Implementation for SEEK
    8. Implementation for EnKF
    9. Implementation for LEnKF
    10. Implementation for NETF
    11. Implementation for LNETF
    12. Implementation for PF
    13. Implementation for 3D-Var
    14. Implementation for 3D Ensemble Var
    15. Implementation for Hybrid 3D-Var
  6. Memory and timing information
  7. Ensemble Generation
  8. Diagnostics
}}} [[PageOutline(2-3,Contents of this page)]] || This page describes the implementation of the analysis step without using PDAF-OMI. Please see the [wiki:ImplementationofAnalysisStep page on the analysis with OMI] for the more modern and efficient implementation variant using PDAF-OMI. || == Overview == With Version 2.0 with introduced 3D variational assimilation methods to PDAF. There are genenerally three different variants: parameterized 3D-Var, 3D Ensemble Var, and hybrid (parameterized + ensemble) 3D-Var. This page describes the implementation of the analysis step for the hybrid 3D-Var in the classical way (without using PDAF-OMI). For the analysis step of 3D-Var we need different operations related to the observations. 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 as this procedure should simplify the implementation. The names of the required routines are specified in the call to the routine `PDAF_assimilate_3dvar` in the fully-parallel implementation (or `PDAF_put_state_3dvar` for the 'flexible' implementation) described below. With regard to the parallelization, all these routines (except `U_collect_state`) are executed by the filter processes (`filterpe=.true.`) only. For Hybrid 3D-Var the background covariance matrix '''B''' is represented by a combination of a parameterized covariance matrix with a covariance matrix part represented by the ensemble. In practive this means that in the square root of '''B''' one concatenates parameterized and ensemble columns. The 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. For completeness we discuss here all user-supplied routines that are specified in the interface to `PDAF_assimilate_hyb3dvar_*` and `PDAF_put_state_hyb3dvar_*`. 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. == Analysis Routines == The 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. There are two variants that either compute the transformataion of the ensemble transformation using the local LESTKF method, or the global ESTKF. === `PDAF_assimilate_hyb3dvar_lestkf` === This routine is called for the case of transforming the ensemble perturbations using the local LESTKF. The interface is: {{{ SUBROUTINE PDAF_assimilate_hyb3dvar_lestkf(U_collect_state, U_distribute_state, & U_init_dim_obs, U_obs_op, U_init_obs, U_prodRinvA, & U_cvt_ens, U_cvt_adj_ens, U_cvt, U_cvt_adj, U_obs_op_lin, U_obs_op_adj, & U_init_dim_obs_f, U_obs_op_f, U_init_obs_f, U_init_obs_l, U_prodRinvA_l, & U_init_n_domains_p, U_init_dim_l, U_init_dim_obs_l, U_g2l_state, U_l2g_state, & U_g2l_obs, U_init_obsvar, U_init_obsvar_l, & U_prepoststep, U_next_observation, outflag) }}} with the following arguments: * [#U_collect_statecollect_state_pdaf.F90 U_collect_state]: The name of the user-supplied routine that initializes a state vector from the array holding the ensemble of model states from the model fields. This is basically the inverse operation to `U_distribute_state` used in `PDAF_get_state` as well as here. * [#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. * [#U_init_dim_obsinit_dim_obs_pdaf.F90 U_init_dim_obs]: The name of the user-supplied routine that provides the size of observation vector * [#U_obs_opobs_op_pdaf.F90 U_obs_op]: The name of the user-supplied routine that acts as the observation operator on some state vector * [#U_init_obsinit_obs_pdaf.F90 U_init_obs]: The name of the user-supplied routine that initializes the vector of observations * [#U_prodRinvAprodrinva_pdaf.F90 U_prodRinvA]: The name of the user-supplied routine that computes the product of the inverse of the observation error covariance matrix with some matrix provided to the routine by PDAF. This operation occurs during the analysis step of the ETKF. * [#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. * [#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. * [#U_cvtcvt_pdaf.F90 U_cvt]: The name of the user-supplied routine that applies the control-vector transformation (square-root of the B-matrix) on some control vector to obtain a state vector. * [#U_cvt_adjcvt_adj_pdaf.F90 U_cvt_adj]: The name of the user-supplied routine that applies the adjoint control-vector transformation (with square-root of the B-matrix) on some state vector to obtain the control vector. * [#U_obs_op_linobs_op_lin_pdaf.F90 U_obs_op_lin]: The name of the user-supplied routine that acts as the linearized observation operator on some state vector * [#U_obs_op_adjobs_op_adj_pdaf.F90 U_obs_op_adj]: The name of the user-supplied routine that acts as the adjoint observation operator on some state vector * [#U_init_dim_obs_finit_dim_obs_f_pdaf.F90 U_init_dim_obs_f]: The name of the user-supplied routine that provides the size of the full observation vector * [#U_obs_op_fobs_op_f_pdaf.F90 U_obs_op_f]: The name of the user-supplied routine that acts as the full observation operator on some state vector * [#U_init_obs_finit_obs_f_pdaf.F90 U_init_obs_f]: The name of the user-supplied routine that initializes the full vector of observations * [#U_init_obs_linit_obs_l_pdaf.F90 U_init_obs_l]: The name of the user-supplied routine that initializes the vector of observations for a local analysis domain * [#U_prepoststepprepoststep_ens_pdaf.F90 U_prepoststep]: The name of the pre/poststep routine as in `PDAF_get_state` * [#U_prodRinvA_lprodrinva_l_pdaf.F90 U_prodRinvA_l]: The name of the user-supplied routine that computes the product of the inverse of the observation error covariance matrix with some matrix provided to the routine by PDAF. * [#U_init_n_domains_pinit_n_domains_pdaf.F90 U_init_n_domains]: The name of the routine that provides the number of local analysis domains * [#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 * [#U_init_dim_obs_linit_dim_obs_l_pdaf.F90 U_init_dim_obs_l]: The name of the routine that initializes the size of the observation vector for a local analysis domain * [#U_g2l_stateg2l_state_pdaf.F90 U_g2l_state]: The name of the routine that initializes a local state vector from the global state vector * [#U_l2g_statel2g_state_pdaf.F90 U_l2g_state]: The name of the routine that initializes the corresponding part of the global state vector from the the provided local state vector * [#U_g2l_obsg2l_obs_pdaf.F90 U_g2l_obs]: The name of the routine that initializes a local observation vector from a full observation vector * [#U_init_obsvarinit_obsvar_pdaf.F90 U_init_obsvar]: The name of the user-supplied routine that provides a global mean observation error variance (This routine will only be executed, if an adaptive forgetting factor is used) * [#U_init_obsvar_linit_obsvar_l_pdaf.F90 U_init_obsvar_l]: The name of the user-supplied routine that provides a mean observation error variance for the local analysis domain (This routine will only be executed, if a local adaptive forgetting factor is used) * [#U_prepoststepprepoststep_ens_pdaf.F90 U_prepoststep]: The name of the pre/poststep routine as in `PDAF_get_state` * [#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`. * `status`: The integer status flag. It is zero, if the routine is exited without errors. === `PDAF_assimilate_hyb3dvar_estkf` === This routine is called for the case of transforming the ensemble perturbations using the global ESTKF. The interface is: {{{ SUBROUTINE PDAF_assimilate_hyb3dvar_estkf(U_collect_state, U_distribute_state, & U_init_dim_obs, U_obs_op, U_init_obs, U_prodRinvA, & U_cvt_ens, U_cvt_adj_ens, U_cvt, U_cvt_adj, U_obs_op_lin, U_obs_op_adj, & U_init_obsvar, U_prepoststep, U_next_observation, outflag) }}} with the following arguments: * [#U_collect_statecollect_state_pdaf.F90 U_collect_state]: The name of the user-supplied routine that initializes a state vector from the array holding the ensemble of model states from the model fields. This is basically the inverse operation to `U_distribute_state` used in `PDAF_get_state` as well as here. * [#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. * [#U_init_dim_obsinit_dim_obs_pdaf.F90 U_init_dim_obs]: The name of the user-supplied routine that provides the size of observation vector * [#U_obs_opobs_op_pdaf.F90 U_obs_op]: The name of the user-supplied routine that acts as the observation operator on some state vector * [#U_init_obsinit_obs_pdaf.F90 U_init_obs]: The name of the user-supplied routine that initializes the vector of observations * [#U_prodRinvAprodrinva_pdaf.F90 U_prodRinvA]: The name of the user-supplied routine that computes the product of the inverse of the observation error covariance matrix with some matrix provided to the routine by PDAF. This operation occurs during the analysis step of the ETKF. * [#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. * [#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. * [#U_cvtcvt_pdaf.F90 U_cvt]: The name of the user-supplied routine that applies the control-vector transformation (square-root of the B-matrix) on some control vector to obtain a state vector. * [#U_cvt_adjcvt_adj_pdaf.F90 U_cvt_adj]: The name of the user-supplied routine that applies the adjoint control-vector transformation (with square-root of the B-matrix) on some state vector to obtain the control vector. * [#U_obs_op_linobs_op_lin_pdaf.F90 U_obs_op_lin]: The name of the user-supplied routine that acts as the linearized observation operator on some state vector * [#U_obs_op_adjobs_op_adj_pdaf.F90 U_obs_op_adj]: The name of the user-supplied routine that acts as the adjoint observation operator on some state vector * [#U_init_obsvarinit_obsvar_pdaf.F90 U_init_obsvar]: The name of the user-supplied routine that provides a mean observation error variance to PDAF (This routine will only be executed, if an adaptive forgetting factor is used) * [#U_prepoststepprepoststep_ens_pdaf.F90 U_prepoststep]: The name of the pre/poststep routine as in `PDAF_get_state` * [#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`. * `status`: The integer status flag. It is zero, if the routine is exited without errors. === `PDAF_put_state_hyb3dvar_lestkf` === When the 'flexible' implementation variant is chosen for the assimilation system, the routine `PDAF_put_state_*` has to be used instead of `PDAF_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. The interface when using one of the global filters is the following: {{{ SUBROUTINE PDAF_put_state_hyb3dvar_lestkf(U_collect_state, & U_init_dim_obs, U_obs_op, U_init_obs, U_prodRinvA, & U_cvt_ens, U_cvt_adj_ens, U_cvt, U_cvt_adj, U_obs_op_lin, U_obs_op_adj, & U_init_dim_obs_f, U_obs_op_f, U_init_obs_f, U_init_obs_l, U_prodRinvA_l, & U_init_n_domains_p, U_init_dim_l, U_init_dim_obs_l, U_g2l_state, U_l2g_state, & U_g2l_obs, U_init_obsvar, U_init_obsvar_l, & U_prepoststep, outflag) }}} === `PDAF_put_state_hyb3dvar_estkf` === The interface of this routine is analogous to that of `PDAF_assimilate_hyb3dvar_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. The interface when using one of the global filters is the following: {{{ SUBROUTINE PDAF_put_state_hyb3dvar_estkf(U_collect_state, & U_init_dim_obs, U_obs_op, U_init_obs, U_prodRinvA, & U_cvt_ens, U_cvt_adj_ens, U_cvt, U_cvt_adj, U_obs_op_lin, U_obs_op_adj, & U_init_obsvar, U_prepoststep, outflag) }}} == User-supplied routines == Here all user-supplied routines are described that are required in the calls to `PDAF_assimilate_hyb3dvar_*` and `PDAF_put_state_hyb3dvar_*`. For some of the generic routines, we link to the page on [ModifyModelforEnsembleIntegration modifying the model code for the ensemble integration]. To indicate user-supplied routines we use the prefix `U_`. In the tutorials in `tutorial/` and in the template directory `templates/` these routines exist without the prefix, but with the extension `_pdaf`. The files are named correspondingly. In the section titles below we provide the name of the template file in parentheses. In the subroutine interfaces some variables appear with the suffix `_p`. This suffix indicates that the variable is particular to a model sub-domain, if a domain decomposed model is used. Thus, the value(s) in the variable will be different for different model sub-domains. === `U_collect_state` (collect_state_pdaf.F90) === This routine is independent of the filter algorithm used. See the page on [InsertAnalysisStep#U_collect_statecollect_state_pdaf.F90 inserting the analysis step] for the description of this routine. === `U_distribute_state` (distribute_state_pdaf.F90) === This routine is independent of the filter algorithm used. See the page on [InsertAnalysisStep#U_distribute_statedistribute_state_pdaf.F90 inserting the analysis step] for the description of this routine. === `U_init_dim_obs` (init_dim_obs_pdaf.F90) === This routine is used by all global filter algorithms (SEEK, SEIK, EnKF, ETKF, NETF, PF) and the 3D-Var methods. The interface for this routine is: {{{ SUBROUTINE init_dim_obs(step, dim_obs_p) INTEGER, INTENT(in) :: step ! Current time step INTEGER, INTENT(out) :: dim_obs_p ! Dimension of observation vector }}} The routine is called at the beginning of each analysis step. It has to initialize the size `dim_obs_p` of the observation vector according to the current time step. Without parallelization `dim_obs_p` will be the size for the full model domain. When a domain-decomposed model is used, `dim_obs_p` will be the size of the observation vector for the sub-domain of the calling process. Some hints: * It can be useful to not only determine the size of the observation vector at this point. One can also already gather information about the locations of the observations, which will be used later, e.g. to implement the observation operator. An array for the locations can be defined in a module like `mod_assimilation` of the example implementation. === `U_obs_op` (obs_op_pdaf.F90) === This routine is used by all global filter algorithms (SEEK, SEIK, EnKF, ETKF, NETF, PF) and the 3D-Var methods. The interface for this routine is: {{{ SUBROUTINE obs_op(step, dim_p, dim_obs_p, state_p, m_state_p) INTEGER, INTENT(in) :: step ! Current time step INTEGER, INTENT(in) :: dim_p ! PE-local dimension of state INTEGER, INTENT(in) :: dim_obs_p ! Dimension of observed state REAL, INTENT(in) :: state_p(dim_p) ! PE-local model state REAL, INTENT(out) :: m_state_p(dim_obs_p) ! PE-local observed state }}} The routine is called during the analysis step. It has to perform the operation of the observation operator acting on a state vector that is provided as `state_p`. The observed state has to be returned in `m_state_p`. For a model using domain decomposition, the operation is on the PE-local sub-domain of the model and has to provide the observed sub-state for the PE-local domain. Hint: * If the observation operator involves a global operation, e.g. some global integration, while using domain-decomposition one has to gather the information from the other model domains using MPI communication. === `U_init_obs` (init_obs_pdaf.F90) === This routine is used by all global filter algorithms (SEEK, SEIK, EnKF, ETKF, NETF, PF) and the 3D-Var methods. The interface for this routine is: {{{ SUBROUTINE init_obs(step, dim_obs_p, observation_p) INTEGER, INTENT(in) :: step ! Current time step INTEGER, INTENT(in) :: dim_obs_p ! PE-local dimension of obs. vector REAL, INTENT(out) :: observation_p(dim_obs_p) ! PE-local observation vector }}} The routine is called during the analysis step. It has to provide the vector of observations in `observation_p` for the current time step. For a model using domain decomposition, the vector of observations that exist on the model sub-domain for the calling process has to be initialized. === `U_prodRinvA` (prodrinva_pdaf.F90) === This routine is used by all filter algorithms that use the inverse of the observation error covariance matrix (SEEK, SEIK, and ETKF). The interface for this routine is: {{{ SUBROUTINE prodRinvA(step, dim_obs_p, dim_ens, obs_p, A_p, C_p) INTEGER, INTENT(in) :: step ! Current time step INTEGER, INTENT(in) :: dim_obs_p ! PE-local dimension of obs. vector INTEGER, INTENT(in) :: dim_ens ! Ensemble size REAL, INTENT(in) :: obs_p(dim_obs_p) ! PE-local vector of observations REAL, INTENT(in) :: A_p(dim_obs_p, dim_ens) ! Input matrix from analysis routine REAL, INTENT(out) :: C_p(dim_obs_p, dim_ens) ! Output matrix }}} The routine is called during the analysis step. In the algorithms the product of the inverse of the observation error covariance matrix with some matrix has to be computed. For the ETKF, this matrix holds the observed part of the ensemble perturbations. The matrix is provided as `A_p`. The product has to be given as `C_p`. For a model with domain decomposition, `A_p` contains the part of the matrix that resides on the model sub-domain of the calling process. The product has to be computed for this sub-domain, too. Hints: * The routine does not require that the product is implemented as a real matrix-matrix product. Rather, the product can be implemented in its most efficient form. For example, if the observation error covariance matrix is diagonal, only the multiplication of the diagonal with matrix `A_p` has to be implemented. * The observation vector `obs_p` is provided through the interface for cases where the observation error variance is relative to the actual value of the observations. * The interface has a difference for ESTKF/SEIK and ETKF: For ETKF the third argument is the ensemble size (`dim_ens`), while for ESTKF and SEIK it is the 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 and SEIK filters. (Practically, one can usually ignore this difference as the fourth argument of the interface can be named arbitrarily in the routine.) === `U_cvt_ens` (cvt_ens_pdaf.F90) === The interface for this routine is: {{{ SUBROUTINE cvt_ens_pdaf(iter, dim_p, dim_ens, dim_cv_ens_p, ens_p, cv_p, Vcv_p) INTEGER, INTENT(in) :: iter ! Iteration of optimization INTEGER, INTENT(in) :: dim_p ! PE-local observation dimension INTEGER, INTENT(in) :: dim_ens ! Ensemble size INTEGER, INTENT(in) :: dim_cv_ens_p ! Dimension of control vector REAL, INTENT(in) :: ens_p(dim_p, dim_ens) ! PE-local ensemble REAL, INTENT(in) :: cv_p(dim_cv_ens_p) ! PE-local control vector REAL, INTENT(inout) :: Vcv_p(dim_p) ! PE-local state increment }}} The routine is called during the analysis step during the iterative minimization of the cost function. It 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 hybrid 3D-Var, a part of this square root is usually expressed through the ensemble. If 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. === `U_cvt_adj_ens` (cvt_adj_ens_pdaf.F90) === The interface for this routine is: {{{ SUBROUTINE cvt_adj_ens_pdaf(iter, dim_p, dim_ens, dim_cv_ens_p, ens_p, Vcv_p, cv_p) INTEGER, INTENT(in) :: iter ! Iteration of optimization INTEGER, INTENT(in) :: dim_p ! PE-local observation dimension INTEGER, INTENT(in) :: dim_ens ! Ensemble size INTEGER, INTENT(in) :: dim_cv_ens_p ! PE-local dimension of control vector REAL, INTENT(in) :: ens_p(dim_p, dim_ens) ! PE-local ensemble REAL, INTENT(in) :: Vcv_p(dim_p) ! PE-local input vector REAL, INTENT(inout) :: cv_p(dim_cv_ens_p) ! PE-local result vector }}} The routine is called during the analysis step during the iterative minimization of the cost function. It 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'''. For the hybrid 3D-Var, a part of this square root is usually expressed through the ensemble. If 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. === `U_cvt` (cvt_pdaf.F90) === The interface for this routine is: {{{ SUBROUTINE cvt_pdaf(iter, dim_p, dim_cvec, cv_p, Vv_p) INTEGER, INTENT(in) :: iter ! Iteration of optimization INTEGER, INTENT(in) :: dim_p ! PE-local observation dimension INTEGER, INTENT(in) :: dim_cvec ! Dimension of control vector REAL, INTENT(in) :: cv_p(dim_cvec) ! PE-local control vector REAL, INTENT(inout) :: Vv_p(dim_p) ! PE-local result vector (state vector increment) }}} The routine is called during the analysis step during the iterative minimization of the cost function. It 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'''. If 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. === `U_cvt_adj` (cvt_adj_pdaf.F90) === The interface for this routine is: {{{ SUBROUTINE cvt_adj_pdaf(iter, dim_p, dim_cvec, Vv_p, cv_p) INTEGER, INTENT(in) :: iter ! Iteration of optimization INTEGER, INTENT(in) :: dim_p ! PE-local observation dimension INTEGER, INTENT(in) :: dim_cvec ! Dimension of control vector REAL, INTENT(in) :: Vv_p(dim_p) ! PE-local result vector (state vector increment) REAL, INTENT(inout) :: cv_p(dim_cvec) ! PE-local control vector }}} The routine is called during the analysis step during the iterative minimization of the cost function. It has to apply the adjoint control vector transformation to a state vector and return the control vector. Usually this transformation is the multiplication with transposed of the square-root of the background error covariance matrix '''B'''. If 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. === `U_obs_op_lin` (obs_op_lin_pdaf.F90) === This routine is used by all 3D-Var methods. The interface for this routine is: {{{ SUBROUTINE obs_op_lin(step, dim_p, dim_obs_p, state_p, m_state_p) INTEGER, INTENT(in) :: step ! Current time step INTEGER, INTENT(in) :: dim_p ! PE-local dimension of state INTEGER, INTENT(in) :: dim_obs_p ! Dimension of observed state REAL, INTENT(in) :: state_p(dim_p) ! PE-local model state REAL, INTENT(out) :: m_state_p(dim_obs_p) ! PE-local observed state }}} The routine is called during the analysis step. It has to perform the operation of the linearized observation operator acting on a state vector increment that is provided as `state_p`. The observed state has to be returned in `m_state_p`. For a model using domain decomposition, the operation is on the PE-local sub-domain of the model and has to provide the observed sub-state for the PE-local domain. Hint: * If the observation operator involves a global operation, e.g. some global integration, while using domain-decomposition one has to gather the information from the other model domains using MPI communication. === `U_obs_op_adj` (obs_op_adj_pdaf.F90) === This routine is used by all 3D-Var methods. The interface for this routine is: {{{ SUBROUTINE obs_op_adj(step, dim_p, dim_obs_p, state_p, m_state_p) INTEGER, INTENT(in) :: step ! Current time step INTEGER, INTENT(in) :: dim_p ! PE-local dimension of state INTEGER, INTENT(in) :: dim_obs_p ! Dimension of observed state REAL, INTENT(in) :: m_state_p(dim_obs_p) ! PE-local observed state REAL, INTENT(out) :: state_p(dim_p) ! PE-local model state }}} The routine is called during the analysis step. It has to perform the operation of the adjoint observation operator acting on a vector in observation space that is provided as m_state_p. The resulting state vector has to be returned in `m_state_p`. For a model using domain decomposition, the operation is on the PE-local sub-domain of the model and has to provide the observed sub-state for the PE-local domain. Hint: * If the observation operator involves a global operation, e.g. some global integration, while using domain-decomposition one has to gather the information from the other model domains using MPI communication. === `U_init_dim_obs_f` (init_dim_obs_f_pdaf.F90) === This routine is used by all filter algorithms with domain-localization (LSEIK, LETKF, LESTKF, LNETF) and is independent of the particular algorithm. The interface for this routine is: {{{ SUBROUTINE init_dim_obs_f(step, dim_obs_f) INTEGER, INTENT(in) :: step ! Current time step INTEGER, INTENT(out) :: dim_obs_f ! Dimension of full observation vector }}} The routine is called at the beginning of each analysis step, before the loop over all local analysis domains is entered. It has to initialize the size `dim_obs_f` of the full observation vector according to the current time step. For simplicity, `dim_obs_f` can be the size for the global model domain. Some hints: * It can be useful to not only determine the size of the observation vector at this point. One can also already gather information about the location of the observations, which can be used later, e.g. to implement the observation operator. In addition, one can already prepare an array that holds the full observation vector. This can be used later by `U_init_obs_l` to initialize a local vector of observations by selecting the relevant parts of the full observation vector. The required arrays can be defined in a module like `mod_assimilation`. * The routine is similar to `init_dim_obs` used in the global filters. However, if the global filter is used with a domain-decomposed model, it only initializes the size of the observation vector for the local model sub-domain. This is different for the local filters, as the local analysis also requires observational data from neighboring model sub-domains. Nonetheless, one can base on an implemented routine `init_dim_obs` to implement `init_dim_obs_f`. === `U_obs_op_f` (obs_op_f_pdaf.F90) === This routine is used by all filter algorithms with domain-localization (LSEIK, LETKF, LESTKF, LNETF) and is independent of the particular algorithm. The interface for this routine is: {{{ SUBROUTINE obs_op_f(step, dim_p, dim_obs_f, state_p, m_state_f) INTEGER, INTENT(in) :: step ! Current time step INTEGER, INTENT(in) :: dim_p ! PE-local dimension of state INTEGER, INTENT(in) :: dim_obs_f ! Dimension of the full observed state REAL, INTENT(in) :: state_p(dim_p) ! PE-local model state REAL, INTENT(out) :: m_state_f(dim_obs_f) ! Full observed state }}} The routine is called during the analysis step, before the loop over the local analysis domain is entered. It has to perform the operation of the observation operator acting on a state vector, which is provided as `state_p`. The observed state has to be returned in `m_state_f`. It is the observed state corresponding to the 'full' observation vector. Hint: * The routine is similar to `init_dim_obs` used for the global filters. However, with a domain-decomposed model `m_state_f` will need to contain parts of the state vector from neighboring model sub-domains. Thus, one needs to collect this information which resides in the memory of other processes. PDAF provides the routine [wiki:PDAF_gather_obs_f PDAF_gather_obs_f] for this task. The example implementation in `tutorial/classical/online_2D_parallelmodel` shows the use of `PDAF_gather_obs_f`. === `U_init_obs_f` (init_obs_f_pdaf.F90) === This routine is used by all filter algorithms with domain-localization (LSEIK, LESTKF, LETKF, LNETF) and is independent of the particular algorithm. The routine is only called if the globally adaptive forgetting factor is used (`type_forget=1` in the example implementation). For the local filters there is also the alternative to use locally adaptive forgetting factors (`type_forget=2` in the example implementation) The interface for this routine is: {{{ SUBROUTINE init_obs_f(step, dim_obs_f, observation_f) INTEGER, INTENT(in) :: step ! Current time step INTEGER, INTENT(in) :: dim_obs_f ! Dimension of full observation vector REAL, INTENT(out) :: observation_f(dim_obs_f) ! Full observation vector }}} The routine is called during the analysis step before the loop over the local analysis domains is entered. It has to provide the full vector of observations in `observation_f` for the current time step. The caller is the routine that computes an adaptive forgetting factor (PDAF_set_forget). Hints: * As for the other 'full' routines: While the global counterpart of this routine (`init_obs`) has to initialize the observation vector only for the local model sub-domain, the 'full' routine has to include observations that spatially belong to neighboring model sub-domains. As an easy choice one can simply initialize a vector of all globally available observations. * If the adaptive forgetting factor is not used, this routine only has to exist. However, no functionality is required. === `U_init_obs_l` (init_obs_l_pdaf.F90) === This routine is used by all filter algorithms with domain-localization (LSEIK, LESTKF, LETKF, LNETF) and is independent of the particular algorithm. The interface for this routine is: {{{ SUBROUTINE init_obs_l(domain_p, step, dim_obs_l, observation_l) INTEGER, INTENT(in) :: domain_p ! Current local analysis domain INTEGER, INTENT(in) :: step ! Current time step INTEGER, INTENT(in) :: dim_obs_l ! Local dimension of observation vector REAL, INTENT(out) :: observation_l(dim_obs_l) ! Local observation vector }}} The routine is called during the analysis step during the loop over the local analysis domain. It has to provide the vector of observations for the analysis in the local analysis domain with index `domain_p` in `observation_l` for the current time step. Hints: * For parallel efficiency, the LETKF algorithm is implemented in a way that first the full vectors are initialized. These are then restricted to the local analysis domain during the loop over all local analysis domains. Thus, if the full vector of observations has been initialized before `U_init_obs_l` is executed (e.g. by `U_init_dim_obs_f`), the operations performed in this routine will be to select the part of the full observation vector that is relevant for the current local analysis domain. * The routine `U_init_dim_obs_l` is executed before this routine. Thus, if that routine already prepares the information which elements of `observation_f` are needed for `observation_l`, this information can be used efficiently here. === `U_prodRinvA_l` (prodrinva_l_pdaf.F90) === This routine is used by the local filters (LSEIK, LESTKF, LETKF, LNETF). There is a slight difference between LSEIK and LETKF for this routine, which is described below. The interface for this routine is: {{{ SUBROUTINE prodRinvA_l(domain_p, step, dim_obs_l, dim_ens, obs_l, A_l, C_l) INTEGER, INTENT(in) :: domain_p ! Current local analysis domain INTEGER, INTENT(in) :: step ! Current time step INTEGER, INTENT(in) :: dim_obs_l ! Dimension of local observation vector INTEGER, INTENT(in) :: dim_ens ! Ensemble size REAL, INTENT(in) :: obs_l(dim_obs_l) ! Local vector of observations REAL, INTENT(inout) :: A_l(dim_obs_l, dim_ens) ! Input matrix from analysis routine REAL, INTENT(out) :: C_l(dim_obs_l, dim_ens) ! Output matrix }}} The 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 ESTKF 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`. This 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. Hints: * 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. * 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. * 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. * The interface has a difference for ESTKF/SEIK and ETKF: For ETKF the third argument is the ensemble size (`dim_ens`), while for ESTKF/SEIK 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 ETKF, while it is `rank` for the ESTKF/SEIK filter. (Practically, one can usually ignore this difference as the fourth argument of the interface can be named arbitrarily in the routine.) === `U_init_n_domains` (init_n_domains_pdaf.F90) === This routine is used by all filter algorithms with domain-localization (LSEIK, LESTKF, LETKF, LNETF) and is independent of the particular algorithm. The interface for this routine is: {{{ SUBROUTINE init_n_domains(step, n_domains_p) INTEGER, INTENT(in) :: step ! Current time step INTEGER, INTENT(out) :: n_domains_p ! Number of analysis domains for local model sub-domain }}} The routine 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. Hints: * 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. === `U_init_dim_l` (init_dim_l_pdaf.F90) === This routine is used by all filter algorithms with domain-localization (LSEIK, LESTKF, LETKF, LNETF) and is independent of the particular algorithm. The interface for this routine is: {{{ SUBROUTINE init_dim_l(step, domain_p, dim_l) INTEGER, INTENT(in) :: step ! Current time step INTEGER, INTENT(in) :: domain_p ! Current local analysis domain INTEGER, INTENT(out) :: dim_l ! Local state dimension }}} The routine is called during the loop over the local analysis domains in the analysis step. It has to provide in `dim_l` the dimension of the state vector for the local analysis domain with index `domain_p`. Hints: * If a local analysis domain is a single vertical column of the model grid, the size of the state in the local analysis domain will be just the number of vertical grid points at this location. === `U_init_dim_obs_l` (init_dim_obs_l_pdaf.F90) === This routine is used by all filter algorithms with domain-localization (LSEIK, LESTKF, LETKF, LNETF) and is independent of the particular algorithm. The interface for this routine is: {{{ SUBROUTINE init_dim_obs_l(domain_p, step, dim_obs_f, dim_obs_l) INTEGER, INTENT(in) :: domain_p ! Current local analysis domain INTEGER, INTENT(in) :: step ! Current time step INTEGER, INTENT(in) :: dim_obs_f ! Full dimension of observation vector INTEGER, INTENT(out) :: dim_obs_l ! Local dimension of observation vector }}} The routine is called during the loop over the local analysis domains in the analysis step. It has to initialize in `dim_obs_l` the size of the observation vector used for the local analysis domain with index `domain_p`. Some hints: * Usually, the observations to be considered for a local analysis are those which reside within some distance from the local analysis domain. Thus, if the local analysis domain is a single vertical column of the model grid and if the model grid is a regular ij-grid, then one could use some range of i/j indices to select the observations and determine the local number of them. More generally, one can compute the physical distance of an observation from the local analysis domain and decide on this basis, if the observation has to be considered. * In the loop over the local analysis domains, the routine is always called before `U_init_obs_l` is executed. Thus, as `U_init_dim_obs_local` has to check which observations should be used for the local analysis domain, one can already initialize an integer array that stores the index of observations to be considered. This index should be the position of the observation in the array `observation_f`. With this, the initialization of the local observation vector in `U_init_obs_l` can be sped up. * For PDAF, we could not join the routines `U_init_dim_obs_l` and `U_init_obs_l`, because the array for the local observations is allocated internally to PDAF after its size has been determined in `U_init_dim_obs_l`. === `U_g2l_state` (g2l_state_pdaf.F90) === This routine is used by all filter algorithms with domain-localization (LSEIK, LESTKF, LETKF, LNETF) and is independent of the particular algorithm. The interface for this routine is: {{{ SUBROUTINE global2local_state(step, domain_p, dim_p, state_p, dim_l, state_l) INTEGER, INTENT(in) :: step ! Current time step INTEGER, INTENT(in) :: domain_p ! Current local analysis domain INTEGER, INTENT(in) :: dim_p ! State dimension for model sub-domain INTEGER, INTENT(in) :: dim_l ! Local state dimension REAL, INTENT(in) :: state_p(dim_p) ! State vector for model sub-domain REAL, INTENT(out) :: state_l(dim_l) ! State vector on local analysis domain }}} The routine is called during the loop over the local analysis domains in the analysis step. It has to provide the local state vector `state_l` that corresponds to the local analysis domain with index `domain_p`. Provided to the routine is the state vector `state_p`. With a domain decomposed model, this is the state for the local model sub-domain. Hints: * In the simple case that a local analysis domain is a single vertical column of the model grid, the operation in this routine would be to take out of `state_p` the data for the vertical column indexed by `domain_p`. === `U_l2g_state` (l2g_state_pdaf.F90) === This routine is used by all filter algorithms with domain-localization (LSEIK, LESTKF, LETKF, LNETF) and is independent of the particular algorithm. The interface for this routine is: {{{ SUBROUTINE l2g_state(step, domain_p, dim_l, state_l, dim_p, state_p) INTEGER, INTENT(in) :: step ! Current time step INTEGER, INTENT(in) :: domain_p ! Current local analysis domain INTEGER, INTENT(in) :: dim_p ! State dimension for model sub-domain INTEGER, INTENT(in) :: dim_l ! Local state dimension REAL, INTENT(in) :: state_p(dim_p) ! State vector for model sub-domain REAL, INTENT(out) :: state_l(dim_l) ! State vector on local analysis domain }}} The routine is called during the loop over the local analysis domains in the analysis step. It has to initialize the part of the global state vector `state_p` that corresponds to the local analysis domain with index `domain_p`. Provided to the routine is the state vector `state_l` for the local analysis domain. Hints: * In the simple case that a local analysis domain is a single vertical column of the model grid, the operation in this routine would be to write into `state_p` the data for the vertical column indexed by `domain_p`. === `U_g2l_obs` (g2l_obs_pdaf.F90) === This routine is used by all filter algorithms with domain-localization (LSEIK, LESTKF, LETKF, LNETF) and is independent of the particular algorithm. The interface for this routine is: {{{ SUBROUTINE g2l_obs(domain_p, step, dim_obs_f, dim_obs_l, mstate_f, mstate_l) INTEGER, INTENT(in) :: domain_p ! Current local analysis domain INTEGER, INTENT(in) :: step ! Current time step INTEGER, INTENT(in) :: dim_obs_f ! Dimension of full observation vector for model sub-domain INTEGER, INTENT(in) :: dim_obs_l ! Dimension of observation vector for local analysis domain REAL, INTENT(in) :: mstate_f(dim_obs_f) ! Full observation vector for model sub-domain REAL, INTENT(out) :: mstate_l(dim_obs_l) ! Observation vector for local analysis domain }}} The routine is called during the loop over the local analysis domains in the analysis step. It has to provide a local observation vector `mstate_l` for the observation domain that corresponds to the local analysis domain with index `domain_p`. Provided to the routine is the full observation vector `mstate_f` from which the local part has to be extracted. Hints: * The vector `mstate_f` that is provided to the routine is one of the observed state vectors that are produced by `U_obs_op_full`. * Some operations performed here are analogous to those required to initialize a local vector of observations in `U_init_obs_l`. If that routine reads first a full vector of observations (e.g. in `U_init_dim_obs_f`), this vector has to be restricted to the relevant observations for the current local analysis domain. For this operation, one can for example initialize an index array when `U_init_dim_obs_l` is executed. (Which happens before `U_global2local_obs`) === `U_init_obsvar` (init_obsvar_pdaf.F90) === This routine is used by the global filter algorithms SEIK, ESTKF, and ETKF as well as the local filters LESTKF, LSEIK, LETKF. The routine is only called if the adaptive forgetting factor is used (`type_forget=1` in the example implementation). The difference in this routine between global and local filters is that the global filters use 'global' while the local filters use 'full' quantities. The interface for this routine is: {{{ SUBROUTINE init_obsvar(step, dim_obs_f, obs_f, meanvar_f) INTEGER, INTENT(in) :: step ! Current time step INTEGER, INTENT(in) :: dim_obs_f ! Full dimension of observation vector REAL, INTENT(in) :: obs_f(dim_obs_f) ! Full observation vector REAL, INTENT(out) :: meanvar_f ! Mean observation error variance }}} The routine is called in the local filters before the loop over all local analysis domains is entered. The call is by the routine that computes an adaptive forgetting factor (`PDAF_set_forget`). The routine has to initialize an average full observation error variance, which should be consistent with the observation vector initialized in `U_init_ob_f`. Hints: * For a model with domain-decomposition one might use the mean variance for the model sub-domain of the calling process. Alternatively one can compute a mean variance for the full model domain using MPI communication (e.g. the function `MPI_allreduce`). * The observation vector `obs_p` is provided to the routine for the case that the observation error variance is relative to the value of the observations. * If the adaptive forgetting factor is not used, this routine has only to exist for the compilation, but it does not need functionality. === `U_init_obsvar_l` (init_obsvar_l_pdaf.F90) === This routine is used by all filter algorithms with domain-localization (LSEIK, LESTKF, LETKF, LNETF) and is independent of the particular algorithm. The routine is only called if the local adaptive forgetting factor is used (`type_forget=2` in the example implementation). The interface for this routine is: {{{ SUBROUTINE init_obsvar_l(domain_p, step, dim_obs_l, obs_l, meanvar_l) INTEGER, INTENT(in) :: domain_p ! Current local analysis domain INTEGER, INTENT(in) :: step ! Current time step INTEGER, INTENT(in) :: dim_obs_l ! Local dimension of observation vector REAL, INTENT(in) :: obs_l(dim_obs_l) ! Local observation vector REAL, INTENT(out) :: meanvar_l ! Mean local observation error variance }}} The routine is called in the local filters during the loop over all local analysis domains by the routine that computes a local adaptive forgetting factor (`PDAF_set_forget_local`). The routine has to initialize a local mean observation error variance for all observations used for the analysis in the current local analysis domain. Hints: * If the local adaptive forgetting factor is not used, this routine has only to exist for the compilation, but it does not need functionality. === `U_prepoststep` (prepoststep_ens_pdaf.F90) === The routine has already been described for modifying the model for the ensemble integration and for inserting the analysis step. See the page on [InsertAnalysisStep#U_prepoststepprepoststep_ens_pdaf.F90 inserting the analysis step] for the description of this routine. === `U_next_observation` (next_observation_pdaf.F90) === This routine is independent of the filter algorithm used. See the page on [InsertAnalysisStep#U_next_observationnext_observation_pdaf.F90 inserting the analysis step] for the description of this routine. == Execution order of user-supplied routines == The user-supplied routines are essentially executed in the order they are listed in the interfaces to the routines. An except is the order of the routines for the local ESTKF. The order can be important as some routines can perform preparatory work for later routines. For example, `U_init_dim_obs` prepares an index array that provides the information for executing the observation operator in `U_obs_op`. Before the analysis step is called the following routine is executed: 1. [#U_collect_statecollect_state_pdaf.F90 U_collect_state] The analysis step is executed when the ensemble integration of the forecast is completed. During the analysis step the following routines are executed in the given order: 1. [#U_prepoststepprepoststep_ens_pdaf.F90 U_prepoststep] (Call to act on the forecast ensemble, called with negative value of the time step) 1. [#U_init_dim_obsinit_dim_obs_pdaf.F90 U_init_dim_obs] 1. [#U_obs_opobs_op_pdaf.F90 U_obs_op] 1. [#U_init_obsinit_obs_pdaf.F90 U_init_obs] Inside the analysis step the interative optimization is computed. This involves the repeated call of the routines: 1. [#U_cvtcvt_pdaf.F90 U_cvt] 1. [#U_cvt_enscvt_ens_pdaf.F90 U_cvt_ens] 1. [#U_obs_op_linobs_op_lin_pdaf.F90 U_obs_op_lin] 1. [#U_prodRinvAprodrinva_pdaf.F90 U_prodRinvA] 1. [#U_obs_op_adjobs_op_adj_pdaf.F90 U_obs_op_adj] 1. [#U_cvt_adjcvt_adj_pdaf.F90 U_cvt_adj] 1. [#U_cvt_adj_enscvt_adj_ens_pdaf.F90 U_cvt_adj_ens] After the iterative optimization the following routines are executes to complte the analysis step: 1. [#U_cvt_enscvt_pdaf.F90 U_cvt] (Call to the parameterized part of the control vector transform to compute the final state vector increment) 1. [#U_cvt_enscvt_ens_pdaf.F90 U_cvt_ens] (Call to the eensemble-part of the control vector transform to compute the final state vector increment) 1. [#U_prepoststepprepoststep_ens_pdaf.F90 U_prepoststep] (Call to act on the analysis ensemble, called with (positive) value of the time step) The 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:ImplementAnalysislestkf page on implementing the analysis step of the LESTKF] and on the [wiki:ImplementAnalysisestkf page on implementing the analysis step of the ESTKF]. In case of the routines `PDAF_assimilate_*`, the following routines are executed after the analysis step: 1. [#U_distribute_statedistribute_state_pdaf.F90 U_distribute_state] 1. [#U_next_observationnext_observation_pdaf.F90 U_next_observation]