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 11468 for branches/UKMO/dev_r5518_obs_oper_update/NEMOGCM/NEMO/OPA_SRC/OBS/obs_write.F90 – NEMO

Ignore:
Timestamp:
2019-08-23T10:37:22+02:00 (5 years ago)
Author:
mattmartin
Message:

Merged changes to allow writing of climatological information to feedback files.

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/UKMO/dev_r5518_obs_oper_update/NEMOGCM/NEMO/OPA_SRC/OBS/obs_write.F90

    r9308 r11468  
    9595      INTEGER :: je 
    9696      INTEGER :: iadd 
     97      INTEGER :: iadd_clm ! 1 if climatology present 
    9798      INTEGER :: iext 
    9899      REAL(wp) :: zpres 
    99100 
     101 
     102      iadd_clm = 0  
     103      IF ( profdata%lclim ) iadd_clm = 1 
     104       
    100105      IF ( PRESENT( padd ) ) THEN 
    101106         iadd = padd%inum 
     
    123128         clfiletype='profb' 
    124129         CALL alloc_obfbdata( fbdata, 2, profdata%nprof, ilevel, & 
    125             &                 1 + iadd, 1 + iext, .TRUE. ) 
     130            &                 1 + iadd_clm + iadd, 1 + iext, .TRUE. ) 
    126131         fbdata%cname(1)      = profdata%cvars(1) 
    127132         fbdata%cname(2)      = profdata%cvars(2) 
     
    137142         fbdata%caddunit(1,1) = 'Degrees centigrade' 
    138143         fbdata%caddunit(1,2) = 'PSU' 
     144         IF ( profdata%lclim ) THEN 
     145            fbdata%caddlong(2,1) = 'Climatology interpolated potential temperature' 
     146            fbdata%caddlong(2,2) = 'Climatology interpolated practical salinity' 
     147            fbdata%caddunit(2,1) = 'Degrees centigrade' 
     148            fbdata%caddunit(2,2) = 'PSU' 
     149         ENDIF 
    139150         fbdata%cgrid(:)      = 'T' 
    140151         DO je = 1, iext 
     
    144155         END DO 
    145156         DO ja = 1, iadd 
    146             fbdata%caddname(1+ja) = padd%cdname(ja) 
     157            fbdata%caddname(1+iadd_clm+ja) = padd%cdname(ja) 
    147158            DO jvar = 1, 2 
    148                fbdata%caddlong(1+ja,jvar) = padd%cdlong(ja,jvar) 
    149                fbdata%caddunit(1+ja,jvar) = padd%cdunit(ja,jvar) 
     159               fbdata%caddlong(1+iadd_clm+ja,jvar) = padd%cdlong(ja,jvar) 
     160               fbdata%caddunit(1+iadd_clm+ja,jvar) = padd%cdunit(ja,jvar) 
    150161            END DO 
    151162         END DO 
     
    154165 
    155166         clfiletype='velfb' 
    156          CALL alloc_obfbdata( fbdata, 2, profdata%nprof, ilevel, 1, 0, .TRUE. ) 
     167         CALL alloc_obfbdata( fbdata, 2, profdata%nprof, ilevel, & 
     168            &                 1 + iadd_clm + iadd, 0, .TRUE. ) 
    157169         fbdata%cname(1)      = profdata%cvars(1) 
    158170         fbdata%cname(2)      = profdata%cvars(2) 
     
    170182         fbdata%caddunit(1,1) = 'm/s' 
    171183         fbdata%caddunit(1,2) = 'm/s' 
     184         IF ( profdata%lclim ) THEN 
     185            fbdata%caddlong(2,1) = 'Climatology interpolated zonal velocity' 
     186            fbdata%caddlong(2,2) = 'Climatology interpolated meridional velocity' 
     187            fbdata%caddunit(2,1) = 'm/s' 
     188            fbdata%caddunit(2,2) = 'm/s' 
     189         ENDIF          
    172190         fbdata%cgrid(1)      = 'U'  
    173191         fbdata%cgrid(2)      = 'V' 
    174192         DO ja = 1, iadd 
    175             fbdata%caddname(1+ja) = padd%cdname(ja) 
    176             fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1) 
    177             fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1) 
     193            fbdata%caddname(1+iadd_clm+ja) = padd%cdname(ja) 
     194            fbdata%caddlong(1+iadd_clm+ja,1) = padd%cdlong(ja,1) 
     195            fbdata%caddunit(1+iadd_clm+ja,1) = padd%cdunit(ja,1) 
    178196         END DO 
    179197 
     
    246264         & ( TRIM(profdata%cvars(1)) /= 'UVEL' ) ) THEN 
    247265         CALL alloc_obfbdata( fbdata, 1, profdata%nprof, ilevel, & 
    248             &                 1 + iadd, iext, .TRUE. ) 
     266            &                 1 + iadd_clm + iadd, iext, .TRUE. ) 
    249267         fbdata%cname(1)      = profdata%cvars(1) 
    250268         fbdata%coblong(1)    = cllongname 
     
    252270         fbdata%caddlong(1,1) = 'Model interpolated ' // TRIM(cllongname) 
    253271         fbdata%caddunit(1,1) = clunits 
     272         IF ( profdata%lclim ) THEN 
     273            fbdata%caddlong(2,1) = 'Climatological interpolated ' // TRIM(cllongname) 
     274            fbdata%caddunit(2,1) = clunits 
     275         ENDIF          
    254276         fbdata%cgrid(:)      = clgrid 
    255277         DO je = 1, iext 
     
    259281         END DO 
    260282         DO ja = 1, iadd 
    261             fbdata%caddname(1+ja) = padd%cdname(ja) 
    262             fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1) 
    263             fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1) 
     283            fbdata%caddname(1+iadd_clm+ja) = padd%cdname(ja) 
     284            fbdata%caddlong(1+iadd_clm+ja,1) = padd%cdlong(ja,1) 
     285            fbdata%caddunit(1+iadd_clm+ja,1) = padd%cdunit(ja,1) 
    264286         END DO 
    265287      ENDIF 
    266288 
    267289      fbdata%caddname(1)   = 'Hx' 
    268  
     290      IF ( profdata%lclim ) fbdata%caddname(1+iadd_clm)   = 'CLM' 
     291       
    269292      WRITE(clfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(clfiletype), nproc 
    270293 
     
    319342            DO jk = profdata%npvsta(jo,jvar), profdata%npvend(jo,jvar) 
    320343               ik = profdata%var(jvar)%nvlidx(jk) 
    321                fbdata%padd(ik,jo,1,jvar) = profdata%var(jvar)%vmod(jk) 
    322344               fbdata%pob(ik,jo,jvar)    = profdata%var(jvar)%vobs(jk) 
    323345               fbdata%pdep(ik,jo)        = profdata%var(jvar)%vdep(jk) 
     
    333355               ENDIF 
    334356               fbdata%iobsk(ik,jo,jvar)  = profdata%var(jvar)%mvk(jk) 
     357                
     358               fbdata%padd(ik,jo,1,jvar) = profdata%var(jvar)%vmod(jk) 
     359               IF ( profdata%lclim ) THEN            
     360                  fbdata%padd(ik,jo,1+iadd_clm,jvar) = profdata%var(jvar)%vclm(jk)      
     361               ENDIF               
    335362               DO ja = 1, iadd 
    336                   fbdata%padd(ik,jo,1+ja,jvar) = & 
     363                  fbdata%padd(ik,jo,1+iadd_clm+ja,jvar) = & 
    337364                     & profdata%var(jvar)%vext(jk,padd%ipoint(ja)) 
    338365               END DO 
     
    420447      INTEGER :: indx_std 
    421448      INTEGER :: iadd_std 
    422  
     449      INTEGER :: iadd_clm      
     450      INTEGER :: iadd_mdt  
     451 
     452      IF ( PRESENT( pext ) ) THEN 
     453         iext = pext%inum 
     454      ELSE 
     455         iext = 0 
     456      ENDIF 
     457 
     458 
     459      ! Set up number of additional variables to be ouput: 
     460      ! Hx, CLM, STD, MDT... 
     461  
    423462      IF ( PRESENT( padd ) ) THEN 
    424463         iadd = padd%inum 
     
    426465         iadd = 0 
    427466      ENDIF 
    428  
    429       IF ( PRESENT( pext ) ) THEN 
    430          iext = pext%inum 
    431       ELSE 
    432          iext = 0 
    433       ENDIF 
    434  
     467       
    435468      iadd_std = 0 
    436469      indx_std = -1 
     
    444477      ENDIF 
    445478       
     479      iadd_clm = 0 
     480      IF ( surfdata%lclim ) iadd_clm = 1 
     481       
     482      iadd_mdt = 0 
     483      IF ( TRIM(surfdata%cvars(1)) == 'SLA' ) iadd_mdt = 1 
     484       
    446485      CALL init_obfbdata( fbdata ) 
    447486 
    448487      SELECT CASE ( TRIM(surfdata%cvars(1)) ) 
    449488      CASE('SLA') 
    450  
     489          
    451490         ! SLA needs special treatment because of MDT, so is all done here 
    452491         ! Other variables are done more generically 
     492         ! No climatology for SLA, MDT is our best estimate of that and is already output. 
    453493 
    454494         CALL alloc_obfbdata( fbdata, 1, surfdata%nsurf, 1, & 
    455             &                 2 + iadd + iadd_std, 1 + iext, .TRUE. ) 
     495            &                 1 + iadd_mdt + iadd_std + iadd, & 
     496            &                 1 + iext, .TRUE. ) 
    456497 
    457498         clfiletype = 'slafb' 
     
    474515         fbdata%cgrid(1)      = 'T' 
    475516         DO ja = 1, iadd 
    476             fbdata%caddname(2+iadd_std+ja) = padd%cdname(ja) 
    477             fbdata%caddlong(2+iadd_std+ja,1) = padd%cdlong(ja,1) 
    478             fbdata%caddunit(2+iadd_std+ja,1) = padd%cdunit(ja,1) 
     517            fbdata%caddname(1+iadd_mdt+iadd_std+ja) = padd%cdname(ja) 
     518            fbdata%caddlong(1+iadd_mdt+iadd_std+ja,1) = padd%cdlong(ja,1) 
     519            fbdata%caddunit(1+iadd_mdt+iadd_std+ja,1) = padd%cdunit(ja,1) 
    479520         END DO 
    480521 
     
    485526         clunits    = 'Degree centigrade' 
    486527         clgrid     = 'T' 
    487  
     528          
    488529      CASE('ICECONC') 
    489530 
     
    499540         clunits    = 'psu' 
    500541         clgrid     = 'T' 
    501  
     542          
    502543      CASE('SLCHLTOT','LOGCHL','LogChl','logchl') 
    503544 
     
    610651       
    611652         CALL alloc_obfbdata( fbdata, 1, surfdata%nsurf, 1, & 
    612             &                 1 + iadd + iadd_std, iext, .TRUE. ) 
     653            &                 1 + iadd_std + iadd_clm + iadd, iext, .TRUE. ) 
    613654 
    614655         fbdata%cname(1)      = surfdata%cvars(1) 
     
    628669         fbdata%cgrid(1)      = clgrid 
    629670         DO ja = 1, iadd 
    630             fbdata%caddname(1+iadd_std+ja) = padd%cdname(ja) 
    631             fbdata%caddlong(1+iadd_std+ja,1) = padd%cdlong(ja,1) 
    632             fbdata%caddunit(1+iadd_std+ja,1) = padd%cdunit(ja,1) 
     671            fbdata%caddname(1+iadd_mdt+iadd_std+iadd_clm+ja) = padd%cdname(ja) 
     672            fbdata%caddlong(1+iadd_mdt+iadd_std+iadd_clm+ja,1) = padd%cdlong(ja,1) 
     673            fbdata%caddunit(1+iadd_mdt+iadd_std+iadd_clm+ja,1) = padd%cdunit(ja,1) 
    633674         END DO 
    634675 
     
    637678      fbdata%caddname(1)   = 'Hx' 
    638679      IF ( indx_std /= -1 ) THEN 
    639          IF ( TRIM(surfdata%cvars(1)) == 'SLA' ) iadd_std = iadd_std + 1 
    640          fbdata%caddname(1+iadd_std)   = surfdata%cext(indx_std) 
    641          fbdata%caddlong(1+iadd_std,1) = 'Obs error standard deviation' 
    642          fbdata%caddunit(1+iadd_std,1) = fbdata%cobunit(1) 
    643       ENDIF 
     680         fbdata%caddname(1+iadd_mdt+iadd_std)   = surfdata%cext(indx_std) 
     681         fbdata%caddlong(1+iadd_mdt+iadd_std,1) = 'Obs error standard deviation' 
     682         fbdata%caddunit(1+iadd_mdt+iadd_std,1) = fbdata%cobunit(1) 
     683      ENDIF 
     684       
     685      IF ( surfdata%lclim ) THEN 
     686         fbdata%caddname(1+iadd_mdt+iadd_std+iadd_clm)   = 'CLM' 
     687         fbdata%caddlong(1+iadd_mdt+iadd_std+iadd_clm,1) = 'Climatology' 
     688         fbdata%caddunit(1+iadd_mdt+iadd_std+iadd_clm,1) = fbdata%cobunit(1) 
     689      ENDIF 
     690       
    644691       
    645692      WRITE(clfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(clfiletype), nproc 
     
    689736            &           fbdata%ptim(jo),   & 
    690737            &           krefdate = 19500101 ) 
    691          fbdata%padd(1,jo,1,1) = surfdata%rmod(jo,1) 
    692          IF ( TRIM(surfdata%cvars(1)) == 'SLA' ) fbdata%padd(1,jo,2,1) = surfdata%rext(jo,1) 
     738                     
    693739         fbdata%pob(1,jo,1)    = surfdata%robs(jo,1)  
    694740         fbdata%pdep(1,jo)     = 0.0 
     
    706752         ENDIF 
    707753         fbdata%iobsk(1,jo,1)  = 0 
    708          IF ( TRIM(surfdata%cvars(1)) == 'SLA' ) fbdata%pext(1,jo,1) = surfdata%rext(jo,2) 
     754  
     755         ! Additional variables. 
     756         ! Hx is always the first additional variable 
     757         fbdata%padd(1,jo,1,1) = surfdata%rmod(jo,1) 
     758         ! MDT is output as an additional variable if SLA obs type 
     759         IF ( TRIM(surfdata%cvars(1)) == 'SLA' ) THEN 
     760            fbdata%padd(1,jo,1+iadd_mdt,1) = surfdata%rext(jo,1) 
     761         ENDIF     
     762         ! STD is output as an additional variable if available 
    709763         IF ( indx_std /= -1 ) THEN 
    710             fbdata%padd(1,jo,1+iadd_std,1) = surfdata%rext(jo,indx_std) 
     764            fbdata%padd(1,jo,1+iadd_mdt+iadd_std,1) = surfdata%rext(jo,indx_std) 
    711765         ENDIF 
     766         ! CLM is output as an additional variable if available 
     767         IF ( surfdata%lclim ) THEN 
     768            fbdata%padd(1,jo,1+iadd_mdt+iadd_std+iadd_clm,1) = surfdata%rclm(jo,1) 
     769         ENDIF 
     770         ! Then other additional variables are output 
     771         DO ja = 1, iadd 
     772            fbdata%padd(1,jo,1+iadd_mdt+iadd_std+iadd_clm+ja,1) = & 
     773               & surfdata%rext(jo,padd%ipoint(ja)) 
     774         END DO 
    712775          
    713          DO ja = 1, iadd 
    714             fbdata%padd(1,jo,2+iadd_std+ja,1) = & 
    715                & surfdata%rext(jo,padd%ipoint(ja)) 
    716          END DO 
     776         ! Extra variables 
     777         IF ( TRIM(surfdata%cvars(1)) == 'SLA' ) fbdata%pext(1,jo,1) = surfdata%rext(jo,2)          
    717778         DO je = 1, iext 
    718779            fbdata%pext(1,jo,1+je) = & 
Note: See TracChangeset for help on using the changeset viewer.