Changeset 6225 for branches/2014/dev_r4704_NOC5_MPP_BDY_UPDATE/NEMOGCM/NEMO/OPA_SRC/OBS/obs_write.F90
- Timestamp:
- 2016-01-08T10:35:19+01:00 (9 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
branches/2014/dev_r4704_NOC5_MPP_BDY_UPDATE/NEMOGCM/NEMO/OPA_SRC/OBS/obs_write.F90
r2287 r6225 6 6 7 7 !!---------------------------------------------------------------------- 8 !! cdf_wri_p3d : Write profile observation diagnostics in NetCDF format 9 !! obs_wri_sla : Write SLA observation related diagnostics 10 !! obs_wri_sst : Write SST observation related diagnostics 11 !! obs_wri_seaice: Write seaice observation related diagnostics 12 !! cdf_wri_vel : Write velocity observation diagnostics in NetCDF format 8 !! obs_wri_prof : Write profile observations in feedback format 9 !! obs_wri_surf : Write surface observations in feedback format 10 !! obs_wri_stats : Print basic statistics on the data being written out 13 11 !!---------------------------------------------------------------------- 14 12 … … 29 27 USE obs_conv ! Conversion between units 30 28 USE obs_const 31 USE obs_ sla_types32 USE obs_rot_vel ! Rotation of velocities29 USE obs_mpp ! MPP support routines for observation diagnostics 30 USE lib_mpp ! MPP routines 33 31 34 32 IMPLICIT NONE … … 36 34 !! * Routine accessibility 37 35 PRIVATE 38 PUBLIC obs_wri_p3d, & ! Write profile observation related diagnostics 39 & obs_wri_sla, & ! Write SLA observation related diagnostics 40 & obs_wri_sst, & ! Write SST observation related diagnostics 41 & obs_wri_sss, & ! Write SSS observation related diagnostics 42 & obs_wri_seaice, & ! Write seaice observation related diagnostics 43 & obs_wri_vel, & ! Write velocity observation related diagnostics 36 PUBLIC obs_wri_prof, & ! Write profile observation files 37 & obs_wri_surf, & ! Write surface observation files 44 38 & obswriinfo 45 39 … … 60 54 CONTAINS 61 55 62 SUBROUTINE obs_wri_p 3d( cprefix,profdata, padd, pext )56 SUBROUTINE obs_wri_prof( profdata, padd, pext ) 63 57 !!----------------------------------------------------------------------- 64 58 !! 65 !! *** ROUTINE obs_wri_p3d *** 66 !! 67 !! ** Purpose : Write temperature and salinity (profile) observation 68 !! related diagnostics 59 !! *** ROUTINE obs_wri_prof *** 60 !! 61 !! ** Purpose : Write profile feedback files 69 62 !! 70 63 !! ** Method : NetCDF … … 79 72 !! ! 07-03 (K. Mogensen) General handling of profiles 80 73 !! ! 09-01 (K. Mogensen) New feedback format 74 !! ! 15-02 (M. Martin) Combined routine for writing profiles 81 75 !!----------------------------------------------------------------------- 82 76 83 !! * Modules used84 85 77 !! * Arguments 86 CHARACTER(LEN=*), INTENT(IN) :: cprefix ! Prefix for output files87 78 TYPE(obs_prof), INTENT(INOUT) :: profdata ! Full set of profile data 88 79 TYPE(obswriinfo), OPTIONAL :: padd ! Additional info for each variable 89 80 TYPE(obswriinfo), OPTIONAL :: pext ! Extra info 90 81 91 82 !! * Local declarations 92 83 TYPE(obfbdata) :: fbdata 93 CHARACTER(LEN=40) :: cfname 84 CHARACTER(LEN=40) :: clfname 85 CHARACTER(LEN=6) :: clfiletype 94 86 INTEGER :: ilevel 95 87 INTEGER :: jvar … … 99 91 INTEGER :: ja 100 92 INTEGER :: je 93 INTEGER :: iadd 94 INTEGER :: iext 101 95 REAL(wp) :: zpres 102 INTEGER :: nadd103 INTEGER :: next104 96 105 97 IF ( PRESENT( padd ) ) THEN 106 nadd = padd%inum98 iadd = padd%inum 107 99 ELSE 108 nadd = 0100 iadd = 0 109 101 ENDIF 110 102 111 103 IF ( PRESENT( pext ) ) THEN 112 next = pext%inum104 iext = pext%inum 113 105 ELSE 114 next = 0115 ENDIF 116 106 iext = 0 107 ENDIF 108 117 109 CALL init_obfbdata( fbdata ) 118 110 … … 122 114 ilevel = MAX( ilevel, MAXVAL( profdata%var(jvar)%nvlidx(:) ) ) 123 115 END DO 124 CALL alloc_obfbdata( fbdata, 2, profdata%nprof, ilevel, & 125 & 1 + nadd, 1 + next, .TRUE. ) 126 127 fbdata%cname(1) = 'POTM' 128 fbdata%cname(2) = 'PSAL' 129 fbdata%coblong(1) = 'Potential temperature' 130 fbdata%coblong(2) = 'Practical salinity' 131 fbdata%cobunit(1) = 'Degrees centigrade' 132 fbdata%cobunit(2) = 'PSU' 133 fbdata%cextname(1) = 'TEMP' 134 fbdata%cextlong(1) = 'Insitu temperature' 135 fbdata%cextunit(1) = 'Degrees centigrade' 136 DO je = 1, next 137 fbdata%cextname(1+je) = pext%cdname(je) 138 fbdata%cextlong(1+je) = pext%cdlong(je,1) 139 fbdata%cextunit(1+je) = pext%cdunit(je,1) 140 END DO 116 117 SELECT CASE ( TRIM(profdata%cvars(1)) ) 118 CASE('POTM') 119 120 clfiletype='profb' 121 CALL alloc_obfbdata( fbdata, 2, profdata%nprof, ilevel, & 122 & 1 + iadd, 1 + iext, .TRUE. ) 123 fbdata%cname(1) = profdata%cvars(1) 124 fbdata%cname(2) = profdata%cvars(2) 125 fbdata%coblong(1) = 'Potential temperature' 126 fbdata%coblong(2) = 'Practical salinity' 127 fbdata%cobunit(1) = 'Degrees centigrade' 128 fbdata%cobunit(2) = 'PSU' 129 fbdata%cextname(1) = 'TEMP' 130 fbdata%cextlong(1) = 'Insitu temperature' 131 fbdata%cextunit(1) = 'Degrees centigrade' 132 fbdata%caddlong(1,1) = 'Model interpolated potential temperature' 133 fbdata%caddlong(1,2) = 'Model interpolated practical salinity' 134 fbdata%caddunit(1,1) = 'Degrees centigrade' 135 fbdata%caddunit(1,2) = 'PSU' 136 fbdata%cgrid(:) = 'T' 137 DO je = 1, iext 138 fbdata%cextname(1+je) = pext%cdname(je) 139 fbdata%cextlong(1+je) = pext%cdlong(je,1) 140 fbdata%cextunit(1+je) = pext%cdunit(je,1) 141 END DO 142 DO ja = 1, iadd 143 fbdata%caddname(1+ja) = padd%cdname(ja) 144 DO jvar = 1, 2 145 fbdata%caddlong(1+ja,jvar) = padd%cdlong(ja,jvar) 146 fbdata%caddunit(1+ja,jvar) = padd%cdunit(ja,jvar) 147 END DO 148 END DO 149 150 CASE('UVEL') 151 152 clfiletype='velfb' 153 CALL alloc_obfbdata( fbdata, 2, profdata%nprof, ilevel, 1, 0, .TRUE. ) 154 fbdata%cname(1) = profdata%cvars(1) 155 fbdata%cname(2) = profdata%cvars(2) 156 fbdata%coblong(1) = 'Zonal velocity' 157 fbdata%coblong(2) = 'Meridional velocity' 158 fbdata%cobunit(1) = 'm/s' 159 fbdata%cobunit(2) = 'm/s' 160 DO je = 1, iext 161 fbdata%cextname(je) = pext%cdname(je) 162 fbdata%cextlong(je) = pext%cdlong(je,1) 163 fbdata%cextunit(je) = pext%cdunit(je,1) 164 END DO 165 fbdata%caddlong(1,1) = 'Model interpolated zonal velocity' 166 fbdata%caddlong(1,2) = 'Model interpolated meridional velocity' 167 fbdata%caddunit(1,1) = 'm/s' 168 fbdata%caddunit(1,2) = 'm/s' 169 fbdata%cgrid(1) = 'U' 170 fbdata%cgrid(2) = 'V' 171 DO ja = 1, iadd 172 fbdata%caddname(1+ja) = padd%cdname(ja) 173 fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1) 174 fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1) 175 END DO 176 177 END SELECT 178 141 179 fbdata%caddname(1) = 'Hx' 142 fbdata%caddlong(1,1) = 'Model interpolated potential temperature' 143 fbdata%caddlong(1,2) = 'Model interpolated practical salinity' 144 fbdata%caddunit(1,1) = 'Degrees centigrade' 145 fbdata%caddunit(1,2) = 'PSU' 146 fbdata%cgrid(:) = 'T' 147 DO ja = 1, nadd 148 fbdata%caddname(1+ja) = padd%cdname(ja) 149 DO jvar = 1, 2 150 fbdata%caddlong(1+ja,jvar) = padd%cdlong(ja,jvar) 151 fbdata%caddunit(1+ja,jvar) = padd%cdunit(ja,jvar) 152 END DO 153 END DO 154 155 WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc 180 181 WRITE(clfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(clfiletype), nproc 156 182 157 183 IF(lwp) THEN 158 184 WRITE(numout,*) 159 WRITE(numout,*)'obs_wri_p 3d:'185 WRITE(numout,*)'obs_wri_prof :' 160 186 WRITE(numout,*)'~~~~~~~~~~~~~' 161 WRITE(numout,*)'Writing profile feedback file : ',TRIM(cfname)162 ENDIF 163 164 ! Transform obs_prof data structure into obfb data structure187 WRITE(numout,*)'Writing '//TRIM(clfiletype)//' feedback file : ',TRIM(clfname) 188 ENDIF 189 190 ! Transform obs_prof data structure into obfb data structure 165 191 fbdata%cdjuldref = '19500101000000' 166 192 DO jo = 1, profdata%nprof … … 219 245 ENDIF 220 246 fbdata%iobsk(ik,jo,jvar) = profdata%var(jvar)%mvk(jk) 221 DO ja = 1, nadd247 DO ja = 1, iadd 222 248 fbdata%padd(ik,jo,1+ja,jvar) = & 223 249 & profdata%var(jvar)%vext(jk,padd%ipoint(ja)) 224 250 END DO 225 DO je = 1, next251 DO je = 1, iext 226 252 fbdata%pext(ik,jo,1+je) = & 227 253 & profdata%var(jvar)%vext(jk,pext%ipoint(je)) 228 254 END DO 229 IF ( jvar == 1 ) THEN 255 IF ( ( jvar == 1 ) .AND. & 256 & ( TRIM(profdata%cvars(1)) == 'POTM' ) ) THEN 230 257 fbdata%pext(ik,jo,1) = profdata%var(jvar)%vext(jk,1) 231 258 ENDIF … … 234 261 END DO 235 262 236 ! Convert insitu temperature to potential temperature using the model 237 ! salinity if no potential temperature 238 DO jo = 1, fbdata%nobs 239 IF ( fbdata%pphi(jo) < 9999.0 ) THEN 240 DO jk = 1, fbdata%nlev 241 IF ( ( fbdata%pob(jk,jo,1) >= 9999.0 ) .AND. & 242 & ( fbdata%pdep(jk,jo) < 9999.0 ) .AND. & 243 & ( fbdata%padd(jk,jo,1,2) < 9999.0 ) .AND. & 244 & ( fbdata%pext(jk,jo,1) < 9999.0 ) ) THEN 245 zpres = dep_to_p( REAL(fbdata%pdep(jk,jo),wp), & 246 & REAL(fbdata%pphi(jo),wp) ) 247 fbdata%pob(jk,jo,1) = potemp( & 248 & REAL(fbdata%padd(jk,jo,1,2), wp), & 249 & REAL(fbdata%pext(jk,jo,1), wp), & 250 & zpres, 0.0_wp ) 251 ENDIF 252 END DO 253 ENDIF 254 END DO 255 263 IF ( TRIM(profdata%cvars(1)) == 'POTM' ) THEN 264 ! Convert insitu temperature to potential temperature using the model 265 ! salinity if no potential temperature 266 DO jo = 1, fbdata%nobs 267 IF ( fbdata%pphi(jo) < 9999.0 ) THEN 268 DO jk = 1, fbdata%nlev 269 IF ( ( fbdata%pob(jk,jo,1) >= 9999.0 ) .AND. & 270 & ( fbdata%pdep(jk,jo) < 9999.0 ) .AND. & 271 & ( fbdata%padd(jk,jo,1,2) < 9999.0 ) .AND. & 272 & ( fbdata%pext(jk,jo,1) < 9999.0 ) ) THEN 273 zpres = dep_to_p( REAL(fbdata%pdep(jk,jo),wp), & 274 & REAL(fbdata%pphi(jo),wp) ) 275 fbdata%pob(jk,jo,1) = potemp( & 276 & REAL(fbdata%padd(jk,jo,1,2), wp), & 277 & REAL(fbdata%pext(jk,jo,1), wp), & 278 & zpres, 0.0_wp ) 279 ENDIF 280 END DO 281 ENDIF 282 END DO 283 ENDIF 284 256 285 ! Write the obfbdata structure 257 CALL write_obfbdata( cfname, fbdata ) 258 286 CALL write_obfbdata( clfname, fbdata ) 287 288 ! Output some basic statistics 289 CALL obs_wri_stats( fbdata ) 290 259 291 CALL dealloc_obfbdata( fbdata ) 260 261 END SUBROUTINE obs_wri_p 3d262 263 SUBROUTINE obs_wri_s la( cprefix, sladata, padd, pext )292 293 END SUBROUTINE obs_wri_prof 294 295 SUBROUTINE obs_wri_surf( surfdata, padd, pext ) 264 296 !!----------------------------------------------------------------------- 265 297 !! 266 !! *** ROUTINE obs_wri_sla *** 267 !! 268 !! ** Purpose : Write SLA observation diagnostics 269 !! related 298 !! *** ROUTINE obs_wri_surf *** 299 !! 300 !! ** Purpose : Write surface observation files 270 301 !! 271 302 !! ** Method : NetCDF … … 275 306 !! ! 07-03 (K. Mogensen) Original 276 307 !! ! 09-01 (K. Mogensen) New feedback format. 308 !! ! 15-02 (M. Martin) Combined surface writing routine. 277 309 !!----------------------------------------------------------------------- 278 310 … … 281 313 282 314 !! * Arguments 283 CHARACTER(LEN=*), INTENT(IN) :: cprefix ! Prefix for output files 284 TYPE(obs_surf), INTENT(INOUT) :: sladata ! Full set of SLAa 315 TYPE(obs_surf), INTENT(INOUT) :: surfdata ! Full set of surface data 285 316 TYPE(obswriinfo), OPTIONAL :: padd ! Additional info for each variable 286 317 TYPE(obswriinfo), OPTIONAL :: pext ! Extra info … … 288 319 !! * Local declarations 289 320 TYPE(obfbdata) :: fbdata 290 CHARACTER(LEN=40) :: cfname ! netCDF filename 291 CHARACTER(LEN=12), PARAMETER :: cpname = 'obs_wri_sla' 321 CHARACTER(LEN=40) :: clfname ! netCDF filename 322 CHARACTER(LEN=6) :: clfiletype 323 CHARACTER(LEN=12), PARAMETER :: cpname = 'obs_wri_surf' 292 324 INTEGER :: jo 293 325 INTEGER :: ja 294 326 INTEGER :: je 295 INTEGER :: nadd296 INTEGER :: next327 INTEGER :: iadd 328 INTEGER :: iext 297 329 298 330 IF ( PRESENT( padd ) ) THEN 299 nadd = padd%inum331 iadd = padd%inum 300 332 ELSE 301 nadd = 0333 iadd = 0 302 334 ENDIF 303 335 304 336 IF ( PRESENT( pext ) ) THEN 305 next = pext%inum337 iext = pext%inum 306 338 ELSE 307 next = 0339 iext = 0 308 340 ENDIF 309 341 310 342 CALL init_obfbdata( fbdata ) 311 343 312 CALL alloc_obfbdata( fbdata, 1, sladata%nsurf, 1, & 313 & 2 + nadd, 1 + next, .TRUE. ) 314 315 fbdata%cname(1) = 'SLA' 316 fbdata%coblong(1) = 'Sea level anomaly' 317 fbdata%cobunit(1) = 'Metres' 318 fbdata%cextname(1) = 'MDT' 319 fbdata%cextlong(1) = 'Mean dynamic topography' 320 fbdata%cextunit(1) = 'Metres' 321 DO je = 1, next 322 fbdata%cextname(1+je) = pext%cdname(je) 323 fbdata%cextlong(1+je) = pext%cdlong(je,1) 324 fbdata%cextunit(1+je) = pext%cdunit(je,1) 325 END DO 344 SELECT CASE ( TRIM(surfdata%cvars(1)) ) 345 CASE('SLA') 346 347 CALL alloc_obfbdata( fbdata, 1, surfdata%nsurf, 1, & 348 & 2 + iadd, 1 + iext, .TRUE. ) 349 350 clfiletype = 'slafb' 351 fbdata%cname(1) = surfdata%cvars(1) 352 fbdata%coblong(1) = 'Sea level anomaly' 353 fbdata%cobunit(1) = 'Metres' 354 fbdata%cextname(1) = 'MDT' 355 fbdata%cextlong(1) = 'Mean dynamic topography' 356 fbdata%cextunit(1) = 'Metres' 357 DO je = 1, iext 358 fbdata%cextname(je) = pext%cdname(je) 359 fbdata%cextlong(je) = pext%cdlong(je,1) 360 fbdata%cextunit(je) = pext%cdunit(je,1) 361 END DO 362 fbdata%caddlong(1,1) = 'Model interpolated SSH - MDT' 363 fbdata%caddunit(1,1) = 'Metres' 364 fbdata%caddname(2) = 'SSH' 365 fbdata%caddlong(2,1) = 'Model Sea surface height' 366 fbdata%caddunit(2,1) = 'Metres' 367 fbdata%cgrid(1) = 'T' 368 DO ja = 1, iadd 369 fbdata%caddname(2+ja) = padd%cdname(ja) 370 fbdata%caddlong(2+ja,1) = padd%cdlong(ja,1) 371 fbdata%caddunit(2+ja,1) = padd%cdunit(ja,1) 372 END DO 373 374 CASE('SST') 375 376 CALL alloc_obfbdata( fbdata, 1, surfdata%nsurf, 1, & 377 & 1 + iadd, iext, .TRUE. ) 378 379 clfiletype = 'sstfb' 380 fbdata%cname(1) = surfdata%cvars(1) 381 fbdata%coblong(1) = 'Sea surface temperature' 382 fbdata%cobunit(1) = 'Degree centigrade' 383 DO je = 1, iext 384 fbdata%cextname(je) = pext%cdname(je) 385 fbdata%cextlong(je) = pext%cdlong(je,1) 386 fbdata%cextunit(je) = pext%cdunit(je,1) 387 END DO 388 fbdata%caddlong(1,1) = 'Model interpolated SST' 389 fbdata%caddunit(1,1) = 'Degree centigrade' 390 fbdata%cgrid(1) = 'T' 391 DO ja = 1, iadd 392 fbdata%caddname(1+ja) = padd%cdname(ja) 393 fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1) 394 fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1) 395 END DO 396 397 CASE('ICECON') 398 399 CALL alloc_obfbdata( fbdata, 1, surfdata%nsurf, 1, & 400 & 1 + iadd, iext, .TRUE. ) 401 402 clfiletype = 'sicfb' 403 fbdata%cname(1) = surfdata%cvars(1) 404 fbdata%coblong(1) = 'Sea ice' 405 fbdata%cobunit(1) = 'Fraction' 406 DO je = 1, iext 407 fbdata%cextname(je) = pext%cdname(je) 408 fbdata%cextlong(je) = pext%cdlong(je,1) 409 fbdata%cextunit(je) = pext%cdunit(je,1) 410 END DO 411 fbdata%caddlong(1,1) = 'Model interpolated ICE' 412 fbdata%caddunit(1,1) = 'Fraction' 413 fbdata%cgrid(1) = 'T' 414 DO ja = 1, iadd 415 fbdata%caddname(1+ja) = padd%cdname(ja) 416 fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1) 417 fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1) 418 END DO 419 420 END SELECT 421 326 422 fbdata%caddname(1) = 'Hx' 327 fbdata%caddlong(1,1) = 'Model interpolated SSH - MDT' 328 fbdata%caddunit(1,1) = 'Metres' 329 fbdata%caddname(2) = 'SSH' 330 fbdata%caddlong(2,1) = 'Model Sea surface height' 331 fbdata%caddunit(2,1) = 'Metres' 332 fbdata%cgrid(1) = 'T' 333 DO ja = 1, nadd 334 fbdata%caddname(2+ja) = padd%cdname(ja) 335 fbdata%caddlong(2+ja,1) = padd%cdlong(ja,1) 336 fbdata%caddunit(2+ja,1) = padd%cdunit(ja,1) 337 END DO 338 339 WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc 423 424 WRITE(clfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(clfiletype), nproc 340 425 341 426 IF(lwp) THEN 342 427 WRITE(numout,*) 343 WRITE(numout,*)'obs_wri_s la:'428 WRITE(numout,*)'obs_wri_surf :' 344 429 WRITE(numout,*)'~~~~~~~~~~~~~' 345 WRITE(numout,*)'Writing SLA feedback file : ',TRIM(cfname)346 ENDIF 347 348 ! Transform obs_prof data structure into obfbdata structure430 WRITE(numout,*)'Writing '//TRIM(surfdata%cvars(1))//' feedback file : ',TRIM(clfname) 431 ENDIF 432 433 ! Transform surf data structure into obfbdata structure 349 434 fbdata%cdjuldref = '19500101000000' 350 DO jo = 1, s ladata%nsurf351 fbdata%plam(jo) = s ladata%rlam(jo)352 fbdata%pphi(jo) = s ladata%rphi(jo)353 WRITE(fbdata%cdtyp(jo),'(I4)') s ladata%ntyp(jo)435 DO jo = 1, surfdata%nsurf 436 fbdata%plam(jo) = surfdata%rlam(jo) 437 fbdata%pphi(jo) = surfdata%rphi(jo) 438 WRITE(fbdata%cdtyp(jo),'(I4)') surfdata%ntyp(jo) 354 439 fbdata%ivqc(jo,:) = 0 355 440 fbdata%ivqcf(:,jo,:) = 0 356 IF ( s ladata%nqc(jo) > 10 ) THEN441 IF ( surfdata%nqc(jo) > 10 ) THEN 357 442 fbdata%ioqc(jo) = 4 358 443 fbdata%ioqcf(1,jo) = 0 359 fbdata%ioqcf(2,jo) = s ladata%nqc(jo) - 10444 fbdata%ioqcf(2,jo) = surfdata%nqc(jo) - 10 360 445 ELSE 361 fbdata%ioqc(jo) = s ladata%nqc(jo)446 fbdata%ioqc(jo) = surfdata%nqc(jo) 362 447 fbdata%ioqcf(:,jo) = 0 363 448 ENDIF … … 366 451 fbdata%itqc(jo) = 0 367 452 fbdata%itqcf(:,jo) = 0 368 fbdata%cdwmo(jo) = s ladata%cwmo(jo)369 fbdata%kindex(jo) = s ladata%nsfil(jo)453 fbdata%cdwmo(jo) = surfdata%cwmo(jo) 454 fbdata%kindex(jo) = surfdata%nsfil(jo) 370 455 IF (ln_grid_global) THEN 371 fbdata%iobsi(jo,1) = s ladata%mi(jo)372 fbdata%iobsj(jo,1) = s ladata%mj(jo)456 fbdata%iobsi(jo,1) = surfdata%mi(jo) 457 fbdata%iobsj(jo,1) = surfdata%mj(jo) 373 458 ELSE 374 fbdata%iobsi(jo,1) = mig(s ladata%mi(jo))375 fbdata%iobsj(jo,1) = mjg(s ladata%mj(jo))459 fbdata%iobsi(jo,1) = mig(surfdata%mi(jo)) 460 fbdata%iobsj(jo,1) = mjg(surfdata%mj(jo)) 376 461 ENDIF 377 462 CALL greg2jul( 0, & 378 & s ladata%nmin(jo), &379 & s ladata%nhou(jo), &380 & s ladata%nday(jo), &381 & s ladata%nmon(jo), &382 & s ladata%nyea(jo), &463 & surfdata%nmin(jo), & 464 & surfdata%nhou(jo), & 465 & surfdata%nday(jo), & 466 & surfdata%nmon(jo), & 467 & surfdata%nyea(jo), & 383 468 & fbdata%ptim(jo), & 384 469 & krefdate = 19500101 ) 385 fbdata%padd(1,jo,1,1) = s ladata%rmod(jo,1)386 fbdata%padd(1,jo,2,1) = sladata%rext(jo,1)387 fbdata%pob(1,jo,1) = s ladata%robs(jo,1)470 fbdata%padd(1,jo,1,1) = surfdata%rmod(jo,1) 471 IF ( TRIM(surfdata%cvars(1)) == 'SLA' ) fbdata%padd(1,jo,2,1) = surfdata%rext(jo,1) 472 fbdata%pob(1,jo,1) = surfdata%robs(jo,1) 388 473 fbdata%pdep(1,jo) = 0.0 389 474 fbdata%idqc(1,jo) = 0 390 475 fbdata%idqcf(:,1,jo) = 0 391 IF ( s ladata%nqc(jo) > 10 ) THEN476 IF ( surfdata%nqc(jo) > 10 ) THEN 392 477 fbdata%ivqc(jo,1) = 4 393 478 fbdata%ivlqc(1,jo,1) = 4 394 479 fbdata%ivlqcf(1,1,jo,1) = 0 395 fbdata%ivlqcf(2,1,jo,1) = s ladata%nqc(jo) - 10480 fbdata%ivlqcf(2,1,jo,1) = surfdata%nqc(jo) - 10 396 481 ELSE 397 fbdata%ivqc(jo,1) = s ladata%nqc(jo)398 fbdata%ivlqc(1,jo,1) = s ladata%nqc(jo)482 fbdata%ivqc(jo,1) = surfdata%nqc(jo) 483 fbdata%ivlqc(1,jo,1) = surfdata%nqc(jo) 399 484 fbdata%ivlqcf(:,1,jo,1) = 0 400 485 ENDIF 401 486 fbdata%iobsk(1,jo,1) = 0 402 fbdata%pext(1,jo,1) = sladata%rext(jo,2)403 DO ja = 1, nadd487 IF ( TRIM(surfdata%cvars(1)) == 'SLA' ) fbdata%pext(1,jo,1) = surfdata%rext(jo,2) 488 DO ja = 1, iadd 404 489 fbdata%padd(1,jo,2+ja,1) = & 405 & s ladata%rext(jo,padd%ipoint(ja))406 END DO 407 DO je = 1, next490 & surfdata%rext(jo,padd%ipoint(ja)) 491 END DO 492 DO je = 1, iext 408 493 fbdata%pext(1,jo,1+je) = & 409 & s ladata%rext(jo,pext%ipoint(je))494 & surfdata%rext(jo,pext%ipoint(je)) 410 495 END DO 411 496 END DO 412 497 413 498 ! Write the obfbdata structure 414 CALL write_obfbdata( cfname, fbdata ) 499 CALL write_obfbdata( clfname, fbdata ) 500 501 ! Output some basic statistics 502 CALL obs_wri_stats( fbdata ) 415 503 416 504 CALL dealloc_obfbdata( fbdata ) 417 505 418 END SUBROUTINE obs_wri_s la419 420 SUBROUTINE obs_wri_s st( cprefix, sstdata, padd, pext)506 END SUBROUTINE obs_wri_surf 507 508 SUBROUTINE obs_wri_stats( fbdata ) 421 509 !!----------------------------------------------------------------------- 422 510 !! 423 !! *** ROUTINE obs_wri_sst *** 424 !! 425 !! ** Purpose : Write SST observation diagnostics 426 !! related 427 !! 428 !! ** Method : NetCDF 511 !! *** ROUTINE obs_wri_stats *** 512 !! 513 !! ** Purpose : Output some basic statistics of the data being written out 514 !! 515 !! ** Method : 429 516 !! 430 517 !! ** Action : 431 518 !! 432 !! ! 07-07 (S. Ricci) Original 433 !! ! 09-01 (K. Mogensen) New feedback format. 519 !! ! 2014-08 (D. Lea) Initial version 434 520 !!----------------------------------------------------------------------- 435 521 436 !! * Modules used437 IMPLICIT NONE438 439 522 !! * Arguments 440 CHARACTER(LEN=*), INTENT(IN) :: cprefix ! Prefix for output files441 TYPE(obs_surf), INTENT(INOUT) :: sstdata ! Full set of SST442 TYPE(obswriinfo), OPTIONAL :: padd ! Additional info for each variable443 TYPE(obswriinfo), OPTIONAL :: pext ! Extra info444 445 !! * Local declarations446 523 TYPE(obfbdata) :: fbdata 447 CHARACTER(LEN=40) :: cfname ! netCDF filename 448 CHARACTER(LEN=12), PARAMETER :: cpname = 'obs_wri_sst' 524 525 !! * Local declarations 526 INTEGER :: jvar 449 527 INTEGER :: jo 450 INTEGER :: ja 451 INTEGER :: je 452 INTEGER :: nadd 453 INTEGER :: next 454 455 IF ( PRESENT( padd ) ) THEN 456 nadd = padd%inum 457 ELSE 458 nadd = 0 459 ENDIF 460 461 IF ( PRESENT( pext ) ) THEN 462 next = pext%inum 463 ELSE 464 next = 0 465 ENDIF 466 467 CALL init_obfbdata( fbdata ) 468 469 CALL alloc_obfbdata( fbdata, 1, sstdata%nsurf, 1, & 470 & 1 + nadd, next, .TRUE. ) 471 472 fbdata%cname(1) = 'SST' 473 fbdata%coblong(1) = 'Sea surface temperature' 474 fbdata%cobunit(1) = 'Degree centigrade' 475 DO je = 1, next 476 fbdata%cextname(je) = pext%cdname(je) 477 fbdata%cextlong(je) = pext%cdlong(je,1) 478 fbdata%cextunit(je) = pext%cdunit(je,1) 479 END DO 480 fbdata%caddname(1) = 'Hx' 481 fbdata%caddlong(1,1) = 'Model interpolated SST' 482 fbdata%caddunit(1,1) = 'Degree centigrade' 483 fbdata%cgrid(1) = 'T' 484 DO ja = 1, nadd 485 fbdata%caddname(1+ja) = padd%cdname(ja) 486 fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1) 487 fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1) 488 END DO 489 490 WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc 491 492 IF(lwp) THEN 493 WRITE(numout,*) 494 WRITE(numout,*)'obs_wri_sst :' 495 WRITE(numout,*)'~~~~~~~~~~~~~' 496 WRITE(numout,*)'Writing SST feedback file : ',TRIM(cfname) 497 ENDIF 498 499 ! Transform obs_prof data structure into obfbdata structure 500 fbdata%cdjuldref = '19500101000000' 501 DO jo = 1, sstdata%nsurf 502 fbdata%plam(jo) = sstdata%rlam(jo) 503 fbdata%pphi(jo) = sstdata%rphi(jo) 504 WRITE(fbdata%cdtyp(jo),'(I4)') sstdata%ntyp(jo) 505 fbdata%ivqc(jo,:) = 0 506 fbdata%ivqcf(:,jo,:) = 0 507 IF ( sstdata%nqc(jo) > 10 ) THEN 508 fbdata%ioqc(jo) = 4 509 fbdata%ioqcf(1,jo) = 0 510 fbdata%ioqcf(2,jo) = sstdata%nqc(jo) - 10 511 ELSE 512 fbdata%ioqc(jo) = MAX(sstdata%nqc(jo),1) 513 fbdata%ioqcf(:,jo) = 0 528 INTEGER :: jk 529 INTEGER :: inumgoodobs 530 INTEGER :: inumgoodobsmpp 531 REAL(wp) :: zsumx 532 REAL(wp) :: zsumx2 533 REAL(wp) :: zomb 534 535 536 IF (lwp) THEN 537 WRITE(numout,*) '' 538 WRITE(numout,*) 'obs_wri_stats :' 539 WRITE(numout,*) '~~~~~~~~~~~~~~~' 540 ENDIF 541 542 DO jvar = 1, fbdata%nvar 543 zsumx=0.0_wp 544 zsumx2=0.0_wp 545 inumgoodobs=0 546 DO jo = 1, fbdata%nobs 547 DO jk = 1, fbdata%nlev 548 IF ( ( fbdata%pob(jk,jo,jvar) < 9999.0 ) .AND. & 549 & ( fbdata%pdep(jk,jo) < 9999.0 ) .AND. & 550 & ( fbdata%padd(jk,jo,1,jvar) < 9999.0 ) ) THEN 551 552 zomb=fbdata%pob(jk, jo, jvar)-fbdata%padd(jk, jo, 1, jvar) 553 zsumx=zsumx+zomb 554 zsumx2=zsumx2+zomb**2 555 inumgoodobs=inumgoodobs+1 556 ENDIF 557 ENDDO 558 ENDDO 559 560 CALL obs_mpp_sum_integer( inumgoodobs, inumgoodobsmpp ) 561 CALL mpp_sum(zsumx) 562 CALL mpp_sum(zsumx2) 563 564 IF (lwp) THEN 565 WRITE(numout,*) 'Type: ',fbdata%cname(jvar),' Total number of good observations: ',inumgoodobsmpp 566 WRITE(numout,*) 'Overall mean obs minus model of the good observations: ',zsumx/inumgoodobsmpp 567 WRITE(numout,*) 'Overall RMS obs minus model of the good observations: ',sqrt( zsumx2/inumgoodobsmpp ) 568 WRITE(numout,*) '' 514 569 ENDIF 515 fbdata%ipqc(jo) = 0 516 fbdata%ipqcf(:,jo) = 0 517 fbdata%itqc(jo) = 0 518 fbdata%itqcf(:,jo) = 0 519 fbdata%cdwmo(jo) = '' 520 fbdata%kindex(jo) = sstdata%nsfil(jo) 521 IF (ln_grid_global) THEN 522 fbdata%iobsi(jo,1) = sstdata%mi(jo) 523 fbdata%iobsj(jo,1) = sstdata%mj(jo) 524 ELSE 525 fbdata%iobsi(jo,1) = mig(sstdata%mi(jo)) 526 fbdata%iobsj(jo,1) = mjg(sstdata%mj(jo)) 527 ENDIF 528 CALL greg2jul( 0, & 529 & sstdata%nmin(jo), & 530 & sstdata%nhou(jo), & 531 & sstdata%nday(jo), & 532 & sstdata%nmon(jo), & 533 & sstdata%nyea(jo), & 534 & fbdata%ptim(jo), & 535 & krefdate = 19500101 ) 536 fbdata%padd(1,jo,1,1) = sstdata%rmod(jo,1) 537 fbdata%pob(1,jo,1) = sstdata%robs(jo,1) 538 fbdata%pdep(1,jo) = 0.0 539 fbdata%idqc(1,jo) = 0 540 fbdata%idqcf(:,1,jo) = 0 541 IF ( sstdata%nqc(jo) > 10 ) THEN 542 fbdata%ivqc(jo,1) = 4 543 fbdata%ivlqc(1,jo,1) = 4 544 fbdata%ivlqcf(1,1,jo,1) = 0 545 fbdata%ivlqcf(2,1,jo,1) = sstdata%nqc(jo) - 10 546 ELSE 547 fbdata%ivqc(jo,1) = MAX(sstdata%nqc(jo),1) 548 fbdata%ivlqc(1,jo,1) = MAX(sstdata%nqc(jo),1) 549 fbdata%ivlqcf(:,1,jo,1) = 0 550 ENDIF 551 fbdata%iobsk(1,jo,1) = 0 552 DO ja = 1, nadd 553 fbdata%padd(1,jo,1+ja,1) = & 554 & sstdata%rext(jo,padd%ipoint(ja)) 555 END DO 556 DO je = 1, next 557 fbdata%pext(1,jo,je) = & 558 & sstdata%rext(jo,pext%ipoint(je)) 559 END DO 560 561 END DO 562 563 ! Write the obfbdata structure 564 565 CALL write_obfbdata( cfname, fbdata ) 566 567 CALL dealloc_obfbdata( fbdata ) 568 569 END SUBROUTINE obs_wri_sst 570 571 SUBROUTINE obs_wri_sss 572 END SUBROUTINE obs_wri_sss 573 574 SUBROUTINE obs_wri_seaice( cprefix, seaicedata, padd, pext ) 575 !!----------------------------------------------------------------------- 576 !! 577 !! *** ROUTINE obs_wri_seaice *** 578 !! 579 !! ** Purpose : Write sea ice observation diagnostics 580 !! related 581 !! 582 !! ** Method : NetCDF 583 !! 584 !! ** Action : 585 !! 586 !! ! 07-07 (S. Ricci) Original 587 !! ! 09-01 (K. Mogensen) New feedback format. 588 !!----------------------------------------------------------------------- 589 590 !! * Modules used 591 IMPLICIT NONE 592 593 !! * Arguments 594 CHARACTER(LEN=*), INTENT(IN) :: cprefix ! Prefix for output files 595 TYPE(obs_surf), INTENT(INOUT) :: seaicedata ! Full set of sea ice 596 TYPE(obswriinfo), OPTIONAL :: padd ! Additional info for each variable 597 TYPE(obswriinfo), OPTIONAL :: pext ! Extra info 598 599 !! * Local declarations 600 TYPE(obfbdata) :: fbdata 601 CHARACTER(LEN=40) :: cfname ! netCDF filename 602 CHARACTER(LEN=12), PARAMETER :: cpname = 'obs_wri_seaice' 603 INTEGER :: jo 604 INTEGER :: ja 605 INTEGER :: je 606 INTEGER :: nadd 607 INTEGER :: next 608 609 IF ( PRESENT( padd ) ) THEN 610 nadd = padd%inum 611 ELSE 612 nadd = 0 613 ENDIF 614 615 IF ( PRESENT( pext ) ) THEN 616 next = pext%inum 617 ELSE 618 next = 0 619 ENDIF 620 621 CALL init_obfbdata( fbdata ) 622 623 CALL alloc_obfbdata( fbdata, 1, seaicedata%nsurf, 1, 1, 0, .TRUE. ) 624 625 fbdata%cname(1) = 'SEAICE' 626 fbdata%coblong(1) = 'Sea ice' 627 fbdata%cobunit(1) = 'Fraction' 628 DO je = 1, next 629 fbdata%cextname(je) = pext%cdname(je) 630 fbdata%cextlong(je) = pext%cdlong(je,1) 631 fbdata%cextunit(je) = pext%cdunit(je,1) 632 END DO 633 fbdata%caddname(1) = 'Hx' 634 fbdata%caddlong(1,1) = 'Model interpolated ICE' 635 fbdata%caddunit(1,1) = 'Fraction' 636 fbdata%cgrid(1) = 'T' 637 DO ja = 1, nadd 638 fbdata%caddname(1+ja) = padd%cdname(ja) 639 fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1) 640 fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1) 641 END DO 642 643 WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc 644 645 IF(lwp) THEN 646 WRITE(numout,*) 647 WRITE(numout,*)'obs_wri_seaice :' 648 WRITE(numout,*)'~~~~~~~~~~~~~~~~' 649 WRITE(numout,*)'Writing SEAICE feedback file : ',TRIM(cfname) 650 ENDIF 651 652 ! Transform obs_prof data structure into obfbdata structure 653 fbdata%cdjuldref = '19500101000000' 654 DO jo = 1, seaicedata%nsurf 655 fbdata%plam(jo) = seaicedata%rlam(jo) 656 fbdata%pphi(jo) = seaicedata%rphi(jo) 657 WRITE(fbdata%cdtyp(jo),'(I4)') seaicedata%ntyp(jo) 658 fbdata%ivqc(jo,:) = 0 659 fbdata%ivqcf(:,jo,:) = 0 660 IF ( seaicedata%nqc(jo) > 10 ) THEN 661 fbdata%ioqc(jo) = 4 662 fbdata%ioqcf(1,jo) = 0 663 fbdata%ioqcf(2,jo) = seaicedata%nqc(jo) - 10 664 ELSE 665 fbdata%ioqc(jo) = MAX(seaicedata%nqc(jo),1) 666 fbdata%ioqcf(:,jo) = 0 667 ENDIF 668 fbdata%ipqc(jo) = 0 669 fbdata%ipqcf(:,jo) = 0 670 fbdata%itqc(jo) = 0 671 fbdata%itqcf(:,jo) = 0 672 fbdata%cdwmo(jo) = '' 673 fbdata%kindex(jo) = seaicedata%nsfil(jo) 674 IF (ln_grid_global) THEN 675 fbdata%iobsi(jo,1) = seaicedata%mi(jo) 676 fbdata%iobsj(jo,1) = seaicedata%mj(jo) 677 ELSE 678 fbdata%iobsi(jo,1) = mig(seaicedata%mi(jo)) 679 fbdata%iobsj(jo,1) = mjg(seaicedata%mj(jo)) 680 ENDIF 681 CALL greg2jul( 0, & 682 & seaicedata%nmin(jo), & 683 & seaicedata%nhou(jo), & 684 & seaicedata%nday(jo), & 685 & seaicedata%nmon(jo), & 686 & seaicedata%nyea(jo), & 687 & fbdata%ptim(jo), & 688 & krefdate = 19500101 ) 689 fbdata%padd(1,jo,1,1) = seaicedata%rmod(jo,1) 690 fbdata%pob(1,jo,1) = seaicedata%robs(jo,1) 691 fbdata%pdep(1,jo) = 0.0 692 fbdata%idqc(1,jo) = 0 693 fbdata%idqcf(:,1,jo) = 0 694 IF ( seaicedata%nqc(jo) > 10 ) THEN 695 fbdata%ivlqc(1,jo,1) = 4 696 fbdata%ivlqcf(1,1,jo,1) = 0 697 fbdata%ivlqcf(2,1,jo,1) = seaicedata%nqc(jo) - 10 698 ELSE 699 fbdata%ivlqc(1,jo,1) = MAX(seaicedata%nqc(jo),1) 700 fbdata%ivlqcf(:,1,jo,1) = 0 701 ENDIF 702 fbdata%iobsk(1,jo,1) = 0 703 DO ja = 1, nadd 704 fbdata%padd(1,jo,1+ja,1) = & 705 & seaicedata%rext(jo,padd%ipoint(ja)) 706 END DO 707 DO je = 1, next 708 fbdata%pext(1,jo,je) = & 709 & seaicedata%rext(jo,pext%ipoint(je)) 710 END DO 711 712 END DO 713 714 ! Write the obfbdata structure 715 CALL write_obfbdata( cfname, fbdata ) 716 717 CALL dealloc_obfbdata( fbdata ) 718 719 END SUBROUTINE obs_wri_seaice 720 721 SUBROUTINE obs_wri_vel( cprefix, profdata, k2dint, padd, pext ) 722 !!----------------------------------------------------------------------- 723 !! 724 !! *** ROUTINE obs_wri_p3d *** 725 !! 726 !! ** Purpose : Write current (profile) observation 727 !! related diagnostics 728 !! 729 !! ** Method : NetCDF 730 !! 731 !! ** Action : 732 !! 733 !! History : 734 !! ! 09-01 (K. Mogensen) New feedback format routine 735 !!----------------------------------------------------------------------- 736 737 !! * Modules used 738 739 !! * Arguments 740 CHARACTER(LEN=*), INTENT(IN) :: cprefix ! Prefix for output files 741 TYPE(obs_prof), INTENT(INOUT) :: profdata ! Full set of profile data 742 INTEGER, INTENT(IN) :: k2dint ! Horizontal interpolation method 743 TYPE(obswriinfo), OPTIONAL :: padd ! Additional info for each variable 744 TYPE(obswriinfo), OPTIONAL :: pext ! Extra info 745 746 !! * Local declarations 747 TYPE(obfbdata) :: fbdata 748 CHARACTER(LEN=40) :: cfname 749 INTEGER :: ilevel 750 INTEGER :: jvar 751 INTEGER :: jk 752 INTEGER :: ik 753 INTEGER :: jo 754 INTEGER :: ja 755 INTEGER :: je 756 INTEGER :: nadd 757 INTEGER :: next 758 REAL(wp) :: zpres 759 REAL(wp), DIMENSION(:), ALLOCATABLE :: & 760 & zu, & 761 & zv 762 763 IF ( PRESENT( padd ) ) THEN 764 nadd = padd%inum 765 ELSE 766 nadd = 0 767 ENDIF 768 769 IF ( PRESENT( pext ) ) THEN 770 next = pext%inum 771 ELSE 772 next = 0 773 ENDIF 774 775 CALL init_obfbdata( fbdata ) 776 777 ! Find maximum level 778 ilevel = 0 779 DO jvar = 1, 2 780 ilevel = MAX( ilevel, MAXVAL( profdata%var(jvar)%nvlidx(:) ) ) 781 END DO 782 CALL alloc_obfbdata( fbdata, 2, profdata%nprof, ilevel, 2, 0, .TRUE. ) 783 784 fbdata%cname(1) = 'UVEL' 785 fbdata%cname(2) = 'VVEL' 786 fbdata%coblong(1) = 'Zonal velocity' 787 fbdata%coblong(2) = 'Meridional velocity' 788 fbdata%cobunit(1) = 'm/s' 789 fbdata%cobunit(2) = 'm/s' 790 DO je = 1, next 791 fbdata%cextname(je) = pext%cdname(je) 792 fbdata%cextlong(je) = pext%cdlong(je,1) 793 fbdata%cextunit(je) = pext%cdunit(je,1) 794 END DO 795 fbdata%caddname(1) = 'Hx' 796 fbdata%caddlong(1,1) = 'Model interpolated zonal velocity' 797 fbdata%caddlong(1,2) = 'Model interpolated meridional velocity' 798 fbdata%caddunit(1,1) = 'm/s' 799 fbdata%caddunit(1,2) = 'm/s' 800 fbdata%caddname(2) = 'HxG' 801 fbdata%caddlong(2,1) = 'Model interpolated zonal velocity (model grid)' 802 fbdata%caddlong(2,2) = 'Model interpolated meridional velocity (model grid)' 803 fbdata%caddunit(2,1) = 'm/s' 804 fbdata%caddunit(2,2) = 'm/s' 805 fbdata%cgrid(1) = 'U' 806 fbdata%cgrid(2) = 'V' 807 DO ja = 1, nadd 808 fbdata%caddname(2+ja) = padd%cdname(ja) 809 fbdata%caddlong(2+ja,1) = padd%cdlong(ja,1) 810 fbdata%caddunit(2+ja,1) = padd%cdunit(ja,1) 811 END DO 812 813 WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc 814 815 IF(lwp) THEN 816 WRITE(numout,*) 817 WRITE(numout,*)'obs_wri_vel :' 818 WRITE(numout,*)'~~~~~~~~~~~~~' 819 WRITE(numout,*)'Writing velocuty feedback file : ',TRIM(cfname) 820 ENDIF 821 822 ALLOCATE( & 823 & zu(profdata%nvprot(1)), & 824 & zv(profdata%nvprot(2)) & 825 & ) 826 CALL obs_rotvel( profdata, k2dint, zu, zv ) 827 828 ! Transform obs_prof data structure into obfbdata structure 829 fbdata%cdjuldref = '19500101000000' 830 DO jo = 1, profdata%nprof 831 fbdata%plam(jo) = profdata%rlam(jo) 832 fbdata%pphi(jo) = profdata%rphi(jo) 833 WRITE(fbdata%cdtyp(jo),'(I4)') profdata%ntyp(jo) 834 fbdata%ivqc(jo,:) = profdata%ivqc(jo,:) 835 fbdata%ivqcf(:,jo,:) = profdata%ivqcf(:,jo,:) 836 IF ( profdata%nqc(jo) > 10 ) THEN 837 fbdata%ioqc(jo) = 4 838 fbdata%ioqcf(1,jo) = profdata%nqcf(1,jo) 839 fbdata%ioqcf(2,jo) = profdata%nqc(jo) - 10 840 ELSE 841 fbdata%ioqc(jo) = profdata%nqc(jo) 842 fbdata%ioqcf(:,jo) = profdata%nqcf(:,jo) 843 ENDIF 844 fbdata%ipqc(jo) = profdata%ipqc(jo) 845 fbdata%ipqcf(:,jo) = profdata%ipqcf(:,jo) 846 fbdata%itqc(jo) = profdata%itqc(jo) 847 fbdata%itqcf(:,jo) = profdata%itqcf(:,jo) 848 fbdata%cdwmo(jo) = profdata%cwmo(jo) 849 fbdata%kindex(jo) = profdata%npfil(jo) 850 DO jvar = 1, profdata%nvar 851 IF (ln_grid_global) THEN 852 fbdata%iobsi(jo,jvar) = profdata%mi(jo,jvar) 853 fbdata%iobsj(jo,jvar) = profdata%mj(jo,jvar) 854 ELSE 855 fbdata%iobsi(jo,jvar) = mig(profdata%mi(jo,jvar)) 856 fbdata%iobsj(jo,jvar) = mjg(profdata%mj(jo,jvar)) 857 ENDIF 858 END DO 859 CALL greg2jul( 0, & 860 & profdata%nmin(jo), & 861 & profdata%nhou(jo), & 862 & profdata%nday(jo), & 863 & profdata%nmon(jo), & 864 & profdata%nyea(jo), & 865 & fbdata%ptim(jo), & 866 & krefdate = 19500101 ) 867 ! Reform the profiles arrays for output 868 DO jvar = 1, 2 869 DO jk = profdata%npvsta(jo,jvar), profdata%npvend(jo,jvar) 870 ik = profdata%var(jvar)%nvlidx(jk) 871 IF ( jvar == 1 ) THEN 872 fbdata%padd(ik,jo,1,jvar) = zu(jk) 873 ELSE 874 fbdata%padd(ik,jo,1,jvar) = zv(jk) 875 ENDIF 876 fbdata%padd(ik,jo,2,jvar) = profdata%var(jvar)%vmod(jk) 877 fbdata%pob(ik,jo,jvar) = profdata%var(jvar)%vobs(jk) 878 fbdata%pdep(ik,jo) = profdata%var(jvar)%vdep(jk) 879 fbdata%idqc(ik,jo) = profdata%var(jvar)%idqc(jk) 880 fbdata%idqcf(:,ik,jo) = profdata%var(jvar)%idqcf(:,jk) 881 IF ( profdata%var(jvar)%nvqc(jk) > 10 ) THEN 882 fbdata%ivlqc(ik,jo,jvar) = 4 883 fbdata%ivlqcf(1,ik,jo,jvar) = profdata%var(jvar)%nvqcf(1,jk) 884 fbdata%ivlqcf(2,ik,jo,jvar) = profdata%var(jvar)%nvqc(jk) - 10 885 ELSE 886 fbdata%ivlqc(ik,jo,jvar) = profdata%var(jvar)%nvqc(jk) 887 fbdata%ivlqcf(:,ik,jo,jvar) = profdata%var(jvar)%nvqcf(:,jk) 888 ENDIF 889 fbdata%iobsk(ik,jo,jvar) = profdata%var(jvar)%mvk(jk) 890 DO ja = 1, nadd 891 fbdata%padd(ik,jo,2+ja,jvar) = & 892 & profdata%var(jvar)%vext(jk,padd%ipoint(ja)) 893 END DO 894 DO je = 1, next 895 fbdata%pext(ik,jo,je) = & 896 & profdata%var(jvar)%vext(jk,pext%ipoint(je)) 897 END DO 898 END DO 899 END DO 900 END DO 901 902 ! Write the obfbdata structure 903 CALL write_obfbdata( cfname, fbdata ) 904 905 CALL dealloc_obfbdata( fbdata ) 906 907 DEALLOCATE( & 908 & zu, & 909 & zv & 910 & ) 911 912 END SUBROUTINE obs_wri_vel 570 571 ENDDO 572 573 END SUBROUTINE obs_wri_stats 913 574 914 575 END MODULE obs_write
Note: See TracChangeset
for help on using the changeset viewer.