New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
Changeset 11449 for branches/UKMO/dev_r5518_obs_oper_update_addclim/NEMOGCM/NEMO/OPA_SRC/OBS/obs_oper.F90 – NEMO

Ignore:
Timestamp:
2019-08-19T10:54:47+02:00 (5 years ago)
Author:
mattmartin
Message:

Committed first version of changes to output climatology values at obs locations in the fdbk files.

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  
    6262   SUBROUTINE obs_prof_opt( prodatqc, kt, kpi, kpj, kpk, & 
    6363      &                     kit000, kdaystp, kvar,       & 
    64       &                     pvar, pgdept, pgdepw,        & 
    65       &                     pmask,                       &   
     64      &                     pvar, pclim,                 & 
     65      &                     pgdept, pgdepw, pmask,       &   
    6666      &                     plam, pphi,                  & 
    6767      &                     k1dint, k2dint, kdailyavtypes ) 
     
    137137      REAL(KIND=wp), INTENT(IN), DIMENSION(kpi,kpj,kpk) :: & 
    138138         & pvar,   &                 ! Model field for variable 
     139         & pclim,  &                 ! Climatology field for variable          
    139140         & pmask                     ! Land-sea mask for variable 
    140141      REAL(KIND=wp), INTENT(IN), DIMENSION(kpi,kpj) :: & 
     
    172173      REAL(KIND=wp), DIMENSION(kpk) :: & 
    173174         & zobsk,    & 
    174          & zobs2k 
     175         & zobs2k,   & 
     176         & zclm2k          
    175177      REAL(KIND=wp), DIMENSION(2,2,1) :: & 
    176178         & zweig1, & 
     
    178180      REAL(wp), DIMENSION(:,:,:,:), ALLOCATABLE :: & 
    179181         & zmask,  & 
     182         & zclim,  &          
    180183         & zint,   & 
    181184         & zinm,   & 
     
    189192 
    190193      LOGICAL :: ld_dailyav 
     194      LOGICAL :: ld_clim 
    191195 
    192196      !------------------------------------------------------------------------ 
     
    196200      inrc = kt - kit000 + 2 
    197201      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 
    198209 
    199210      ! Daily average types 
     
    262273         & ) 
    263274 
     275      IF ( ld_clim ) ALLOCATE( zclim(2,2,kpk,ipro) ) 
     276 
    264277      DO jobs = prodatqc%nprofup + 1, prodatqc%nprofup + ipro 
    265278         iobs = jobs - prodatqc%nprofup 
     
    286299      CALL obs_int_comm_3d( 2, 2, ipro, kpi, kpj, kpk, igrdi, igrdj, pgdepw, zgdepw )  
    287300 
     301      IF ( ld_clim ) THEN 
     302         CALL obs_int_comm_3d( 2, 2, ipro, kpi, kpj, kpk, igrdi, igrdj, pclim, zclim )             
     303      ENDIF  
     304       
    288305      ! At the end of the day also get interpolated means 
    289306      IF ( ld_dailyav .AND. idayend == 0 ) THEN 
     
    349366                  inum_obs = iend - ista + 1  
    350367                  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                   
    352370                  DO iin=1,2  
    353371                     DO ijn=1,2  
     
    358376                              &     zobs2k, zgdept(iin,ijn,:,iobs), &  
    359377                              &     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 
    360386                        ENDIF  
    361387        
     
    371397                           &    zgdept(iin,ijn,:,iobs), &  
    372398                           &    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                         
    374409                     ENDDO  
    375410                  ENDDO  
     
    385420               iend = prodatqc%npvend(jobs,kvar)  
    386421               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) )                   
    388425               DO iin=1,2   
    389426                  DO ijn=1,2  
     
    394431                           &    zobs2k, zgdept(iin,ijn,:,iobs), &  
    395432                           &    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 
    396440   
    397441                     ENDIF  
     
    408452                         &          zgdept(iin,ijn,:,iobs),         &  
    409453                         &          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    
    410463          
    411464                  ENDDO  
     
    451504                  &              prodatqc%var(kvar)%vmod(iend:iend) )  
    452505 
     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 
    453511                  ! Set QC flag for any observations found below the bottom 
    454512                  ! needed as the check here is more strict than that in obs_prep 
     
    458516  
    459517            DEALLOCATE(interp_corner,iv_indic)  
    460            
     518            IF ( ld_clim ) DEALLOCATE( interp_corner_clim )          
     519              
    461520         ENDIF 
    462521 
     
    475534         & ) 
    476535 
     536      IF ( ld_clim ) DEALLOCATE( zclim ) 
     537       
    477538      ! At the end of the day also get interpolated means 
    478539      IF ( ld_dailyav .AND. idayend == 0 ) THEN 
     
    487548 
    488549   SUBROUTINE obs_surf_opt( surfdataqc, kt, kpi, kpj,            & 
    489       &                     kit000, kdaystp, psurf, psurfmask,   & 
     550      &                     kit000, kdaystp, psurf, pclim, psurfmask,   & 
    490551      &                     k2dint, ldnightav, plamscl, pphiscl, & 
    491552      &                     lindegrees ) 
     
    541602      REAL(wp), INTENT(IN), DIMENSION(kpi,kpj) :: & 
    542603         & psurf,  &                   ! Model surface field 
     604         & pclim,  &                   ! Climatological surface field          
    543605         & psurfmask                   ! Land-sea mask 
    544606      LOGICAL, INTENT(IN) :: ldnightav ! Logical for averaging night-time data 
     
    569631      REAL(wp) :: zlam 
    570632      REAL(wp) :: zphi 
    571       REAL(wp), DIMENSION(1) :: zext, zobsmask 
     633      REAL(wp), DIMENSION(1) :: zext, zobsmask, zclm 
    572634      REAL(wp) :: zdaystp 
    573635      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: & 
     
    577639         & zsurfm, & 
    578640         & zsurftmp, & 
     641         & zclim,  & 
    579642         & zglam,  & 
    580643         & zgphi,  & 
     
    586649         & zouttmp, & 
    587650         & zmeanday    ! to compute model sst in region of 24h daylight (pole) 
    588  
     651          
     652      LOGICAL :: ld_clim ! T => climatological data is available 
    589653      !------------------------------------------------------------------------ 
    590654      ! Local initialization  
     
    593657      inrc = kt - kit000 + 2 
    594658      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 
    595666 
    596667      ! Work out the maximum footprint size for the  
     
    679750         & ) 
    680751 
     752      IF ( ld_clim ) ALLOCATE( zclim(imaxifp,imaxjfp,isurf) ) 
     753 
    681754      DO jobs = surfdataqc%nsurfup + 1, surfdataqc%nsurfup + isurf 
    682755         iobs = jobs - surfdataqc%nsurfup 
     
    720793         &                  igrdip1, igrdjp1, gphif, zgphif ) 
    721794 
     795      IF ( ld_clim ) THEN  
     796         CALL obs_int_comm_2d( imaxifp, imaxjfp, isurf, kpi, kpj, & 
     797            &                  igrdi, igrdj, pclim, zclim ) 
     798      ENDIF 
     799 
    722800      ! At the end of the day get interpolated means 
    723801      IF ( idayend == 0 .AND. ldnightav ) THEN 
     
    775853            CALL obs_int_h2d( 1, 1, zweig, zsurftmp(:,:,iobs), zext ) 
    776854 
     855            IF ( ld_clim ) THEN   
     856               CALL obs_int_h2d( 1, 1, zweig, zclim(:,:,iobs), zclm ) 
     857            ENDIF 
     858 
     859 
    777860         ELSE 
    778861 
     
    788871               &              zweig, zsurftmp(:,:,iobs),  zext ) 
    789872 
     873            IF ( ld_clim ) THEN   
     874               CALL obs_avg_h2d( 1, 1, imaxifp, imaxjfp, & 
     875                  &              zweig, zsurftmp(:,:,iobs),  zclm ) 
     876            ENDIF 
     877 
    790878         ENDIF 
    791879 
     
    797885            surfdataqc%rmod(jobs,1) = zext(1) 
    798886         ENDIF 
     887          
     888         IF ( ld_clim ) surfdataqc%rclm(jobs,1) = zclm(1) 
    799889          
    800890         IF ( zext(1) == obfillflt ) THEN 
     
    821911         & ) 
    822912 
     913      IF ( ld_clim ) DEALLOCATE( zclim ) 
     914 
    823915      ! At the end of the day also deallocate night-time mean array 
    824916      IF ( idayend == 0 .AND. ldnightav ) THEN 
Note: See TracChangeset for help on using the changeset viewer.