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 11455 for branches/UKMO/dev_r5518_obs_oper_update_addclim/NEMOGCM/NEMO/OPA_SRC/OBS/obs_oper.F90 – NEMO

Ignore:
Timestamp:
2019-08-19T17:36:23+02:00 (5 years ago)
Author:
mattmartin
Message:

Commit version which compiles and runs. Not fully tested that it is producing the correct answer yet though.

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/UKMO/dev_r5518_obs_oper_update_addclim/NEMOGCM/NEMO/OPA_SRC/OBS/obs_oper.F90

    r11449 r11455  
    190190      REAL(KIND=wp), DIMENSION(1) :: zmsk 
    191191      REAL(KIND=wp), DIMENSION(:,:,:), ALLOCATABLE :: interp_corner 
    192  
     192      REAL(KIND=wp), DIMENSION(:,:,:), ALLOCATABLE :: interp_corner_clim 
     193       
    193194      LOGICAL :: ld_dailyav 
    194       LOGICAL :: ld_clim 
    195195 
    196196      !------------------------------------------------------------------------ 
     
    200200      inrc = kt - kit000 + 2 
    201201      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 
    209202 
    210203      ! Daily average types 
     
    273266         & ) 
    274267 
    275       IF ( ld_clim ) ALLOCATE( zclim(2,2,kpk,ipro) ) 
     268      IF ( prodatqc%lclim ) ALLOCATE( zclim(2,2,kpk,ipro) ) 
    276269 
    277270      DO jobs = prodatqc%nprofup + 1, prodatqc%nprofup + ipro 
     
    299292      CALL obs_int_comm_3d( 2, 2, ipro, kpi, kpj, kpk, igrdi, igrdj, pgdepw, zgdepw )  
    300293 
    301       IF ( ld_clim ) THEN 
     294      IF ( prodatqc%lclim ) THEN 
    302295         CALL obs_int_comm_3d( 2, 2, ipro, kpi, kpj, kpk, igrdi, igrdj, pclim, zclim )             
    303296      ENDIF  
     
    366359                  inum_obs = iend - ista + 1  
    367360                  ALLOCATE(interp_corner(2,2,inum_obs),iv_indic(inum_obs))  
    368                   IF ( ld_clim ) ALLOCATE( interp_corner_clim(2,2,inum_obs) ) 
     361                  IF ( prodatqc%lclim ) ALLOCATE( interp_corner_clim(2,2,inum_obs) ) 
    369362                   
    370363                  DO iin=1,2  
     
    377370                              &     zmask(iin,ijn,:,iobs))  
    378371 
    379                            IF ( ld_clim ) THEN 
     372                           IF ( prodatqc%lclim ) THEN 
    380373                              CALL obs_int_z1d_spl( kpk, &  
    381374                                 &     zclim(iin,ijn,:,iobs), &  
     
    398391                           &    zmask(iin,ijn,:,iobs))  
    399392 
    400                         IF ( ld_clim ) THEN 
     393                        IF ( prodatqc%lclim ) THEN 
    401394                           CALL obs_int_z1d(kpk, iv_indic, k1dint, inum_obs, &  
    402395                              &    prodatqc%var(kvar)%vdep(ista:iend), &  
     
    422415               ALLOCATE( interp_corner(2,2,inum_obs),      & 
    423416                  &      iv_indic(inum_obs) )  
    424                IF ( ld_clim ) ALLOCATE( interp_corner_clim(2,2,inum_obs) )                   
     417               IF ( prodatqc%lclim ) ALLOCATE( interp_corner_clim(2,2,inum_obs) )                   
    425418               DO iin=1,2   
    426419                  DO ijn=1,2  
     
    432425                           &    zmask(iin,ijn,:,iobs))  
    433426 
    434                         IF ( ld_clim ) THEN 
     427                        IF ( prodatqc%lclim ) THEN 
    435428                           CALL obs_int_z1d_spl( kpk, &  
    436429                              &    zclim(iin,ijn,:,iobs),&  
     
    453446                         &          zmask(iin,ijn,:,iobs) )       
    454447 
    455                      IF ( ld_clim ) THEN 
     448                     IF ( prodatqc%lclim ) THEN 
    456449                        CALL obs_int_z1d(kpk, iv_indic, k1dint, inum_obs,     &  
    457450                            &          prodatqc%var(kvar)%vdep(ista:iend),     &  
     
    504497                  &              prodatqc%var(kvar)%vmod(iend:iend) )  
    505498 
    506                IF ( ld_clim ) THEN 
     499               IF ( prodatqc%lclim ) THEN 
    507500                  CALL obs_int_h2d( 1, 1, zweig, interp_corner_clim(:,:,ikn), &  
    508501                     &              prodatqc%var(kvar)%vclm(iend:iend) ) 
     
    516509  
    517510            DEALLOCATE(interp_corner,iv_indic)  
    518             IF ( ld_clim ) DEALLOCATE( interp_corner_clim )          
     511            IF ( prodatqc%lclim ) DEALLOCATE( interp_corner_clim )          
    519512              
    520513         ENDIF 
     
    534527         & ) 
    535528 
    536       IF ( ld_clim ) DEALLOCATE( zclim ) 
     529      IF ( prodatqc%lclim ) DEALLOCATE( zclim ) 
    537530       
    538531      ! At the end of the day also get interpolated means 
     
    650643         & zmeanday    ! to compute model sst in region of 24h daylight (pole) 
    651644          
    652       LOGICAL :: ld_clim ! T => climatological data is available 
    653645      !------------------------------------------------------------------------ 
    654646      ! Local initialization  
     
    657649      inrc = kt - kit000 + 2 
    658650      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 
    666651 
    667652      ! Work out the maximum footprint size for the  
     
    750735         & ) 
    751736 
    752       IF ( ld_clim ) ALLOCATE( zclim(imaxifp,imaxjfp,isurf) ) 
     737      IF ( surfdataqc%lclim ) ALLOCATE( zclim(imaxifp,imaxjfp,isurf) ) 
    753738 
    754739      DO jobs = surfdataqc%nsurfup + 1, surfdataqc%nsurfup + isurf 
     
    793778         &                  igrdip1, igrdjp1, gphif, zgphif ) 
    794779 
    795       IF ( ld_clim ) THEN  
     780      IF ( surfdataqc%lclim ) THEN  
    796781         CALL obs_int_comm_2d( imaxifp, imaxjfp, isurf, kpi, kpj, & 
    797782            &                  igrdi, igrdj, pclim, zclim ) 
     
    853838            CALL obs_int_h2d( 1, 1, zweig, zsurftmp(:,:,iobs), zext ) 
    854839 
    855             IF ( ld_clim ) THEN   
     840            IF ( surfdataqc%lclim ) THEN   
    856841               CALL obs_int_h2d( 1, 1, zweig, zclim(:,:,iobs), zclm ) 
     842               IF (lwp) THEN 
     843                  WRITE(numout,*)'zclim: ', iobs, zclim(:,:,iobs), zclm 
     844               ENDIF 
    857845            ENDIF 
    858846 
     
    871859               &              zweig, zsurftmp(:,:,iobs),  zext ) 
    872860 
    873             IF ( ld_clim ) THEN   
     861            IF ( surfdataqc%lclim ) THEN   
    874862               CALL obs_avg_h2d( 1, 1, imaxifp, imaxjfp, & 
    875863                  &              zweig, zsurftmp(:,:,iobs),  zclm ) 
     
    886874         ENDIF 
    887875          
    888          IF ( ld_clim ) surfdataqc%rclm(jobs,1) = zclm(1) 
     876         IF ( surfdataqc%lclim ) surfdataqc%rclm(jobs,1) = zclm(1) 
    889877          
    890878         IF ( zext(1) == obfillflt ) THEN 
     
    911899         & ) 
    912900 
    913       IF ( ld_clim ) DEALLOCATE( zclim ) 
     901      IF ( surfdataqc%lclim ) DEALLOCATE( zclim ) 
    914902 
    915903      ! At the end of the day also deallocate night-time mean array 
Note: See TracChangeset for help on using the changeset viewer.