Changeset 11449 for branches/UKMO/dev_r5518_obs_oper_update_addclim/NEMOGCM/NEMO/OPA_SRC/OBS/obs_oper.F90
- Timestamp:
- 2019-08-19T10:54:47+02:00 (5 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
branches/UKMO/dev_r5518_obs_oper_update_addclim/NEMOGCM/NEMO/OPA_SRC/OBS/obs_oper.F90
r9306 r11449 62 62 SUBROUTINE obs_prof_opt( prodatqc, kt, kpi, kpj, kpk, & 63 63 & kit000, kdaystp, kvar, & 64 & pvar, p gdept, pgdepw,&65 & p mask,&64 & pvar, pclim, & 65 & pgdept, pgdepw, pmask, & 66 66 & plam, pphi, & 67 67 & k1dint, k2dint, kdailyavtypes ) … … 137 137 REAL(KIND=wp), INTENT(IN), DIMENSION(kpi,kpj,kpk) :: & 138 138 & pvar, & ! Model field for variable 139 & pclim, & ! Climatology field for variable 139 140 & pmask ! Land-sea mask for variable 140 141 REAL(KIND=wp), INTENT(IN), DIMENSION(kpi,kpj) :: & … … 172 173 REAL(KIND=wp), DIMENSION(kpk) :: & 173 174 & zobsk, & 174 & zobs2k 175 & zobs2k, & 176 & zclm2k 175 177 REAL(KIND=wp), DIMENSION(2,2,1) :: & 176 178 & zweig1, & … … 178 180 REAL(wp), DIMENSION(:,:,:,:), ALLOCATABLE :: & 179 181 & zmask, & 182 & zclim, & 180 183 & zint, & 181 184 & zinm, & … … 189 192 190 193 LOGICAL :: ld_dailyav 194 LOGICAL :: ld_clim 191 195 192 196 !------------------------------------------------------------------------ … … 196 200 inrc = kt - kit000 + 2 197 201 ipro = prodatqc%npstp(inrc) 202 203 ! Check if climatology is available and set flag 204 IF ( SUM( pclim(:,:,:) ) == 0. ) THEN 205 ld_clim = .FALSE. 206 ELSE 207 ld_clim = .TRUE. 208 ENDIF 198 209 199 210 ! Daily average types … … 262 273 & ) 263 274 275 IF ( ld_clim ) ALLOCATE( zclim(2,2,kpk,ipro) ) 276 264 277 DO jobs = prodatqc%nprofup + 1, prodatqc%nprofup + ipro 265 278 iobs = jobs - prodatqc%nprofup … … 286 299 CALL obs_int_comm_3d( 2, 2, ipro, kpi, kpj, kpk, igrdi, igrdj, pgdepw, zgdepw ) 287 300 301 IF ( ld_clim ) THEN 302 CALL obs_int_comm_3d( 2, 2, ipro, kpi, kpj, kpk, igrdi, igrdj, pclim, zclim ) 303 ENDIF 304 288 305 ! At the end of the day also get interpolated means 289 306 IF ( ld_dailyav .AND. idayend == 0 ) THEN … … 349 366 inum_obs = iend - ista + 1 350 367 ALLOCATE(interp_corner(2,2,inum_obs),iv_indic(inum_obs)) 351 368 IF ( ld_clim ) ALLOCATE( interp_corner_clim(2,2,inum_obs) ) 369 352 370 DO iin=1,2 353 371 DO ijn=1,2 … … 358 376 & zobs2k, zgdept(iin,ijn,:,iobs), & 359 377 & zmask(iin,ijn,:,iobs)) 378 379 IF ( ld_clim ) THEN 380 CALL obs_int_z1d_spl( kpk, & 381 & zclim(iin,ijn,:,iobs), & 382 & zclm2k, zgdept(iin,ijn,:,iobs), & 383 & zmask(iin,ijn,:,iobs)) 384 ENDIF 385 360 386 ENDIF 361 387 … … 371 397 & zgdept(iin,ijn,:,iobs), & 372 398 & zmask(iin,ijn,:,iobs)) 373 399 400 IF ( ld_clim ) THEN 401 CALL obs_int_z1d(kpk, iv_indic, k1dint, inum_obs, & 402 & prodatqc%var(kvar)%vdep(ista:iend), & 403 & zclim(iin,ijn,:,iobs), & 404 & zclm2k, interp_corner_clim(iin,ijn,:), & 405 & zgdept(iin,ijn,:,iobs), & 406 & zmask(iin,ijn,:,iobs)) 407 ENDIF 408 374 409 ENDDO 375 410 ENDDO … … 385 420 iend = prodatqc%npvend(jobs,kvar) 386 421 inum_obs = iend - ista + 1 387 ALLOCATE(interp_corner(2,2,inum_obs), iv_indic(inum_obs)) 422 ALLOCATE( interp_corner(2,2,inum_obs), & 423 & iv_indic(inum_obs) ) 424 IF ( ld_clim ) ALLOCATE( interp_corner_clim(2,2,inum_obs) ) 388 425 DO iin=1,2 389 426 DO ijn=1,2 … … 394 431 & zobs2k, zgdept(iin,ijn,:,iobs), & 395 432 & zmask(iin,ijn,:,iobs)) 433 434 IF ( ld_clim ) THEN 435 CALL obs_int_z1d_spl( kpk, & 436 & zclim(iin,ijn,:,iobs),& 437 & zclm2k, zgdept(iin,ijn,:,iobs), & 438 & zmask(iin,ijn,:,iobs)) 439 ENDIF 396 440 397 441 ENDIF … … 408 452 & zgdept(iin,ijn,:,iobs), & 409 453 & zmask(iin,ijn,:,iobs) ) 454 455 IF ( ld_clim ) THEN 456 CALL obs_int_z1d(kpk, iv_indic, k1dint, inum_obs, & 457 & prodatqc%var(kvar)%vdep(ista:iend), & 458 & zclim(iin,ijn,:,iobs), & 459 & zclm2k,interp_corner_clim(iin,ijn,:), & 460 & zgdept(iin,ijn,:,iobs), & 461 & zmask(iin,ijn,:,iobs) ) 462 ENDIF 410 463 411 464 ENDDO … … 451 504 & prodatqc%var(kvar)%vmod(iend:iend) ) 452 505 506 IF ( ld_clim ) THEN 507 CALL obs_int_h2d( 1, 1, zweig, interp_corner_clim(:,:,ikn), & 508 & prodatqc%var(kvar)%vclm(iend:iend) ) 509 ENDIF 510 453 511 ! Set QC flag for any observations found below the bottom 454 512 ! needed as the check here is more strict than that in obs_prep … … 458 516 459 517 DEALLOCATE(interp_corner,iv_indic) 460 518 IF ( ld_clim ) DEALLOCATE( interp_corner_clim ) 519 461 520 ENDIF 462 521 … … 475 534 & ) 476 535 536 IF ( ld_clim ) DEALLOCATE( zclim ) 537 477 538 ! At the end of the day also get interpolated means 478 539 IF ( ld_dailyav .AND. idayend == 0 ) THEN … … 487 548 488 549 SUBROUTINE obs_surf_opt( surfdataqc, kt, kpi, kpj, & 489 & kit000, kdaystp, psurf, p surfmask, &550 & kit000, kdaystp, psurf, pclim, psurfmask, & 490 551 & k2dint, ldnightav, plamscl, pphiscl, & 491 552 & lindegrees ) … … 541 602 REAL(wp), INTENT(IN), DIMENSION(kpi,kpj) :: & 542 603 & psurf, & ! Model surface field 604 & pclim, & ! Climatological surface field 543 605 & psurfmask ! Land-sea mask 544 606 LOGICAL, INTENT(IN) :: ldnightav ! Logical for averaging night-time data … … 569 631 REAL(wp) :: zlam 570 632 REAL(wp) :: zphi 571 REAL(wp), DIMENSION(1) :: zext, zobsmask 633 REAL(wp), DIMENSION(1) :: zext, zobsmask, zclm 572 634 REAL(wp) :: zdaystp 573 635 REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: & … … 577 639 & zsurfm, & 578 640 & zsurftmp, & 641 & zclim, & 579 642 & zglam, & 580 643 & zgphi, & … … 586 649 & zouttmp, & 587 650 & zmeanday ! to compute model sst in region of 24h daylight (pole) 588 651 652 LOGICAL :: ld_clim ! T => climatological data is available 589 653 !------------------------------------------------------------------------ 590 654 ! Local initialization … … 593 657 inrc = kt - kit000 + 2 594 658 isurf = surfdataqc%nsstp(inrc) 659 660 ! Check if climatological information is available 661 IF ( SUM(pclim(:,:)) == 0._wp ) THEN 662 ld_clim = .FALSE. 663 ELSE 664 ld_clim = .TRUE. 665 ENDIF 595 666 596 667 ! Work out the maximum footprint size for the … … 679 750 & ) 680 751 752 IF ( ld_clim ) ALLOCATE( zclim(imaxifp,imaxjfp,isurf) ) 753 681 754 DO jobs = surfdataqc%nsurfup + 1, surfdataqc%nsurfup + isurf 682 755 iobs = jobs - surfdataqc%nsurfup … … 720 793 & igrdip1, igrdjp1, gphif, zgphif ) 721 794 795 IF ( ld_clim ) THEN 796 CALL obs_int_comm_2d( imaxifp, imaxjfp, isurf, kpi, kpj, & 797 & igrdi, igrdj, pclim, zclim ) 798 ENDIF 799 722 800 ! At the end of the day get interpolated means 723 801 IF ( idayend == 0 .AND. ldnightav ) THEN … … 775 853 CALL obs_int_h2d( 1, 1, zweig, zsurftmp(:,:,iobs), zext ) 776 854 855 IF ( ld_clim ) THEN 856 CALL obs_int_h2d( 1, 1, zweig, zclim(:,:,iobs), zclm ) 857 ENDIF 858 859 777 860 ELSE 778 861 … … 788 871 & zweig, zsurftmp(:,:,iobs), zext ) 789 872 873 IF ( ld_clim ) THEN 874 CALL obs_avg_h2d( 1, 1, imaxifp, imaxjfp, & 875 & zweig, zsurftmp(:,:,iobs), zclm ) 876 ENDIF 877 790 878 ENDIF 791 879 … … 797 885 surfdataqc%rmod(jobs,1) = zext(1) 798 886 ENDIF 887 888 IF ( ld_clim ) surfdataqc%rclm(jobs,1) = zclm(1) 799 889 800 890 IF ( zext(1) == obfillflt ) THEN … … 821 911 & ) 822 912 913 IF ( ld_clim ) DEALLOCATE( zclim ) 914 823 915 ! At the end of the day also deallocate night-time mean array 824 916 IF ( idayend == 0 .AND. ldnightav ) THEN
Note: See TracChangeset
for help on using the changeset viewer.