[2128] | 1 | MODULE obs_read_vel |
---|
| 2 | !!====================================================================== |
---|
| 3 | !! *** MODULE obs_read_vel *** |
---|
| 4 | !! Observation diagnostics: Read the velocity profile observations |
---|
| 5 | !!====================================================================== |
---|
| 6 | |
---|
| 7 | !!---------------------------------------------------------------------- |
---|
| 8 | !! obs_rea_vel_dri : Driver for reading profile obs |
---|
| 9 | !!---------------------------------------------------------------------- |
---|
| 10 | |
---|
| 11 | !! * Modules used |
---|
| 12 | USE par_kind ! Precision variables |
---|
| 13 | USE par_oce ! Ocean parameters |
---|
| 14 | USE in_out_manager ! I/O manager |
---|
| 15 | USE dom_oce ! Ocean space and time domain variables |
---|
| 16 | USE obs_mpp ! MPP support routines for observation diagnostics |
---|
| 17 | USE julian ! Julian date routines |
---|
| 18 | USE obs_utils ! Observation operator utility functions |
---|
| 19 | USE obs_prep ! Prepare observation arrays |
---|
| 20 | USE obs_grid ! Grid search |
---|
| 21 | USE obs_sort ! Sorting observation arrays |
---|
| 22 | USE obs_profiles_def ! Profile definitions |
---|
| 23 | USE obs_conv ! Various conversion routines |
---|
| 24 | USE obs_types ! Observation type definitions |
---|
| 25 | USE netcdf ! NetCDF library |
---|
| 26 | USE obs_oper ! Observation operators |
---|
| 27 | USE obs_vel_io ! Velocity profile files I/O (non-FB files) |
---|
[2715] | 28 | USE lib_mpp ! For ctl_warn/stop |
---|
[2128] | 29 | |
---|
| 30 | IMPLICIT NONE |
---|
| 31 | |
---|
| 32 | !! * Routine accessibility |
---|
| 33 | PRIVATE |
---|
| 34 | |
---|
| 35 | PUBLIC obs_rea_vel_dri ! Read the profile observations |
---|
| 36 | |
---|
[2287] | 37 | !!---------------------------------------------------------------------- |
---|
| 38 | !! NEMO/OPA 3.3 , NEMO Consortium (2010) |
---|
| 39 | !! $Id$ |
---|
| 40 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
| 41 | !!---------------------------------------------------------------------- |
---|
| 42 | |
---|
[2128] | 43 | CONTAINS |
---|
| 44 | |
---|
| 45 | SUBROUTINE obs_rea_vel_dri( kformat, & |
---|
| 46 | & profdata, knumfiles, cfilenames, & |
---|
| 47 | & kvars, kextr, kstp, ddobsini, ddobsend, & |
---|
| 48 | & ldignmis, ldavtimset, ldmod ) |
---|
| 49 | !!--------------------------------------------------------------------- |
---|
| 50 | !! |
---|
| 51 | !! *** ROUTINE obs_rea_pro_dri *** |
---|
| 52 | !! |
---|
| 53 | !! ** Purpose : Read from file the profile observations |
---|
| 54 | !! |
---|
| 55 | !! ** Method : Depending on kformat either ENACT, CORIOLIS or |
---|
| 56 | !! feedback data files are read |
---|
| 57 | !! |
---|
| 58 | !! ** Action : |
---|
| 59 | !! |
---|
| 60 | !! References : |
---|
| 61 | !! |
---|
| 62 | !! History : |
---|
| 63 | !! ! : 2009-01 (K. Mogensen) : New merged version of old routines |
---|
| 64 | !!---------------------------------------------------------------------- |
---|
| 65 | !! * Modules used |
---|
| 66 | |
---|
| 67 | !! * Arguments |
---|
| 68 | INTEGER :: kformat ! Format of input data |
---|
| 69 | ! ! 1: ENACT |
---|
| 70 | ! ! 2: Coriolis |
---|
| 71 | TYPE(obs_prof), INTENT(OUT) :: profdata ! Profile data to be read |
---|
| 72 | INTEGER, INTENT(IN) :: knumfiles ! Number of files to read in |
---|
| 73 | CHARACTER(LEN=128), INTENT(IN) :: cfilenames(knumfiles) ! File names to read in |
---|
| 74 | INTEGER, INTENT(IN) :: kvars ! Number of variables in profdata |
---|
| 75 | INTEGER, INTENT(IN) :: kextr ! Number of extra fields for each var in profdata |
---|
| 76 | INTEGER, INTENT(IN) :: kstp ! Ocean time-step index |
---|
| 77 | LOGICAL, INTENT(IN) :: ldignmis ! Ignore missing files |
---|
| 78 | LOGICAL, INTENT(IN) :: ldavtimset ! Set time to be equal to the end of the day |
---|
| 79 | LOGICAL, INTENT(IN) :: ldmod ! Initialize model from input data |
---|
| 80 | REAL(KIND=dp), INTENT(IN) :: ddobsini ! Obs. ini time in YYYYMMDD.HHMMSS |
---|
| 81 | REAL(KIND=dp), INTENT(IN) :: ddobsend ! Obs. end time in YYYYMMDD.HHMMSS |
---|
| 82 | |
---|
| 83 | !! * Local declarations |
---|
| 84 | CHARACTER(LEN=15), PARAMETER :: cpname='obs_rea_vel_dri' |
---|
| 85 | INTEGER :: jvar |
---|
| 86 | INTEGER :: ji |
---|
| 87 | INTEGER :: jj |
---|
| 88 | INTEGER :: jk |
---|
| 89 | INTEGER :: ij |
---|
| 90 | INTEGER :: iflag |
---|
| 91 | INTEGER :: inobf |
---|
| 92 | INTEGER :: i_file_id |
---|
| 93 | INTEGER :: inowin |
---|
| 94 | INTEGER :: iyea |
---|
| 95 | INTEGER :: imon |
---|
| 96 | INTEGER :: iday |
---|
| 97 | INTEGER :: ihou |
---|
| 98 | INTEGER :: imin |
---|
| 99 | INTEGER :: isec |
---|
| 100 | INTEGER, DIMENSION(knumfiles) :: & |
---|
| 101 | & irefdate |
---|
| 102 | INTEGER, DIMENSION(ntyp1770+1) :: & |
---|
| 103 | & itypuv, & |
---|
| 104 | & itypuvmpp |
---|
| 105 | INTEGER :: iuv3dtmpp |
---|
| 106 | INTEGER, DIMENSION(:), ALLOCATABLE :: & |
---|
| 107 | & iobsiu, & |
---|
| 108 | & iobsju, & |
---|
| 109 | & iprocu, & |
---|
| 110 | & iobsiv, & |
---|
| 111 | & iobsjv, & |
---|
| 112 | & iprocv, & |
---|
| 113 | & iindx, & |
---|
| 114 | & ifileidx, & |
---|
| 115 | & iprofidx |
---|
| 116 | INTEGER :: itype |
---|
| 117 | REAL(wp), DIMENSION(:), ALLOCATABLE :: & |
---|
| 118 | & zphi, & |
---|
| 119 | & zlam |
---|
| 120 | REAL(dp), DIMENSION(:), ALLOCATABLE :: & |
---|
| 121 | & zdat |
---|
| 122 | LOGICAL :: & |
---|
| 123 | & llvalprof |
---|
| 124 | TYPE(obfbdata), POINTER, DIMENSION(:) :: & |
---|
| 125 | & inpfiles |
---|
| 126 | REAL(dp), DIMENSION(knumfiles) :: & |
---|
| 127 | & djulini, & |
---|
| 128 | & djulend |
---|
| 129 | INTEGER :: iprof |
---|
| 130 | INTEGER :: iproftot |
---|
| 131 | INTEGER :: iuv3dt |
---|
| 132 | INTEGER, DIMENSION(kvars) :: iv3dt |
---|
| 133 | CHARACTER(len=8) :: cl_refdate |
---|
| 134 | |
---|
| 135 | ! Local initialization |
---|
| 136 | iprof = 0 |
---|
| 137 | iuv3dt = 0 |
---|
| 138 | |
---|
| 139 | !----------------------------------------------------------------------- |
---|
| 140 | ! Check data the model part is just with feedback data files |
---|
| 141 | !----------------------------------------------------------------------- |
---|
| 142 | IF ( ldmod .AND. ( kformat /= 0 ) ) THEN |
---|
| 143 | CALL ctl_stop( 'Model can only be read from feedback data' ) |
---|
| 144 | RETURN |
---|
| 145 | ENDIF |
---|
| 146 | |
---|
| 147 | !----------------------------------------------------------------------- |
---|
| 148 | ! Count the number of files needed and allocate the obfbdata type |
---|
| 149 | !----------------------------------------------------------------------- |
---|
| 150 | |
---|
| 151 | inobf = knumfiles |
---|
| 152 | |
---|
| 153 | ALLOCATE( inpfiles(inobf) ) |
---|
| 154 | |
---|
| 155 | prof_files : DO jj = 1, inobf |
---|
| 156 | |
---|
| 157 | !--------------------------------------------------------------------- |
---|
| 158 | ! Prints |
---|
| 159 | !--------------------------------------------------------------------- |
---|
| 160 | IF(lwp) THEN |
---|
| 161 | WRITE(numout,*) |
---|
| 162 | WRITE(numout,*) ' obs_rea_vel_dri : Reading from file = ', & |
---|
| 163 | & TRIM( TRIM( cfilenames(jj) ) ) |
---|
| 164 | WRITE(numout,*) ' ~~~~~~~~~~~~~~~' |
---|
| 165 | WRITE(numout,*) |
---|
| 166 | ENDIF |
---|
| 167 | |
---|
| 168 | !--------------------------------------------------------------------- |
---|
| 169 | ! Initialization: Open file and get dimensions only |
---|
| 170 | !--------------------------------------------------------------------- |
---|
| 171 | |
---|
| 172 | iflag = nf90_open( TRIM( TRIM( cfilenames(jj) ) ), nf90_nowrite, & |
---|
| 173 | & i_file_id ) |
---|
| 174 | |
---|
| 175 | IF ( iflag /= nf90_noerr ) THEN |
---|
| 176 | |
---|
| 177 | IF ( ldignmis ) THEN |
---|
| 178 | inpfiles(jj)%nobs = 0 |
---|
| 179 | CALL ctl_warn( 'File ' // TRIM( TRIM( cfilenames(jj) ) ) // & |
---|
| 180 | & ' not found' ) |
---|
| 181 | ELSE |
---|
| 182 | CALL ctl_stop( 'File ' // TRIM( TRIM( cfilenames(jj) ) ) // & |
---|
| 183 | & ' not found' ) |
---|
| 184 | ENDIF |
---|
| 185 | |
---|
| 186 | ELSE |
---|
| 187 | |
---|
| 188 | !------------------------------------------------------------------ |
---|
| 189 | ! Close the file since it is opened in read_proffile |
---|
| 190 | !------------------------------------------------------------------ |
---|
| 191 | |
---|
| 192 | iflag = nf90_close( i_file_id ) |
---|
| 193 | |
---|
| 194 | !------------------------------------------------------------------ |
---|
| 195 | ! Read the profile file into inpfiles |
---|
| 196 | !------------------------------------------------------------------ |
---|
| 197 | IF ( kformat == 0 ) THEN |
---|
| 198 | CALL init_obfbdata( inpfiles(jj) ) |
---|
| 199 | IF(lwp) THEN |
---|
| 200 | WRITE(numout,*) |
---|
| 201 | WRITE(numout,*)'Reading from feedback file :', & |
---|
| 202 | & TRIM( cfilenames(jj) ) |
---|
| 203 | ENDIF |
---|
| 204 | CALL read_obfbdata( TRIM( cfilenames(jj) ), inpfiles(jj), & |
---|
| 205 | & ldgrid = .TRUE. ) |
---|
| 206 | IF ( ldmod .AND. ( inpfiles(jj)%nadd == 0 ) ) THEN |
---|
| 207 | CALL ctl_stop( 'Model not in input data' ) |
---|
| 208 | RETURN |
---|
| 209 | ENDIF |
---|
| 210 | ELSEIF ( kformat == 1 ) THEN |
---|
| 211 | CALL read_taondbc( TRIM( cfilenames(jj) ), inpfiles(jj), & |
---|
| 212 | & numout, lwp, .TRUE. ) |
---|
| 213 | ELSE |
---|
| 214 | CALL ctl_stop( 'File format unknown' ) |
---|
| 215 | ENDIF |
---|
| 216 | |
---|
| 217 | !------------------------------------------------------------------ |
---|
| 218 | ! Change longitude (-180,180) |
---|
| 219 | !------------------------------------------------------------------ |
---|
| 220 | |
---|
| 221 | DO ji = 1, inpfiles(jj)%nobs |
---|
| 222 | |
---|
| 223 | IF ( inpfiles(jj)%plam(ji) < -180. ) & |
---|
| 224 | & inpfiles(jj)%plam(ji) = inpfiles(jj)%plam(ji) + 360. |
---|
| 225 | |
---|
| 226 | IF ( inpfiles(jj)%plam(ji) > 180. ) & |
---|
| 227 | & inpfiles(jj)%plam(ji) = inpfiles(jj)%plam(ji) - 360. |
---|
| 228 | |
---|
| 229 | END DO |
---|
| 230 | |
---|
| 231 | !------------------------------------------------------------------ |
---|
| 232 | ! Calculate the date (change eventually) |
---|
| 233 | !------------------------------------------------------------------ |
---|
| 234 | cl_refdate=inpfiles(jj)%cdjuldref(1:8) |
---|
| 235 | READ(cl_refdate,'(I8)') irefdate(jj) |
---|
| 236 | |
---|
| 237 | CALL ddatetoymdhms( ddobsini, iyea, imon, iday, ihou, imin, isec ) |
---|
| 238 | CALL greg2jul( isec, imin, ihou, iday, imon, iyea, djulini(jj), & |
---|
| 239 | & krefdate = irefdate(jj) ) |
---|
| 240 | CALL ddatetoymdhms( ddobsend, iyea, imon, iday, ihou, imin, isec ) |
---|
| 241 | CALL greg2jul( isec, imin, ihou, iday, imon, iyea, djulend(jj), & |
---|
| 242 | & krefdate = irefdate(jj) ) |
---|
| 243 | |
---|
| 244 | IF ( ldavtimset ) THEN |
---|
| 245 | DO ji = 1, inpfiles(jj)%nobs |
---|
| 246 | ! |
---|
| 247 | ! for daily averaged data force the time |
---|
| 248 | ! to be the end of the day |
---|
| 249 | ! |
---|
| 250 | inpfiles(jj)%ptim(ji) = & |
---|
| 251 | & INT(inpfiles(jj)%ptim(ji)) + 1 |
---|
| 252 | END DO |
---|
| 253 | ENDIF |
---|
| 254 | |
---|
| 255 | IF ( inpfiles(jj)%nobs > 0 ) THEN |
---|
| 256 | inpfiles(jj)%iproc = -1 |
---|
| 257 | inpfiles(jj)%iobsi = -1 |
---|
| 258 | inpfiles(jj)%iobsj = -1 |
---|
| 259 | ENDIF |
---|
| 260 | inowin = 0 |
---|
| 261 | DO ji = 1, inpfiles(jj)%nobs |
---|
| 262 | IF ( ( inpfiles(jj)%ptim(ji) > djulini(jj) ) .AND. & |
---|
| 263 | & ( inpfiles(jj)%ptim(ji) <= djulend(jj) ) ) THEN |
---|
| 264 | inowin = inowin + 1 |
---|
| 265 | ENDIF |
---|
| 266 | END DO |
---|
| 267 | ALLOCATE( zlam(inowin) ) |
---|
| 268 | ALLOCATE( zphi(inowin) ) |
---|
| 269 | ALLOCATE( iobsiu(inowin) ) |
---|
| 270 | ALLOCATE( iobsju(inowin) ) |
---|
| 271 | ALLOCATE( iprocu(inowin) ) |
---|
| 272 | ALLOCATE( iobsiv(inowin) ) |
---|
| 273 | ALLOCATE( iobsjv(inowin) ) |
---|
| 274 | ALLOCATE( iprocv(inowin) ) |
---|
| 275 | inowin = 0 |
---|
| 276 | DO ji = 1, inpfiles(jj)%nobs |
---|
| 277 | IF ( ( inpfiles(jj)%ptim(ji) > djulini(jj) ) .AND. & |
---|
| 278 | & ( inpfiles(jj)%ptim(ji) <= djulend(jj) ) ) THEN |
---|
| 279 | inowin = inowin + 1 |
---|
| 280 | zlam(inowin) = inpfiles(jj)%plam(ji) |
---|
| 281 | zphi(inowin) = inpfiles(jj)%pphi(ji) |
---|
| 282 | ENDIF |
---|
| 283 | END DO |
---|
| 284 | |
---|
| 285 | CALL obs_grid_search( inowin, zlam, zphi, iobsiu, iobsju, iprocu, & |
---|
| 286 | & 'U' ) |
---|
| 287 | CALL obs_grid_search( inowin, zlam, zphi, iobsiv, iobsjv, iprocv, & |
---|
| 288 | & 'V' ) |
---|
| 289 | |
---|
| 290 | inowin = 0 |
---|
| 291 | DO ji = 1, inpfiles(jj)%nobs |
---|
| 292 | IF ( ( inpfiles(jj)%ptim(ji) > djulini(jj) ) .AND. & |
---|
| 293 | & ( inpfiles(jj)%ptim(ji) <= djulend(jj) ) ) THEN |
---|
| 294 | inowin = inowin + 1 |
---|
| 295 | inpfiles(jj)%iproc(ji,1) = iprocu(inowin) |
---|
| 296 | inpfiles(jj)%iobsi(ji,1) = iobsiu(inowin) |
---|
| 297 | inpfiles(jj)%iobsj(ji,1) = iobsju(inowin) |
---|
| 298 | inpfiles(jj)%iproc(ji,2) = iprocv(inowin) |
---|
| 299 | inpfiles(jj)%iobsi(ji,2) = iobsiv(inowin) |
---|
| 300 | inpfiles(jj)%iobsj(ji,2) = iobsjv(inowin) |
---|
| 301 | IF ( inpfiles(jj)%iproc(ji,1) /= & |
---|
| 302 | & inpfiles(jj)%iproc(ji,2) ) THEN |
---|
| 303 | CALL ctl_stop( 'Error in obs_read_vel:', & |
---|
| 304 | & 'U and V observation on different processors') |
---|
| 305 | ENDIF |
---|
| 306 | ENDIF |
---|
| 307 | END DO |
---|
| 308 | DEALLOCATE( zlam, zphi, iobsiu, iobsju, iprocu, & |
---|
| 309 | & iobsiv, iobsjv, iprocv ) |
---|
| 310 | |
---|
| 311 | DO ji = 1, inpfiles(jj)%nobs |
---|
| 312 | IF ( ( inpfiles(jj)%ptim(ji) > djulini(jj) ) .AND. & |
---|
| 313 | & ( inpfiles(jj)%ptim(ji) <= djulend(jj) ) ) THEN |
---|
| 314 | IF ( nproc == 0 ) THEN |
---|
| 315 | IF ( inpfiles(jj)%iproc(ji,1) > nproc ) CYCLE |
---|
| 316 | ELSE |
---|
| 317 | IF ( inpfiles(jj)%iproc(ji,1) /= nproc ) CYCLE |
---|
| 318 | ENDIF |
---|
| 319 | llvalprof = .FALSE. |
---|
| 320 | IF ( ( ( nproc == 0 ) .AND. & |
---|
| 321 | & ( inpfiles(jj)%iproc(ji,1) <= nproc ) ) .OR. & |
---|
| 322 | & ( inpfiles(jj)%iproc(ji,1) == nproc ) ) THEN |
---|
| 323 | loop_uv_count : DO ij = 1,inpfiles(jj)%nlev |
---|
| 324 | IF ( inpfiles(jj)%pdep(ij,ji) >= 6000. ) & |
---|
| 325 | & CYCLE |
---|
| 326 | IF ( ( inpfiles(jj)%ivlqc(ij,ji,1) <= 2 ) .AND. & |
---|
| 327 | & ( inpfiles(jj)%ivlqc(ij,ji,2) <= 2 ) .AND. & |
---|
| 328 | & ( inpfiles(jj)%idqc(ij,ji) <= 2 ) ) THEN |
---|
| 329 | iuv3dt = iuv3dt + 1 |
---|
| 330 | llvalprof = .TRUE. |
---|
| 331 | ENDIF |
---|
| 332 | END DO loop_uv_count |
---|
| 333 | ENDIF |
---|
| 334 | IF ( llvalprof ) iprof = iprof + 1 |
---|
| 335 | ENDIF |
---|
| 336 | END DO |
---|
| 337 | |
---|
| 338 | ENDIF |
---|
| 339 | |
---|
| 340 | END DO prof_files |
---|
| 341 | |
---|
| 342 | !----------------------------------------------------------------------- |
---|
| 343 | ! Get the time ordered indices of the input data |
---|
| 344 | !----------------------------------------------------------------------- |
---|
| 345 | |
---|
| 346 | !--------------------------------------------------------------------- |
---|
| 347 | ! Loop over input data files to count total number of profiles |
---|
| 348 | !--------------------------------------------------------------------- |
---|
| 349 | iproftot = 0 |
---|
| 350 | DO jj = 1, inobf |
---|
| 351 | DO ji = 1, inpfiles(jj)%nobs |
---|
| 352 | IF ( ( inpfiles(jj)%ptim(ji) > djulini(jj) ) .AND. & |
---|
| 353 | & ( inpfiles(jj)%ptim(ji) <= djulend(jj) ) ) THEN |
---|
| 354 | iproftot = iproftot + 1 |
---|
| 355 | ENDIF |
---|
| 356 | END DO |
---|
| 357 | END DO |
---|
| 358 | |
---|
| 359 | ALLOCATE( iindx(iproftot), ifileidx(iproftot), & |
---|
| 360 | & iprofidx(iproftot), zdat(iproftot) ) |
---|
| 361 | jk = 0 |
---|
| 362 | DO jj = 1, inobf |
---|
| 363 | DO ji = 1, inpfiles(jj)%nobs |
---|
| 364 | IF ( ( inpfiles(jj)%ptim(ji) > djulini(jj) ) .AND. & |
---|
| 365 | & ( inpfiles(jj)%ptim(ji) <= djulend(jj) ) ) THEN |
---|
| 366 | jk = jk + 1 |
---|
| 367 | ifileidx(jk) = jj |
---|
| 368 | iprofidx(jk) = ji |
---|
| 369 | zdat(jk) = inpfiles(jj)%ptim(ji) |
---|
| 370 | ENDIF |
---|
| 371 | END DO |
---|
| 372 | END DO |
---|
| 373 | CALL sort_dp_indx( iproftot, & |
---|
| 374 | & zdat, & |
---|
| 375 | & iindx ) |
---|
| 376 | |
---|
| 377 | iv3dt(:) = -1 |
---|
| 378 | iv3dt(1:2) = iuv3dt |
---|
| 379 | CALL obs_prof_alloc( profdata, kvars, kextr, iprof, iv3dt, & |
---|
| 380 | & kstp, jpi, jpj, jpk ) |
---|
| 381 | |
---|
| 382 | ! * Read obs/positions, QC, all variable and assign to profdata |
---|
| 383 | |
---|
| 384 | profdata%nprof = 0 |
---|
| 385 | profdata%nvprot(:) = 0 |
---|
| 386 | |
---|
| 387 | iprof = 0 |
---|
| 388 | |
---|
| 389 | iuv3dt = 0 |
---|
| 390 | itypuv (:) = 0 |
---|
| 391 | itypuvmpp(:) = 0 |
---|
| 392 | |
---|
| 393 | DO jk = 1, iproftot |
---|
| 394 | |
---|
| 395 | jj = ifileidx(iindx(jk)) |
---|
| 396 | ji = iprofidx(iindx(jk)) |
---|
| 397 | IF ( ( inpfiles(jj)%ptim(ji) > djulini(jj) ) .AND. & |
---|
| 398 | & ( inpfiles(jj)%ptim(ji) <= djulend(jj) ) ) THEN |
---|
| 399 | |
---|
| 400 | IF ( nproc == 0 ) THEN |
---|
| 401 | IF ( inpfiles(jj)%iproc(ji,1) > nproc ) CYCLE |
---|
| 402 | ELSE |
---|
| 403 | IF ( inpfiles(jj)%iproc(ji,1) /= nproc ) CYCLE |
---|
| 404 | ENDIF |
---|
| 405 | |
---|
| 406 | llvalprof = .FALSE. |
---|
| 407 | |
---|
| 408 | loop_prof : DO ij = 1, inpfiles(jj)%nlev |
---|
| 409 | |
---|
| 410 | IF ( inpfiles(jj)%pdep(ij,ji) >= 6000. ) & |
---|
| 411 | & CYCLE |
---|
| 412 | |
---|
| 413 | IF ( ( inpfiles(jj)%ivlqc(ij,ji,1) <= 2 ) .AND. & |
---|
| 414 | & ( inpfiles(jj)%ivlqc(ij,ji,2) <= 2 ) .AND. & |
---|
| 415 | & ( inpfiles(jj)%idqc(ij,ji) <= 2 ) ) THEN |
---|
| 416 | |
---|
| 417 | llvalprof = .TRUE. |
---|
| 418 | EXIT loop_prof |
---|
| 419 | |
---|
| 420 | ENDIF |
---|
| 421 | |
---|
| 422 | END DO loop_prof |
---|
| 423 | |
---|
| 424 | ! Set profile information |
---|
| 425 | |
---|
| 426 | IF ( llvalprof ) THEN |
---|
| 427 | |
---|
| 428 | iprof = iprof + 1 |
---|
| 429 | |
---|
| 430 | CALL jul2greg( isec, & |
---|
| 431 | & imin, & |
---|
| 432 | & ihou, & |
---|
| 433 | & iday, & |
---|
| 434 | & imon, & |
---|
| 435 | & iyea, & |
---|
| 436 | & inpfiles(jj)%ptim(ji), & |
---|
| 437 | & irefdate(jj) ) |
---|
| 438 | |
---|
| 439 | |
---|
| 440 | ! Profile time coordinates |
---|
| 441 | profdata%nyea(iprof) = iyea |
---|
| 442 | profdata%nmon(iprof) = imon |
---|
| 443 | profdata%nday(iprof) = iday |
---|
| 444 | profdata%nhou(iprof) = ihou |
---|
| 445 | profdata%nmin(iprof) = imin |
---|
| 446 | |
---|
| 447 | ! Profile space coordinates |
---|
| 448 | profdata%rlam(iprof) = inpfiles(jj)%plam(ji) |
---|
| 449 | profdata%rphi(iprof) = inpfiles(jj)%pphi(ji) |
---|
| 450 | |
---|
| 451 | ! Coordinate search parameters |
---|
| 452 | profdata%mi (iprof,1) = inpfiles(jj)%iobsi(ji,1) |
---|
| 453 | profdata%mj (iprof,1) = inpfiles(jj)%iobsj(ji,1) |
---|
| 454 | profdata%mi (iprof,2) = inpfiles(jj)%iobsi(ji,2) |
---|
| 455 | profdata%mj (iprof,2) = inpfiles(jj)%iobsj(ji,2) |
---|
| 456 | |
---|
| 457 | ! Profile WMO number |
---|
| 458 | profdata%cwmo(iprof) = inpfiles(jj)%cdwmo(ji) |
---|
| 459 | |
---|
| 460 | ! Instrument type |
---|
| 461 | READ( inpfiles(jj)%cdtyp(ji), '(I4)' ) itype |
---|
| 462 | profdata%ntyp(iprof) = itype |
---|
| 463 | |
---|
| 464 | ! QC stuff |
---|
| 465 | |
---|
| 466 | profdata%nqc(iprof) = inpfiles(jj)%ioqc(ji) |
---|
| 467 | profdata%nqcf(:,iprof) = inpfiles(jj)%ioqcf(:,ji) |
---|
| 468 | profdata%ipqc(iprof) = inpfiles(jj)%ipqc(ji) |
---|
| 469 | profdata%ipqcf(:,iprof) = inpfiles(jj)%ipqcf(:,ji) |
---|
| 470 | profdata%itqc(iprof) = inpfiles(jj)%itqc(ji) |
---|
| 471 | profdata%itqcf(:,iprof) = inpfiles(jj)%itqcf(:,ji) |
---|
| 472 | profdata%ivqc(iprof,:) = inpfiles(jj)%ivqc(ji,:) |
---|
| 473 | profdata%ivqcf(:,iprof,:) = inpfiles(jj)%ivqcf(:,ji,:) |
---|
| 474 | |
---|
| 475 | ! Bookkeeping data to match profiles |
---|
| 476 | profdata%npidx(iprof) = iprof |
---|
| 477 | profdata%npfil(iprof) = iindx(jk) |
---|
| 478 | |
---|
| 479 | ! Observation QC flag (whole profile) |
---|
| 480 | profdata%nqc(iprof) = 0 !TODO |
---|
| 481 | |
---|
| 482 | loop_uv : DO ij = 1, inpfiles(jj)%nlev |
---|
| 483 | |
---|
| 484 | IF ( inpfiles(jj)%pdep(ij,ji) >= 6000. ) & |
---|
| 485 | & CYCLE |
---|
| 486 | |
---|
| 487 | IF ( ( inpfiles(jj)%ivlqc(ij,ji,1) <= 2 ) .AND. & |
---|
| 488 | & ( inpfiles(jj)%ivlqc(ij,ji,2) <= 2 ) .AND. & |
---|
| 489 | & ( inpfiles(jj)%idqc(ij,ji) <= 2 ) ) THEN |
---|
| 490 | iuv3dt = iuv3dt + 1 |
---|
| 491 | ELSE |
---|
| 492 | CYCLE |
---|
| 493 | ENDIF |
---|
| 494 | |
---|
| 495 | ! Depth of U observation |
---|
| 496 | profdata%var(1)%vdep(iuv3dt) = & |
---|
| 497 | & inpfiles(jj)%pdep(ij,ji) |
---|
| 498 | |
---|
| 499 | ! Depth of U observation QC |
---|
| 500 | profdata%var(1)%idqc(iuv3dt) = & |
---|
| 501 | & inpfiles(jj)%idqc(ij,ji) |
---|
| 502 | |
---|
| 503 | ! Depth of U observation QC flags |
---|
| 504 | profdata%var(1)%idqcf(:,iuv3dt) = & |
---|
| 505 | & inpfiles(jj)%idqcf(:,ij,ji) |
---|
| 506 | |
---|
| 507 | ! Profile index |
---|
| 508 | profdata%var(1)%nvpidx(iuv3dt) = iprof |
---|
| 509 | |
---|
| 510 | ! Vertical index in original profile |
---|
| 511 | profdata%var(1)%nvlidx(iuv3dt) = ij |
---|
| 512 | |
---|
| 513 | ! Profile U value |
---|
| 514 | |
---|
| 515 | profdata%var(1)%vobs(iuv3dt) = & |
---|
| 516 | & inpfiles(jj)%pob(ij,ji,1) |
---|
| 517 | IF ( ldmod ) THEN |
---|
| 518 | profdata%var(1)%vmod(iuv3dt) = & |
---|
| 519 | & inpfiles(jj)%padd(ij,ji,1,1) |
---|
| 520 | ENDIF |
---|
| 521 | |
---|
| 522 | ! Profile U qc |
---|
| 523 | profdata%var(1)%nvqc(iuv3dt) = & |
---|
| 524 | & inpfiles(jj)%ivlqc(ij,ji,1) |
---|
| 525 | |
---|
| 526 | ! Profile U qc flags |
---|
| 527 | profdata%var(1)%nvqcf(:,iuv3dt) = & |
---|
| 528 | & inpfiles(jj)%ivlqcf(:,ij,ji,1) |
---|
| 529 | |
---|
| 530 | |
---|
| 531 | ! Depth of V observation |
---|
| 532 | profdata%var(2)%vdep(iuv3dt) = & |
---|
| 533 | & inpfiles(jj)%pdep(ij,ji) |
---|
| 534 | |
---|
| 535 | ! Depth of V observation QC |
---|
| 536 | profdata%var(2)%idqc(iuv3dt) = & |
---|
| 537 | & inpfiles(jj)%idqc(ij,ji) |
---|
| 538 | |
---|
| 539 | ! Depth of V observation QC flags |
---|
| 540 | profdata%var(2)%idqcf(:,iuv3dt) = & |
---|
| 541 | & inpfiles(jj)%idqcf(:,ij,ji) |
---|
| 542 | |
---|
| 543 | ! Profile index |
---|
| 544 | profdata%var(2)%nvpidx(iuv3dt) = iprof |
---|
| 545 | |
---|
| 546 | ! Vertical index in original profile |
---|
| 547 | profdata%var(2)%nvlidx(iuv3dt) = ij |
---|
| 548 | |
---|
| 549 | ! Profile V value |
---|
| 550 | profdata%var(2)%vobs(iuv3dt) = & |
---|
| 551 | & inpfiles(jj)%pob(ij,ji,2) |
---|
| 552 | IF ( ldmod ) THEN |
---|
| 553 | profdata%var(2)%vmod(iuv3dt) = & |
---|
| 554 | & inpfiles(jj)%padd(ij,ji,1,2) |
---|
| 555 | ENDIF |
---|
| 556 | |
---|
| 557 | ! Profile V qc |
---|
| 558 | profdata%var(2)%nvqc(iuv3dt) = & |
---|
| 559 | & inpfiles(jj)%ivlqc(ij,ji,2) |
---|
| 560 | |
---|
| 561 | ! Profile V qc flags |
---|
| 562 | profdata%var(2)%nvqcf(:,iuv3dt) = & |
---|
| 563 | & inpfiles(jj)%ivlqcf(:,ij,ji,2) |
---|
| 564 | |
---|
| 565 | ! Observation type |
---|
| 566 | itypuv( profdata%ntyp(iprof) + 1 ) = & |
---|
| 567 | & itypuv( profdata%ntyp(iprof) + 1 ) + 1 |
---|
| 568 | |
---|
| 569 | |
---|
| 570 | END DO loop_uv |
---|
| 571 | |
---|
| 572 | ENDIF |
---|
| 573 | |
---|
| 574 | ENDIF |
---|
| 575 | |
---|
| 576 | END DO |
---|
| 577 | |
---|
| 578 | !----------------------------------------------------------------------- |
---|
| 579 | ! Sum up over processors |
---|
| 580 | !----------------------------------------------------------------------- |
---|
| 581 | |
---|
| 582 | CALL obs_mpp_sum_integer ( iuv3dt, iuv3dtmpp ) |
---|
| 583 | |
---|
| 584 | CALL obs_mpp_sum_integers( itypuv, itypuvmpp, ntyp1770 + 1 ) |
---|
| 585 | |
---|
| 586 | !----------------------------------------------------------------------- |
---|
| 587 | ! Output number of observations. |
---|
| 588 | !----------------------------------------------------------------------- |
---|
| 589 | IF(lwp) THEN |
---|
| 590 | WRITE(numout,*) |
---|
| 591 | WRITE(numout,'(1X,A)') 'Profile U,V velocity data' |
---|
| 592 | WRITE(numout,'(1X,A)') '-------------------------' |
---|
| 593 | WRITE(numout,*) |
---|
| 594 | DO ji = 0, ntyp1770 |
---|
| 595 | IF ( itypuvmpp(ji+1) > 0 ) THEN |
---|
| 596 | WRITE(numout,'(1X,A3,A3,I8)') ctypshort(ji), ' = ', & |
---|
| 597 | & itypuvmpp(ji+1) |
---|
| 598 | ENDIF |
---|
| 599 | END DO |
---|
| 600 | WRITE(numout,'(1X,A)') '--------------' |
---|
| 601 | WRITE(numout,'(1X,A6,I8)') & |
---|
| 602 | & 'Total profile UV data = ',& |
---|
| 603 | & iuv3dtmpp |
---|
| 604 | WRITE(numout,'(1X,A)') '--------------' |
---|
| 605 | ENDIF |
---|
| 606 | |
---|
| 607 | profdata%nvprot(1) = iuv3dt |
---|
| 608 | profdata%nvprot(2) = iuv3dt |
---|
| 609 | profdata%nvprotmpp(1) = iuv3dtmpp |
---|
| 610 | profdata%nvprotmpp(2) = iuv3dtmpp |
---|
| 611 | profdata%nprof = iprof |
---|
| 612 | |
---|
| 613 | !----------------------------------------------------------------------- |
---|
| 614 | ! Model level search |
---|
| 615 | !----------------------------------------------------------------------- |
---|
| 616 | CALL obs_level_search( jpk, gdept_0, & |
---|
| 617 | & profdata%nvprot(1), profdata%var(1)%vdep, & |
---|
| 618 | & profdata%var(1)%mvk ) |
---|
| 619 | CALL obs_level_search( jpk, gdept_0, & |
---|
| 620 | & profdata%nvprot(2), profdata%var(2)%vdep, & |
---|
| 621 | & profdata%var(2)%mvk ) |
---|
| 622 | |
---|
| 623 | !----------------------------------------------------------------------- |
---|
| 624 | ! Set model equivalent to 99999 |
---|
| 625 | !----------------------------------------------------------------------- |
---|
| 626 | IF ( .NOT. ldmod ) THEN |
---|
| 627 | DO jvar = 1, kvars |
---|
| 628 | profdata%var(jvar)%vmod(:) = fbrmdi |
---|
| 629 | END DO |
---|
| 630 | ENDIF |
---|
| 631 | !----------------------------------------------------------------------- |
---|
| 632 | ! Deallocate temporary data |
---|
| 633 | !----------------------------------------------------------------------- |
---|
| 634 | DEALLOCATE( ifileidx, iprofidx, zdat ) |
---|
| 635 | |
---|
| 636 | !----------------------------------------------------------------------- |
---|
| 637 | ! Deallocate input data |
---|
| 638 | !----------------------------------------------------------------------- |
---|
| 639 | DO jj = 1, inobf |
---|
| 640 | CALL dealloc_obfbdata( inpfiles(jj) ) |
---|
| 641 | END DO |
---|
| 642 | DEALLOCATE( inpfiles ) |
---|
| 643 | |
---|
| 644 | END SUBROUTINE obs_rea_vel_dri |
---|
| 645 | |
---|
| 646 | END MODULE obs_read_vel |
---|