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 15246 for NEMO/branches/UKMO/NEMO_4.0.4_generic_obs – NEMO

Ignore:
Timestamp:
2021-09-10T13:08:02+02:00 (3 years ago)
Author:
dford
Message:

Add option to output T&S climatology.

Location:
NEMO/branches/UKMO/NEMO_4.0.4_generic_obs
Files:
5 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/UKMO/NEMO_4.0.4_generic_obs/cfgs/SHARED/namelist_ref

    r15240 r15246  
    12921292   rn_mdtcorr             = 1.61           !    MDT correction 
    12931293   rn_mdtcutoff           = 65.0           !    MDT cutoff for computed correction 
    1294 !!! NOT YET IMPLEMENTED: 
    1295 !!!   OUTPUT CLIMATOLOGY 
     1294                                           ! Relevant if 'POTM', 'PSAL', 'SST', or 'SSS' in cn_obstypes: 
     1295   ln_output_clim         = .false.        !    Logical switch to output climatological temperature/salinity (if ln_tradmp) 
    12961296/ 
    12971297!----------------------------------------------------------------------- 
  • NEMO/branches/UKMO/NEMO_4.0.4_generic_obs/src/OCE/OBS/diaobs.F90

    r15240 r15246  
    8484      !!---------------------------------------------------------------------- 
    8585#if defined key_si3 
    86    USE ice, ONLY : &     ! Sea ice variables 
    87       & hm_s             ! Snow depth for freeboard conversion 
     86      USE ice, ONLY : &     ! Sea ice variables 
     87         & hm_s             ! Snow depth for freeboard conversion 
    8888#elif defined key_cice 
    89    USE sbc_oce, ONLY : & ! Sea ice variables 
    90       & thick_s          ! Snow depth for freeboard conversion 
     89      USE sbc_oce, ONLY : & ! Sea ice variables 
     90         & thick_s          ! Snow depth for freeboard conversion 
    9191#endif 
     92      USE obs_fbm, ONLY : & 
     93         & fbrmdi           ! Real missing data indicator 
    9294 
    9395      IMPLICIT NONE 
     
    232234               ENDIF 
    233235               ! 
     236               IF( sobsgroups(jgroup)%loutput_clim ) THEN 
     237                  sobsgroups(jgroup)%sprofdata%caddvars(sobsgroups(jgroup)%nadd_clm)  = 'CLM' 
     238                  DO jvar = 1, sobsgroups(jgroup)%nobstypes 
     239                     sobsgroups(jgroup)%sprofdata%var(jvar)%vadd(:,sobsgroups(jgroup)%nadd_clm) = fbrmdi 
     240                     sobsgroups(jgroup)%sprofdata%caddlong(sobsgroups(jgroup)%nadd_clm,jvar) = 'Climatology' 
     241                     sobsgroups(jgroup)%sprofdata%caddunit(sobsgroups(jgroup)%nadd_clm,jvar) = sobsgroups(jgroup)%sprofdata%cunit(jvar) 
     242                  END DO 
     243               ENDIF 
     244               ! 
    234245               CALL obs_pre_prof( sobsgroups(jgroup)%sprofdata,     & 
    235246                  &               sobsgroups(jgroup)%sprofdataqc,   & 
     
    284295                     sobsgroups(jgroup)%ssurfdata%caddlong(sobsgroups(jgroup)%nadd_fbd,jvar) = 'Freeboard' 
    285296                     sobsgroups(jgroup)%ssurfdata%caddunit(sobsgroups(jgroup)%nadd_fbd,jvar) = 'Metres' 
     297                  END DO 
     298               ENDIF 
     299               ! 
     300               IF( sobsgroups(jgroup)%loutput_clim ) THEN 
     301                  sobsgroups(jgroup)%ssurfdata%caddvars(sobsgroups(jgroup)%nadd_clm)  = 'CLM' 
     302                  DO jvar = 1, sobsgroups(jgroup)%nobstypes 
     303                     sobsgroups(jgroup)%ssurfdata%radd(:,:,jvar) = fbrmdi 
     304                     sobsgroups(jgroup)%ssurfdata%caddlong(sobsgroups(jgroup)%nadd_clm,jvar) = 'Climatology' 
     305                     sobsgroups(jgroup)%ssurfdata%caddunit(sobsgroups(jgroup)%nadd_clm,jvar) = sobsgroups(jgroup)%ssurfdata%cunit(jvar) 
    286306                  END DO 
    287307               ENDIF 
     
    368388      USE sbc_oce, ONLY : fr_i, thick_i       ! CICE Ice model variables 
    369389#endif 
     390      USE tradmp,  ONLY : tclim, sclim        ! T&S climatology 
     391      USE obs_fbm, ONLY : fbrmdi              ! Real missing data indicator 
    370392 
    371393      IMPLICIT NONE 
     
    382404      LOGICAL :: lstp0             ! Flag special treatment on zeroth time step 
    383405      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: & 
    384          & zprofvar                ! Model values for variables in a prof ob 
     406         & zprofvar, &             ! Model values for variables in a prof ob 
     407         & zprofclim               ! Climatology values for variables in a prof ob 
    385408      REAL(wp), DIMENSION(:,:), ALLOCATABLE :: & 
    386          & zsurfvar                ! Model values for variables in a surf ob 
     409         & zsurfvar, &             ! Model values for variables in a surf ob 
     410         & zsurfclim               ! Climatology values for variables in a surf ob 
    387411 
    388412      !----------------------------------------------------------------------- 
     
    402426      !----------------------------------------------------------------------- 
    403427 
     428      ALLOCATE( zprofvar(jpi,jpj,jpk),  & 
     429         &      zprofclim(jpi,jpj,jpk), & 
     430         &      zsurfvar(jpi,jpj),      & 
     431         &      zsurfclim(jpi,jpj) ) 
     432 
    404433      DO jgroup = 1, nn_obsgroups 
    405434         IF (sobsgroups(jgroup)%lenabled) THEN 
     
    407436            IF (sobsgroups(jgroup)%lprof) THEN 
    408437 
    409                ALLOCATE( zprofvar(jpi,jpj,jpk) ) 
     438               zprofclim(:,:,:) = fbrmdi 
    410439 
    411440               DO jvar = 1, sobsgroups(jgroup)%nobstypes 
     
    414443                  CASE('POTM') 
    415444                     zprofvar(:,:,:) = tsn(:,:,:,jp_tem) 
     445                     IF (sobsgroups(jgroup)%loutput_clim) THEN 
     446                        zprofclim(:,:,:) = tclim(:,:,:) 
     447                     ENDIF 
    416448                  CASE('PSAL') 
    417449                     zprofvar(:,:,:) = tsn(:,:,:,jp_sal) 
     450                     IF (sobsgroups(jgroup)%loutput_clim) THEN 
     451                        zprofclim(:,:,:) = sclim(:,:,:) 
     452                     ENDIF 
    418453                  CASE('UVEL') 
    419454                     zprofvar(:,:,:) = un(:,:,:) 
     
    428463                     &               nit000, idaystp, jvar,                & 
    429464                     &               zprofvar,                             & 
     465                     &               sobsgroups(jgroup)%loutput_clim,      & 
     466                     &               sobsgroups(jgroup)%nadd_clm,          & 
     467                     &               zprofclim,                            & 
    430468                     &               gdept_n,                              & 
    431469                     &               gdepw_n,                              &  
     
    439477               END DO 
    440478 
    441                DEALLOCATE( zprofvar ) 
    442  
    443479            ELSEIF (sobsgroups(jgroup)%lsurf) THEN 
    444480 
    445                ALLOCATE( zsurfvar(jpi,jpj) ) 
     481               zsurfclim(:,:) = fbrmdi 
    446482 
    447483               DO jvar = 1, sobsgroups(jgroup)%nobstypes 
     
    451487                  CASE('SST') 
    452488                     zsurfvar(:,:) = tsn(:,:,1,jp_tem) 
     489                     IF (sobsgroups(jgroup)%loutput_clim) THEN 
     490                        zsurfclim(:,:) = tclim(:,:,1) 
     491                     ENDIF 
    453492                  CASE('SLA') 
    454493                     zsurfvar(:,:) = sshn(:,:) 
    455494                  CASE('SSS') 
    456495                     zsurfvar(:,:) = tsn(:,:,1,jp_sal) 
     496                     IF (sobsgroups(jgroup)%loutput_clim) THEN 
     497                        zsurfclim(:,:) = sclim(:,:,1) 
     498                     ENDIF 
    457499                  CASE('UVEL') 
    458500                     zsurfvar(:,:) = un(:,:,1) 
     
    515557                        &               nit000, idaystp,                      & 
    516558                        &               jvar, zsurfvar,                       & 
     559                        &               sobsgroups(jgroup)%loutput_clim,      & 
     560                        &               sobsgroups(jgroup)%nadd_clm,          & 
     561                        &               zsurfclim,                            & 
    517562                        &               sobsgroups(jgroup)%rmask(:,:,1,jvar), & 
    518563                        &               sobsgroups(jgroup)%n2dint,            & 
     
    531576               END DO 
    532577 
    533                DEALLOCATE( zsurfvar ) 
    534  
    535578            ENDIF 
    536579 
    537580         ENDIF 
    538581      END DO 
     582 
     583      DEALLOCATE( zprofvar, zprofclim, & 
     584         &        zsurfvar, zsurfclim ) 
    539585 
    540586      IF( ln_timing )   CALL timing_stop('dia_obs') 
  • NEMO/branches/UKMO/NEMO_4.0.4_generic_obs/src/OCE/OBS/obs_group_def.F90

    r15240 r15246  
    1717   USE obs_surf_def       ! Surface data definitions 
    1818   USE obs_profiles_def   ! Profile data definitions 
     19   USE tradmp, ONLY : &   ! Damping to climatology 
     20      & ln_tradmp 
    1921 
    2022   IMPLICIT NONE 
     
    4143   CHARACTER(LEN=8), PUBLIC :: cobsname_sla  = 'SLA'  ! Expected variable name for sea level anomaly 
    4244   CHARACTER(LEN=8), PUBLIC :: cobsname_fbd  = 'FBD'  ! Expected variable name for sea ice freeboard 
     45   CHARACTER(LEN=8), PUBLIC :: cobsname_t3d  = 'POTM' ! Expected variable name for 3D temperature 
     46   CHARACTER(LEN=8), PUBLIC :: cobsname_s3d  = 'PSAL' ! Expected variable name for 3D salinity 
     47   CHARACTER(LEN=8), PUBLIC :: cobsname_t2d  = 'SST'  ! Expected variable name for 2D temperature 
     48   CHARACTER(LEN=8), PUBLIC :: cobsname_s2d  = 'SSS'  ! Expected variable name for 2D salinity 
    4349 
    4450   !! * Type definition for observation groups 
     
    6874      INTEGER  :: nadd_fbd           !: Index of additional variable representing ice freeboard 
    6975      INTEGER  :: next_snow          !: Index of extra variable representing snow thickness 
     76      INTEGER  :: nadd_clm           !: Index of additional variable representing climatology 
    7077      ! 
    7178      LOGICAL  :: lenabled           !: Logical switch for group being processed and not ignored 
     
    8390      LOGICAL  :: lnight             !: Logical switch for calculating night-time average 
    8491      LOGICAL  :: ltime_mean_bkg     !: Logical switch for applying time mean of background (e.g. to remove tidal signal) 
     92      LOGICAL  :: loutput_clim       !: Logical switch for outputting climatological temperature/salinity 
    8593      LOGICAL  :: lfp_indegs         !: Logical: T=> averaging footprint is in degrees, F=> in metres 
    8694      ! 
     
    196204      LOGICAL                                    :: ln_night 
    197205      LOGICAL                                    :: ln_time_mean_bkg 
     206      LOGICAL                                    :: ln_output_clim 
    198207      LOGICAL                                    :: ln_fp_indegs 
    199208      REAL(wp)                                   :: rn_avglamscl 
     
    210219         &                cn_obsbiasfile_varname, ln_night, ln_time_mean_bkg,   & 
    211220         &                rn_time_mean_period, ln_altbias, cn_altbiasfile,      & 
    212          &                nn_msshc, rn_mdtcorr, rn_mdtcutoff, ln_all_at_all 
     221         &                nn_msshc, rn_mdtcorr, rn_mdtcutoff, ln_all_at_all,    & 
     222         &                ln_output_clim 
    213223      !!---------------------------------------------------------------------- 
    214224 
     
    292302                  sdobsgroup%nadd_fbd  = sdobsgroup%naddvars 
    293303                  sdobsgroup%next_snow = sdobsgroup%nextvars 
     304               ENDIF 
     305               ! 
     306               IF ( ln_output_clim .AND. ln_tradmp .AND. (.NOT. sdobsgroup%loutput_clim) ) THEN 
     307                  IF ( (TRIM(sdobsgroup%cobstypes(itype)) == cobsname_t3d) .OR. & 
     308                     & (TRIM(sdobsgroup%cobstypes(itype)) == cobsname_s3d) .OR. & 
     309                     & (TRIM(sdobsgroup%cobstypes(itype)) == cobsname_t2d) .OR. & 
     310                     & (TRIM(sdobsgroup%cobstypes(itype)) == cobsname_s2d) ) THEN 
     311                     sdobsgroup%loutput_clim = .TRUE. 
     312                     sdobsgroup%naddvars = sdobsgroup%naddvars + 1 
     313                     sdobsgroup%nadd_clm = sdobsgroup%naddvars 
     314                  ENDIF 
    294315               ENDIF 
    295316               ! 
     
    438459            WRITE(numout,*) '             MDT  correction                                      rn_mdtcorr = ', sdobsgroup%rmdtcorr 
    439460            WRITE(numout,*) '             MDT cutoff for computed correction                 rn_mdtcutoff = ', sdobsgroup%rmdtcutoff 
     461            WRITE(numout,*) '          Settings only for temperature/salinity data' 
     462            WRITE(numout,*) '             Logical switch for outputting climatology        ln_output_clim = ', sdobsgroup%loutput_clim 
    440463         ENDIF 
    441464 
  • NEMO/branches/UKMO/NEMO_4.0.4_generic_obs/src/OCE/OBS/obs_oper.F90

    r15240 r15246  
    4242   SUBROUTINE obs_prof_opt( prodatqc, kt, kpi, kpj, kpk, & 
    4343      &                     kit000, kdaystp, kvar,       & 
    44       &                     pvar, pgdept, pgdepw,        & 
     44      &                     pvar,                        & 
     45      &                     ldclim, kclim, pclim,        & 
     46      &                     pgdept, pgdepw,              & 
    4547      &                     pmask,                       &   
    4648      &                     plam, pphi,                  & 
     
    106108      INTEGER       , INTENT(in   ) ::   kdaystp         ! Number of time steps per day 
    107109      INTEGER       , INTENT(in   ) ::   kvar            ! Index of variable in prodatqc 
     110      INTEGER       , INTENT(in   ) ::   kclim           ! Index of climatology in prodatqc 
     111      LOGICAL       , INTENT(in   ) ::   ldclim          ! Switch to interpolate climatology 
    108112      REAL(KIND=wp) , INTENT(in   ), DIMENSION(kpi,kpj,kpk) ::   pvar             ! Model field 
     113      REAL(KIND=wp) , INTENT(in   ), DIMENSION(kpi,kpj,kpk) ::   pclim            ! Climatology field 
    109114      REAL(KIND=wp) , INTENT(in   ), DIMENSION(kpi,kpj,kpk) ::   pmask            ! Land-sea mask 
    110115      REAL(KIND=wp) , INTENT(in   ), DIMENSION(kpi,kpj)     ::   plam             ! Model longitude 
     
    138143      REAL(KIND=wp), DIMENSION(kpk) :: & 
    139144         & zobsk,  & 
    140          & zobs2k 
     145         & zobs2k, & 
     146         & zclm2k 
    141147      REAL(KIND=wp), DIMENSION(2,2,1) :: & 
    142148         & zweig1, & 
     
    144150      REAL(wp), DIMENSION(:,:,:,:), ALLOCATABLE :: & 
    145151         & zmask,  & 
     152         & zclim,  & 
    146153         & zint,   & 
    147154         & zinm,   & 
     
    153160      REAL(KIND=wp), DIMENSION(1) :: zmsk 
    154161      REAL(KIND=wp), DIMENSION(:,:,:), ALLOCATABLE :: interp_corner 
     162      REAL(KIND=wp), DIMENSION(:,:,:), ALLOCATABLE :: interp_corner_clim 
    155163 
    156164      LOGICAL :: ld_dailyav 
     
    228236         & ) 
    229237 
     238      IF ( ldclim ) THEN 
     239         ALLOCATE( zclim(2,2,kpk,ipro) ) 
     240      ENDIF 
     241 
    230242      DO jobs = prodatqc%nprofup + 1, prodatqc%nprofup + ipro 
    231243         iobs = jobs - prodatqc%nprofup 
     
    251263      CALL obs_int_comm_3d( 2, 2, ipro, kpi, kpj, kpk, igrdi, igrdj, pgdept, zgdept )  
    252264      CALL obs_int_comm_3d( 2, 2, ipro, kpi, kpj, kpk, igrdi, igrdj, pgdepw, zgdepw )  
     265 
     266      IF ( ldclim ) THEN 
     267         CALL obs_int_comm_3d( 2, 2, ipro, kpi, kpj, kpk, igrdi, igrdj, pclim, zclim ) 
     268      ENDIF 
    253269 
    254270      ! At the end of the day also get interpolated means 
     
    315331                  inum_obs = iend - ista + 1  
    316332                  ALLOCATE(interp_corner(2,2,inum_obs),iv_indic(inum_obs))  
     333                  IF ( ldclim ) THEN 
     334                     ALLOCATE( interp_corner_clim(2,2,inum_obs) ) 
     335                  ENDIF 
    317336 
    318337                  DO iin=1,2  
     
    324343                              &     zobs2k, zgdept(iin,ijn,:,iobs), &  
    325344                              &     zmask(iin,ijn,:,iobs))  
     345 
     346                           IF ( ldclim ) THEN 
     347                              CALL obs_int_z1d_spl( kpk, & 
     348                                 &     zclim(iin,ijn,:,iobs), & 
     349                                 &     zclm2k, zgdept(iin,ijn,:,iobs), & 
     350                                 &     zmask(iin,ijn,:,iobs)) 
     351                           ENDIF 
    326352                        ENDIF  
    327353        
     
    337363                           &    zgdept(iin,ijn,:,iobs), &  
    338364                           &    zmask(iin,ijn,:,iobs))  
     365 
     366                        IF ( ldclim ) THEN 
     367                           CALL obs_int_z1d(kpk, iv_indic, k1dint, inum_obs, & 
     368                              &    prodatqc%var(kvar)%vdep(ista:iend), & 
     369                              &    zclim(iin,ijn,:,iobs), & 
     370                              &    zclm2k, interp_corner_clim(iin,ijn,:), & 
     371                              &    zgdept(iin,ijn,:,iobs), & 
     372                              &    zmask(iin,ijn,:,iobs)) 
     373                        ENDIF 
    339374        
    340375                     ENDDO  
     
    352387               inum_obs = iend - ista + 1  
    353388               ALLOCATE(interp_corner(2,2,inum_obs), iv_indic(inum_obs))  
     389               IF ( ldclim ) THEN 
     390                  ALLOCATE( interp_corner_clim(2,2,inum_obs) ) 
     391               ENDIF 
    354392               DO iin=1,2   
    355393                  DO ijn=1,2  
     
    360398                           &    zobs2k, zgdept(iin,ijn,:,iobs), &  
    361399                           &    zmask(iin,ijn,:,iobs))  
    362    
     400 
     401                        IF ( ldclim ) THEN 
     402                           CALL obs_int_z1d_spl( kpk, & 
     403                              &    zclim(iin,ijn,:,iobs),& 
     404                              &    zclm2k, zgdept(iin,ijn,:,iobs), & 
     405                              &    zmask(iin,ijn,:,iobs)) 
     406                        ENDIF 
    363407                     ENDIF  
    364408        
     
    374418                         &          zgdept(iin,ijn,:,iobs),         &  
    375419                         &          zmask(iin,ijn,:,iobs) )       
     420 
     421                     IF ( ldclim ) THEN 
     422                        CALL obs_int_z1d(kpk, iv_indic, k1dint, inum_obs,     & 
     423                            &          prodatqc%var(kvar)%vdep(ista:iend),     & 
     424                            &          zclim(iin,ijn,:,iobs),            & 
     425                            &          zclm2k,interp_corner_clim(iin,ijn,:), & 
     426                            &          zgdept(iin,ijn,:,iobs),         & 
     427                            &          zmask(iin,ijn,:,iobs) ) 
     428                     ENDIF 
    376429          
    377430                  ENDDO  
     
    417470                  &              prodatqc%var(kvar)%vmod(iend:iend) )  
    418471 
    419                   ! Set QC flag for any observations found below the bottom 
    420                   ! needed as the check here is more strict than that in obs_prep 
     472               IF ( ldclim ) THEN 
     473                  CALL obs_int_h2d( 1, 1, zweig, interp_corner_clim(:,:,ikn), & 
     474                     &              prodatqc%var(kvar)%vadd(iend:iend,kclim) ) 
     475               ENDIF 
     476 
     477               ! Set QC flag for any observations found below the bottom 
     478               ! needed as the check here is more strict than that in obs_prep 
    421479               IF (sum(zweig) == 0.0_wp) prodatqc%var(kvar)%nvqc(iend:iend)=4 
    422480  
     
    424482  
    425483            DEALLOCATE(interp_corner,iv_indic)  
     484            IF ( ldclim ) THEN 
     485               DEALLOCATE( interp_corner_clim ) 
     486            ENDIF 
    426487           
    427488         ENDIF 
     
    441502         & ) 
    442503 
     504      IF ( ldclim ) THEN 
     505         DEALLOCATE( zclim ) 
     506      ENDIF 
     507 
    443508      ! At the end of the day also get interpolated means 
    444509      IF ( ld_dailyav .AND. idayend == 0 ) THEN 
     
    453518 
    454519   SUBROUTINE obs_surf_opt( surfdataqc, kt, kpi, kpj,                & 
    455       &                     kit000, kdaystp, kvar, psurf, psurfmask, & 
     520      &                     kit000, kdaystp, kvar, psurf,            & 
     521      &                     ldclim, kclim, pclim, psurfmask,         & 
    456522      &                     k2dint, ldnightav, plamscl, pphiscl,     & 
    457523      &                     lindegrees, ldtime_mean, kmeanstp,       & 
     
    503569      INTEGER, INTENT(IN) :: kdaystp   ! Number of time steps per day 
    504570      INTEGER, INTENT(IN) :: kvar      ! Index of variable in surfdataqc   
     571      INTEGER, INTENT(IN) :: kclim     ! Index of climatology in surfdataqc 
    505572      INTEGER, INTENT(IN) :: k2dint    ! Horizontal interpolation type (see header) 
     573      LOGICAL, INTENT(IN) :: ldclim    ! Switch to interpolate climatology 
    506574      REAL(wp), INTENT(IN), DIMENSION(kpi,kpj) :: & 
    507575         & psurf,  &                   ! Model surface field 
     576         & pclim,  &                   ! Climatology surface field 
    508577         & psurfmask                   ! Land-sea mask 
    509578      LOGICAL, INTENT(IN) :: ldnightav ! Logical for averaging night-time data 
     
    546615      REAL(wp) :: zlam 
    547616      REAL(wp) :: zphi 
    548       REAL(wp), DIMENSION(1) :: zext, zobsmask 
     617      REAL(wp), DIMENSION(1) :: zext, zobsmask, zclm 
    549618      REAL(wp) :: zdaystp 
    550619      REAL(wp) :: zmeanstp 
     
    555624         & zsurfm, & 
    556625         & zsurftmp, & 
     626         & zclim,  & 
    557627         & zglam,  & 
    558628         & zgphi,  & 
     
    696766         & igrdjp1(imaxifp+1,imaxjfp+1,isurf) & 
    697767         & ) 
     768 
     769      IF ( ldclim ) THEN 
     770         ALLOCATE( zclim(imaxifp,imaxjfp,isurf) ) 
     771      ENDIF 
    698772 
    699773      DO jobs = surfdataqc%nsurfup + 1, surfdataqc%nsurfup + isurf 
     
    745819      ENDIF 
    746820 
    747 ! IF ( k2dint > 4 ) THEN   AROUND THIS? - SEEMS TO BE IN V3.6 
    748821      CALL obs_int_comm_2d( imaxifp+1, imaxjfp+1, isurf, kpi, kpj, & 
    749822         &                  igrdip1, igrdjp1, glamf, zglamf ) 
    750823      CALL obs_int_comm_2d( imaxifp+1, imaxjfp+1, isurf, kpi, kpj, & 
    751824         &                  igrdip1, igrdjp1, gphif, zgphif ) 
     825       
     826      IF ( ldclim ) THEN 
     827         CALL obs_int_comm_2d( imaxifp, imaxjfp, isurf, kpi, kpj, & 
     828            &                  igrdi, igrdj, pclim, zclim ) 
     829      ENDIF 
    752830 
    753831      ! At the end of the day get interpolated means 
     
    808886               CALL obs_int_h2d( 1, 1, zweig, zsurftmp(:,:,iobs), zext ) 
    809887 
     888               IF ( ldclim ) THEN 
     889                  CALL obs_int_h2d( 1, 1, zweig, zclim(:,:,iobs), zclm ) 
     890               ENDIF 
     891 
    810892            ELSE 
    811893 
     
    820902               CALL obs_avg_h2d( 1, 1, imaxifp, imaxjfp, & 
    821903                  &              zweig, zsurftmp(:,:,iobs),  zext ) 
     904 
     905               IF ( ldclim ) THEN 
     906                  CALL obs_avg_h2d( 1, 1, imaxifp, imaxjfp, & 
     907                     &              zweig, zclim(:,:,iobs),  zclm ) 
     908               ENDIF 
    822909 
    823910            ENDIF 
     
    846933            ENDIF 
    847934 
     935            IF ( ldclim ) THEN 
     936               surfdataqc%radd(jobs,kclim,kvar) = zclm(1) 
     937            ENDIF 
     938 
    848939            IF ( zext(1) == obfillflt ) THEN 
    849940               ! If the observation value is a fill value, set QC flag to bad 
     
    870961         & igrdjp1 & 
    871962         & ) 
     963             
     964      IF ( ldclim ) THEN 
     965         DEALLOCATE( zclim ) 
     966      ENDIF 
    872967 
    873968      ! At the end of the day also deallocate time mean array 
  • NEMO/branches/UKMO/NEMO_4.0.4_generic_obs/src/OCE/TRA/tradmp.F90

    r14075 r15246  
    5050   ! 
    5151   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   resto    !: restoring coeff. on T and S (s-1) 
     52   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   tclim    !: temperature climatology on each time step(Celcius) 
     53   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   sclim    !: salinity climatology on each time step (psu) 
    5254 
    5355   !! * Substitutions 
     
    6466      !!                ***  FUNCTION tra_dmp_alloc  *** 
    6567      !!---------------------------------------------------------------------- 
    66       ALLOCATE( resto(jpi,jpj,jpk), STAT= tra_dmp_alloc ) 
     68      ALLOCATE( resto(jpi,jpj,jpk), & 
     69         &      tclim(jpi,jpj,jpk), & 
     70         &      sclim(jpi,jpj,jpk), STAT= tra_dmp_alloc ) 
    6771      ! 
    6872      CALL mpp_sum ( 'tradmp', tra_dmp_alloc ) 
     
    106110      CALL dta_tsd( kt, zts_dta )            ! read and interpolates T-S data at kt 
    107111      ! 
     112      tclim(:,:,:) = zts_dta(:,:,:,jp_tem) 
     113      sclim(:,:,:) = zts_dta(:,:,:,jp_sal) 
     114      ! 
    108115      SELECT CASE ( nn_zdmp )     !==  type of damping  ==! 
    109116      ! 
Note: See TracChangeset for help on using the changeset viewer.