| | 1 | = PDAF-OMI Debugging Information = |
| | 2 | |
| | 3 | [[PageOutline(2-3,Contents of this page)]] |
| | 4 | |
| | 5 | == Overview == |
| | 6 | |
| | 7 | When implementing an observation with PDAF, or when performing the very first implementation of PDAF with a new model, it is useful to check whether the inputs to the PDAF-routines are correctly used. For this porpose, PDAF-OMI provides a debugging functionality. It allows you to activate debugging output e.g. for a single local analysis domain on a single process of a complex application of a local filter like LEKSTF. |
| | 8 | |
| | 9 | |
| | 10 | == Activating Debugging Output == |
| | 11 | |
| | 12 | Debugging output is activated by the routine `PDAFomi_set_debug_flag`. In particular this call can be inserted in any routines contained in `callback_obs_pdafomi.F90`. |
| | 13 | |
| | 14 | For example to activate debugging in `init_dim_obs_l_pdafomi` for the local analysis domain `domain_p=10` and filter process 0, one uses |
| | 15 | {{{ |
| | 16 | SUBROUTINE init_dim_obs_l_pdafomi(domain_p, step, dim_obs, dim_obs_l) |
| | 17 | |
| | 18 | USE PDAFomi, ONLY: PDAFomi_set_debug_flag |
| | 19 | USE mod_parallel_pdaf, ONLY: mype_filter |
| | 20 | |
| | 21 | ... |
| | 22 | |
| | 23 | IF (domain_p==10 .AND. mype_filter==0) THEN |
| | 24 | CALL PDAFomi_set_debug_flag(domain_p) |
| | 25 | ELSE |
| | 26 | CALL PDAFomi_set_debug_flag(0) |
| | 27 | ENDIF |
| | 28 | |
| | 29 | ... |
| | 30 | |
| | 31 | }}} |
| | 32 | |
| | 33 | == Understanding the Debugging Output == |
| | 34 | |
| | 35 | The debugging output mainly writes information about the different variables contained in the full data type `obs_f` allocated as `thisobs` and the local type `obs_l` allocate as `thisobs_l`. For reference we list the full declaration of these types. When reading the debugging output one can check for the meaning of the variables. |
| | 36 | |
| | 37 | {{{ |
| | 38 | TYPE obs_f |
| | 39 | ! ---- Mandatory variables to be set in INIT_DIM_OBS ---- |
| | 40 | INTEGER :: doassim=0 !< Whether to assimilate this observation type |
| | 41 | INTEGER :: disttype !< Type of distance computation to use for localization |
| | 42 | INTEGER :: ncoord !< Number of coordinates use for distance computation |
| | 43 | INTEGER, ALLOCATABLE :: id_obs_p(:,:) !< Indices of process-local observed field in state vector |
| | 44 | |
| | 45 | ! ---- Optional variables - they can be set in INIT_DIM_OBS ---- |
| | 46 | REAL, ALLOCATABLE :: icoeff_p(:,:) !< Interpolation coefficients for obs. operator (optional) |
| | 47 | REAL, ALLOCATABLE :: domainsize(:) !< Size of domain for periodicity (<=0 for no periodicity) (optional) |
| | 48 | |
| | 49 | ! ---- Variables with predefined values - they can be changed in INIT_DIM_OBS ---- |
| | 50 | INTEGER :: obs_err_type=0 !< Type of observation error: (0) Gauss, (1) Laplace |
| | 51 | INTEGER :: use_global_obs=1 !< Whether to use (1) global full obs. |
| | 52 | !< or (0) obs. restricted to those relevant for a process domain |
| | 53 | |
| | 54 | ! ---- The following variables are set in the routine PDAFomi_gather_obs --- |
| | 55 | INTEGER :: dim_obs_p !< number of PE-local observations |
| | 56 | INTEGER :: dim_obs_f !< number of full observations |
| | 57 | INTEGER :: dim_obs_g !< global number of observations |
| | 58 | INTEGER :: off_obs_f !< Offset of this observation in overall full obs. vector |
| | 59 | INTEGER :: off_obs_g !< Offset of this observation in overall global obs. vector |
| | 60 | INTEGER :: obsid !< Index of observation over all assimilated observations |
| | 61 | REAL, ALLOCATABLE :: obs_f(:) !< Full observed field |
| | 62 | REAL, ALLOCATABLE :: ocoord_f(:,:) !< Coordinates of full observation vector |
| | 63 | REAL, ALLOCATABLE :: ivar_obs_f(:) !< Inverse variance of full observations |
| | 64 | INTEGER, ALLOCATABLE :: id_obs_f_lim(:) !< Indices of domain-relevant full obs. in global vector of obs. |
| | 65 | END TYPE obs_f |
| | 66 | }}} |
| | 67 | |
| | 68 | {{{ |
| | 69 | TYPE obs_l |
| | 70 | INTEGER :: dim_obs_l !< number of local observations |
| | 71 | INTEGER :: off_obs_l !< Offset of this observation in overall local obs. vector |
| | 72 | INTEGER, ALLOCATABLE :: id_obs_l(:) !< Indices of local observations in full obs. vector |
| | 73 | REAL, ALLOCATABLE :: distance_l(:) !< Distances of local observations |
| | 74 | REAL, ALLOCATABLE :: ivar_obs_l(:) !< Inverse variance of local observations |
| | 75 | INTEGER :: locweight !< Specify localization function |
| | 76 | REAL :: lradius !< localization radius |
| | 77 | REAL :: sradius !< support radius for localization function |
| | 78 | END TYPE obs_l |
| | 79 | }}} |
| | 80 | |
| | 81 | === Example output === |
| | 82 | |
| | 83 | For illustration we insert the above call to PDAFomi_set_debug flag into `init_dim_obs_l_pdafomi` in `tutorial/online_2D_serialmodel_omi/`. then we compile and execute the tutorial program with: |
| | 84 | {{{ |
| | 85 | mpirun -np 4 ./model_pdaf -dim_ens 4 -filtertype 7 -local_range 5 -locweight 2 -assim_B .true. -assim_A .false. |
| | 86 | }}} |
| | 87 | With these settings only observation type B is activated, which are just 3 values in the model domain |
| | 88 | |
| | 89 | The first lines of the debugging output looks like this: |
| | 90 | {{{ |
| | 91 | ++ OMI-debug set_debug_flag: mype_filter 0 activate 10 |
| | 92 | ++ OMI-debug: 10 PDAFomi_init_dim_obs_l -- START |
| | 93 | ++ OMI-debug: 10 PDAFomi_init_dim_obs_l -- count local observations |
| | 94 | ++ OMI-debug init_dim_obs_l: 10 Re-init dim_obs_l=0 |
| | 95 | ++ OMI-debug init_dim_obs_l: 10 coords_l 1.0000000000000000 10.000000000000000 |
| | 96 | ++ OMI-debug cnt_dim_obs_l: 10 thisobs%ncoord 2 |
| | 97 | ++ OMI-debug cnt_dim_obs_l: 10 thisobs_l%lradius 5.0000000000000000 |
| | 98 | ++ OMI-debug cnt_dim_obs_l: 10 Check for observations within radius |
| | 99 | ++ OMI-debug comp_dist2: 10 compute Cartesian distance |
| | 100 | ++ OMI-debug cnt_dim_obs_l: 10 valid observation with coordinates 5.0000000000000000 8.0000000000000000 |
| | 101 | ++ OMI-debug: 10 PDAFomi_init_dim_obs_l -- initialize local observation arrays |
| | 102 | ++ OMI-debug comp_dist2: 10 compute Cartesian distance |
| | 103 | ++ OMI-debug init_dim_obs_l: 10 thisobs_l%dim_obs_l 1 |
| | 104 | ++ OMI-debug init_dim_obs_l: 10 thisobs_l%id_obs_l 1 |
| | 105 | ++ OMI-debug init_dim_obs_l: 10 thisobs_l%distance_l 4.4721359549995796 |
| | 106 | ++ OMI-debug: 10 PDAFomi_init_dim_obs_l -- END |
| | 107 | }}} |
| | 108 | |
| | 109 | The first line |
| | 110 | {{{ |
| | 111 | ++ OMI-debug set_debug_flag: mype_filter 0 activate 10 |
| | 112 | }}} |
| | 113 | is from PDAFomi_set_debug_flag showing that debugging is activates with value 10 (which is values fo `domain_p` specified in the call) |
| | 114 | |
| | 115 | The next lines are |
| | 116 | {{{ |
| | 117 | ++ OMI-debug: 10 PDAFomi_init_dim_obs_l -- START |
| | 118 | ++ OMI-debug: 10 PDAFomi_init_dim_obs_l -- count local observations |
| | 119 | ++ OMI-debug init_dim_obs_l: 10 Re-init dim_obs_l=0 |
| | 120 | }}} |
| | 121 | show that debugging output for PDAFomi_init_dim_obs_l is shown. Most routines show such a 'START' line. The second line shows that a segment of the orutine started, the counting of local observations. The third line states that dim_obs_l=0 is set. This line, as many others shows the name of the subroutine in short form without `PDAFomi` at the beginning of the line. |
| | 122 | |
| | 123 | The following lines show variable values |
| | 124 | {{{ |
| | 125 | ++ OMI-debug init_dim_obs_l: 10 coords_l 1.0000000000000000 10.000000000000000 |
| | 126 | ++ OMI-debug cnt_dim_obs_l: 10 thisobs%ncoord 2 |
| | 127 | ++ OMI-debug cnt_dim_obs_l: 10 thisobs_l%lradius 5.0000000000000000 |
| | 128 | }}} |
| | 129 | First we see that the coordinates `coords_l` of the grid point corresponding to domain_p=10 are (1.0, 10.0). Further we have two dimensions (thisobs%ncoord=2) and the localization radius is set to 5.0. |
| | 130 | In the following lines |
| | 131 | {{{ |
| | 132 | ++ OMI-debug cnt_dim_obs_l: 10 Check for observations within radius |
| | 133 | ++ OMI-debug comp_dist2: 10 compute Cartesian distance |
| | 134 | ++ OMI-debug cnt_dim_obs_l: 10 valid observation with coordinates 5.0000000000000000 8.0000000000000000 |
| | 135 | }}} |
| | 136 | the it is checked which observations lie within the distance thisobs_l%lradius=5.0 from coords_l=(1.0, 5.0). One observation with coordinates (5.0, 8.0) is found. The followign lines |
| | 137 | {{{ |
| | 138 | ++ OMI-debug: 10 PDAFomi_init_dim_obs_l -- initialize local observation arrays |
| | 139 | ++ OMI-debug comp_dist2: 10 compute Cartesian distance |
| | 140 | ++ OMI-debug init_dim_obs_l: 10 thisobs_l%dim_obs_l 1 |
| | 141 | ++ OMI-debug init_dim_obs_l: 10 thisobs_l%id_obs_l 1 |
| | 142 | ++ OMI-debug init_dim_obs_l: 10 thisobs_l%distance_l 4.4721359549995796 |
| | 143 | }}} |
| | 144 | show that observation arrays are initialized. A Cartesian distance is computerd. `thisobs_l%dim_obs_l` confirms that one local observation was found and it's the first element of the full observation vector (thisobs_l%id_obs_l=1) and the distance is thisobs_l%distance_l=4.4721359549995796. |
| | 145 | |
| | 146 | One could now compare this with the input information. Are the values of coords_l correct? Are the coordinates of the first observation (5.0, 8.0)? |
| | 147 | |
| | 148 | Later in the debugging output we find e.g. |
| | 149 | {{{ |
| | 150 | ++ OMI-debug: 10 PDAFomi_init_obs_l -- Get local vector of observations |
| | 151 | ++ OMI-debug g2l_obs: 10 thisobs%id_obs_l 1 |
| | 152 | ++ OMI-debug g2l_obs: 10 obs_l -0.98253999999999997 |
| | 153 | ++ OMI-debug: 10 PDAFomi_init_obs_l -- Get local vector of inverse obs. variances |
| | 154 | ++ OMI-debug g2l_obs: 10 thisobs%id_obs_l 1 |
| | 155 | ++ OMI-debug g2l_obs: 10 obs_l 4.0000000000000000 |
| | 156 | }}} |
| | 157 | These lines are from the internal routine `PDAFomi_init_obs_l`, which initializes the local vector of observations. (This is functionality also existing in the 'traditional' form of implementing the observation. The documention page about [wiki:init_obs_l_pdaf init_obs_l_pdaf] explains what is done here.) Important is that not only the local observation vector (a single value) is initialized, but also the inverse observation variance (4.0 in this example). |
| | 158 | |
| | 159 | Similarly |
| | 160 | {{{ |
| | 161 | ++ OMI-debug: 10 PDAFomi_g2l_obs -- START Get local observed ensemble member 1 |
| | 162 | ++ OMI-debug g2l_obs: 10 thisobs%id_obs_l 1 |
| | 163 | ++ OMI-debug g2l_obs: 10 obs_l 0.12725110911115362 |
| | 164 | ++ OMI-debug: 10 PDAFomi_g2l_obs -- END |
| | 165 | }}} |
| | 166 | shows the localization of ensemble member 1. This is performed in the interval routine `PDAFomi_g2l_obs` (see the page about [wiki:g2l_obs_pdaf g2l_obs_pdaf]). There is analogous output for all of the 4 ensemble states. |
| | 167 | |
| | 168 | |