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

Ignore:
Timestamp:
2021-08-13T11:34:58+02:00 (3 years ago)
Author:
dford
Message:

Update treatment of SLA and POTM additional/extra variables.

Location:
NEMO/branches/UKMO/NEMO_4.0.4_generic_obs/src/OCE/OBS
Files:
9 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/UKMO/NEMO_4.0.4_generic_obs/src/OCE/OBS/diaobs.F90

    r15180 r15187  
    249249                  &               sobsgroups(jgroup)%cobstypes ) 
    250250                  ! 
     251               IF( sobsgroups(jgroup)%lsla ) THEN 
     252                  sobsgroups(jgroup)%ssurfdata%cextvars(sobsgroups(jgroup)%next_mdt) = 'MDT' 
     253                  sobsgroups(jgroup)%ssurfdata%cextlong(sobsgroups(jgroup)%next_mdt) = 'Mean dynamic topography' 
     254                  sobsgroups(jgroup)%ssurfdata%cextunit(sobsgroups(jgroup)%next_mdt) = 'Metres' 
     255                  sobsgroups(jgroup)%ssurfdata%caddvars(sobsgroups(jgroup)%nadd_ssh) = 'SSH' 
     256                  DO jvar = 1, sobsgroups(jgroup)%nobstypes 
     257                     sobsgroups(jgroup)%ssurfdata%caddlong(sobsgroups(jgroup)%nadd_ssh,jvar) = 'Model Sea surface height' 
     258                     sobsgroups(jgroup)%ssurfdata%caddunit(sobsgroups(jgroup)%nadd_ssh,jvar) = 'Metres' 
     259                  END DO 
     260               ENDIF 
    251261 
    252262               CALL obs_pre_surf( sobsgroups(jgroup)%ssurfdata,      & 
     
    261271               IF( sobsgroups(jgroup)%lsla ) THEN 
    262272                  CALL obs_rea_mdt( sobsgroups(jgroup)%ssurfdataqc, & 
    263                      &              sobsgroups(jgroup)%n2dint ) 
     273                     &              sobsgroups(jgroup)%n2dint,      & 
     274                     &              sobsgroups(jgroup)%next_mdt ) 
    264275                  IF( sobsgroups(jgroup)%laltbias ) THEN 
    265                      CALL obs_rea_altbias( sobsgroups(jgroup)%ssurfdataqc, & 
    266                         &                  sobsgroups(jgroup)%n2dint,      & 
    267                         &                  sobsgroups(jgroup)%caltbiasfile ) 
     276                     !CALL obs_rea_altbias( sobsgroups(jgroup)%ssurfdataqc, & 
     277                     !   &                  sobsgroups(jgroup)%n2dint,      & 
     278                     !   &                  sobsgroups(jgroup)%caltbiasfile ) 
     279                     CALL obs_app_bias( sobsgroups(jgroup)%ssurfdataqc,   & 
     280                        &               sobsgroups(jgroup)%next_mdt,      &  
     281                        &               sobsgroups(jgroup)%n2dint,        &  
     282                        &               1,                                & 
     283                        &               sobsgroups(jgroup)%caltbiasfile,  & 
     284                        &               'altbias',                        & 
     285                        &               ld_extvar=.TRUE. )  
    268286                  ENDIF 
    269287               ENDIF 
     
    423441                     &               sobsgroups(jgroup)%ravglamscl,        & 
    424442                     &               sobsgroups(jgroup)%ravgphiscl,        & 
    425                      &               sobsgroups(jgroup)%lfp_indegs ) 
     443                     &               sobsgroups(jgroup)%lfp_indegs,        & 
     444                     &               kssh=sobsgroups(jgroup)%nadd_ssh,     & 
     445                     &               kmdt=sobsgroups(jgroup)%next_mdt ) 
    426446 
    427447               END DO 
     
    463483      !! * Local declarations 
    464484      INTEGER :: jgroup                   ! Data set loop variable 
    465       INTEGER :: jo, jvar, jk, jadd, jext 
     485      INTEGER :: jo, jvar, jk, jadd, jext, jadd2, jext2 
    466486      REAL(wp), DIMENSION(:), ALLOCATABLE :: & 
    467487         & zu, & 
    468488         & zv 
     489      LOGICAL, DIMENSION(:), ALLOCATABLE :: ll_write 
    469490      TYPE(obswriinfo) :: sladd, slext 
    470491 
     
    513534                  &                      sobsgroups(jgroup)%sprofdata, .TRUE., numout ) 
    514535 
     536               ! Put additional and extra variable information into obswriinfo structure 
     537               ! used by obs_write. 
     538               ! add/ext variables generated by the OBS code (1...sobsgroups(jgroup)%naddvars) 
     539               ! may duplicate ones read in (%naddvars+1...sobsgroups(jgroup)%sprofdata%nadd) 
     540               ! Check for this, and if so only write out the version generated by the OBS code 
    515541               sladd%inum = sobsgroups(jgroup)%sprofdata%nadd 
     542               ALLOCATE( ll_write(sobsgroups(jgroup)%sprofdata%nadd) ) 
     543               ll_write(:) = .TRUE. 
     544               IF ( (sobsgroups(jgroup)%naddvars > 0) .AND. & 
     545                  & (sobsgroups(jgroup)%sprofdata%nadd > sobsgroups(jgroup)%naddvars) ) THEN 
     546                  DO jadd = sobsgroups(jgroup)%naddvars + 1, sobsgroups(jgroup)%sprofdata%nadd 
     547                     DO jadd2 = 1, sobsgroups(jgroup)%naddvars 
     548                        IF ( TRIM(sobsgroups(jgroup)%sprofdata%caddvars(jadd )) == & 
     549                           & TRIM(sobsgroups(jgroup)%sprofdata%caddvars(jadd2)) ) THEN 
     550                           sladd%inum = sladd%inum - 1 
     551                           ll_write(jadd) = .FALSE. 
     552                        ENDIF 
     553                     END DO 
     554                  END DO 
     555               ENDIF 
    516556               IF ( sladd%inum > 0 ) THEN 
    517557                  ALLOCATE( sladd%ipoint(sladd%inum),                                   & 
     
    519559                     &      sladd%cdlong(sladd%inum,sobsgroups(jgroup)%sprofdata%nvar), & 
    520560                     &      sladd%cdunit(sladd%inum,sobsgroups(jgroup)%sprofdata%nvar) ) 
    521                   DO jadd = 1, sladd%inum 
    522                      sladd%ipoint(jadd) = jadd 
    523                      sladd%cdname(jadd) = sobsgroups(jgroup)%sprofdata%caddvars(jadd) 
    524                      DO jvar = 1, sobsgroups(jgroup)%sprofdata%nvar 
    525                         sladd%cdlong(jadd,jvar) = sobsgroups(jgroup)%sprofdata%caddlong(jadd,jvar) 
    526                         sladd%cdunit(jadd,jvar) = sobsgroups(jgroup)%sprofdata%caddunit(jadd,jvar) 
     561                  jadd2 = 0 
     562                  DO jadd = 1, sobsgroups(jgroup)%sprofdata%nadd 
     563                     IF ( ll_write(jadd) ) THEN 
     564                        jadd2 = jadd2 + 1 
     565                        sladd%ipoint(jadd2) = jadd 
     566                        sladd%cdname(jadd2) = sobsgroups(jgroup)%sprofdata%caddvars(jadd) 
     567                        DO jvar = 1, sobsgroups(jgroup)%sprofdata%nvar 
     568                           sladd%cdlong(jadd2,jvar) = sobsgroups(jgroup)%sprofdata%caddlong(jadd,jvar) 
     569                           sladd%cdunit(jadd2,jvar) = sobsgroups(jgroup)%sprofdata%caddunit(jadd,jvar) 
     570                        END DO 
     571                     ENDIF 
     572                  END DO 
     573               ENDIF 
     574               DEALLOCATE( ll_write ) 
     575                
     576               slext%inum = sobsgroups(jgroup)%sprofdata%next 
     577               ALLOCATE( ll_write(sobsgroups(jgroup)%sprofdata%next) ) 
     578               ll_write(:) = .TRUE. 
     579               IF ( (sobsgroups(jgroup)%nextvars > 0) .AND. & 
     580                  & (sobsgroups(jgroup)%sprofdata%next > sobsgroups(jgroup)%nextvars) ) THEN 
     581                  DO jext = sobsgroups(jgroup)%nextvars + 1, sobsgroups(jgroup)%sprofdata%next 
     582                     DO jext2 = 1, sobsgroups(jgroup)%nextvars 
     583                        IF ( TRIM(sobsgroups(jgroup)%sprofdata%cextvars(jext )) == & 
     584                           & TRIM(sobsgroups(jgroup)%sprofdata%cextvars(jext2)) ) THEN 
     585                           slext%inum = slext%inum - 1 
     586                           ll_write(jext) = .FALSE. 
     587                        ENDIF 
    527588                     END DO 
    528589                  END DO 
    529590               ENDIF 
    530                slext%inum = sobsgroups(jgroup)%sprofdata%next 
    531591               IF ( slext%inum > 0 ) THEN 
    532592                  ALLOCATE( slext%ipoint(slext%inum),   & 
     
    534594                     &      slext%cdlong(slext%inum,1), & 
    535595                     &      slext%cdunit(slext%inum,1) ) 
    536                   DO jext = 1, slext%inum 
    537                      slext%ipoint(jext)   = jext 
    538                      slext%cdname(jext)   = sobsgroups(jgroup)%sprofdata%cextvars(jext) 
    539                      slext%cdlong(jext,1) = sobsgroups(jgroup)%sprofdata%cextlong(jext) 
    540                      slext%cdunit(jext,1) = sobsgroups(jgroup)%sprofdata%cextunit(jext) 
    541                   END DO 
    542                ENDIF 
     596                  jext2 = 0 
     597                  DO jext = 1, sobsgroups(jgroup)%sprofdata%next 
     598                     IF ( ll_write(jext) ) THEN 
     599                        jext2 = jext2 + 1 
     600                        slext%ipoint(jext2)   = jext 
     601                        slext%cdname(jext2)   = sobsgroups(jgroup)%sprofdata%cextvars(jext) 
     602                        slext%cdlong(jext2,1) = sobsgroups(jgroup)%sprofdata%cextlong(jext) 
     603                        slext%cdunit(jext2,1) = sobsgroups(jgroup)%sprofdata%cextunit(jext) 
     604                     ENDIF 
     605                  END DO 
     606               ENDIF 
     607               DEALLOCATE( ll_write ) 
    543608 
    544609               CALL obs_wri_prof( sobsgroups(jgroup)%sprofdata, sobsgroups(jgroup)%cgroupname, sladd, slext ) 
     
    556621                  &                      sobsgroups(jgroup)%ssurfdata, .TRUE., numout ) 
    557622 
     623               ! Put additional and extra variable information into obswriinfo structure 
     624               ! used by obs_write. 
     625               ! add/ext variables generated by the OBS code (1...sobsgroups(jgroup)%naddvars) 
     626               ! may duplicate ones read in (%naddvars+1...sobsgroups(jgroup)%ssurfdata%nadd) 
     627               ! Check for this, and if so only write out the version generated by the OBS code 
    558628               sladd%inum = sobsgroups(jgroup)%ssurfdata%nadd 
     629               ALLOCATE( ll_write(sobsgroups(jgroup)%ssurfdata%nadd) ) 
     630               ll_write(:) = .TRUE. 
     631               IF ( (sobsgroups(jgroup)%naddvars > 0) .AND. & 
     632                  & (sobsgroups(jgroup)%ssurfdata%nadd > sobsgroups(jgroup)%naddvars) ) THEN 
     633                  DO jadd = sobsgroups(jgroup)%naddvars + 1, sobsgroups(jgroup)%ssurfdata%nadd 
     634                     DO jadd2 = 1, sobsgroups(jgroup)%naddvars 
     635                        IF ( TRIM(sobsgroups(jgroup)%ssurfdata%caddvars(jadd )) == & 
     636                           & TRIM(sobsgroups(jgroup)%ssurfdata%caddvars(jadd2)) ) THEN 
     637                           sladd%inum = sladd%inum - 1 
     638                           ll_write(jadd) = .FALSE. 
     639                        ENDIF 
     640                     END DO 
     641                  END DO 
     642               ENDIF 
    559643               IF ( sladd%inum > 0 ) THEN 
    560644                  ALLOCATE( sladd%ipoint(sladd%inum),                                   & 
     
    562646                     &      sladd%cdlong(sladd%inum,sobsgroups(jgroup)%ssurfdata%nvar), & 
    563647                     &      sladd%cdunit(sladd%inum,sobsgroups(jgroup)%ssurfdata%nvar) ) 
    564                   DO jadd = 1, sladd%inum 
    565                      sladd%ipoint(jadd) = jadd 
    566                      sladd%cdname(jadd) = sobsgroups(jgroup)%ssurfdata%caddvars(jadd) 
    567                      DO jvar = 1, sobsgroups(jgroup)%ssurfdata%nvar 
    568                         sladd%cdlong(jadd,jvar) = sobsgroups(jgroup)%ssurfdata%caddlong(jadd,jvar) 
    569                         sladd%cdunit(jadd,jvar) = sobsgroups(jgroup)%ssurfdata%caddunit(jadd,jvar) 
     648                  jadd2 = 0 
     649                  DO jadd = 1, sobsgroups(jgroup)%ssurfdata%nadd 
     650                     IF ( ll_write(jadd) ) THEN 
     651                        jadd2 = jadd2 + 1 
     652                        sladd%ipoint(jadd2) = jadd 
     653                        sladd%cdname(jadd2) = sobsgroups(jgroup)%ssurfdata%caddvars(jadd) 
     654                        DO jvar = 1, sobsgroups(jgroup)%ssurfdata%nvar 
     655                           sladd%cdlong(jadd2,jvar) = sobsgroups(jgroup)%ssurfdata%caddlong(jadd,jvar) 
     656                           sladd%cdunit(jadd2,jvar) = sobsgroups(jgroup)%ssurfdata%caddunit(jadd,jvar) 
     657                        END DO 
     658                     ENDIF 
     659                  END DO 
     660               ENDIF 
     661               DEALLOCATE( ll_write ) 
     662                
     663               slext%inum = sobsgroups(jgroup)%ssurfdata%nextra 
     664               ALLOCATE( ll_write(sobsgroups(jgroup)%ssurfdata%nextra) ) 
     665               ll_write(:) = .TRUE. 
     666               IF ( (sobsgroups(jgroup)%nextvars > 0) .AND. & 
     667                  & (sobsgroups(jgroup)%ssurfdata%nextra > sobsgroups(jgroup)%nextvars) ) THEN 
     668                  DO jext = sobsgroups(jgroup)%nextvars + 1, sobsgroups(jgroup)%ssurfdata%nextra 
     669                     DO jext2 = 1, sobsgroups(jgroup)%nextvars 
     670                        IF ( TRIM(sobsgroups(jgroup)%ssurfdata%cextvars(jext )) == & 
     671                           & TRIM(sobsgroups(jgroup)%ssurfdata%cextvars(jext2)) ) THEN 
     672                           slext%inum = slext%inum - 1 
     673                           ll_write(jext) = .FALSE. 
     674                        ENDIF 
    570675                     END DO 
    571676                  END DO 
    572677               ENDIF 
    573                slext%inum = sobsgroups(jgroup)%ssurfdata%nextra 
    574678               IF ( slext%inum > 0 ) THEN 
    575679                  ALLOCATE( slext%ipoint(slext%inum),   & 
     
    577681                     &      slext%cdlong(slext%inum,1), & 
    578682                     &      slext%cdunit(slext%inum,1) ) 
    579                   DO jext = 1, slext%inum 
    580                      slext%ipoint(jext)   = jext 
    581                      slext%cdname(jext)   = sobsgroups(jgroup)%ssurfdata%cextvars(jext) 
    582                      slext%cdlong(jext,1) = sobsgroups(jgroup)%ssurfdata%cextlong(jext) 
    583                      slext%cdunit(jext,1) = sobsgroups(jgroup)%ssurfdata%cextunit(jext) 
    584                   END DO 
    585                ENDIF 
     683                  jext2 = 0 
     684                  DO jext = 1, sobsgroups(jgroup)%ssurfdata%nextra 
     685                     IF ( ll_write(jext) ) THEN 
     686                        jext2 = jext2 + 1 
     687                        slext%ipoint(jext2)   = jext 
     688                        slext%cdname(jext2)   = sobsgroups(jgroup)%ssurfdata%cextvars(jext) 
     689                        slext%cdlong(jext2,1) = sobsgroups(jgroup)%ssurfdata%cextlong(jext) 
     690                        slext%cdunit(jext2,1) = sobsgroups(jgroup)%ssurfdata%cextunit(jext) 
     691                     ENDIF 
     692                  END DO 
     693               ENDIF 
     694               DEALLOCATE( ll_write ) 
    586695 
    587696               CALL obs_wri_surf( sobsgroups(jgroup)%ssurfdata, sobsgroups(jgroup)%cgroupname, sladd, slext ) 
  • NEMO/branches/UKMO/NEMO_4.0.4_generic_obs/src/OCE/OBS/obs_bias.F90

    r15180 r15187  
    3030CONTAINS 
    3131   SUBROUTINE obs_app_bias( obsdata, kvar, k2dint, knumtypes, & 
    32                             cl_bias_files, cd_biasname ) 
     32                            cl_bias_files, cd_biasname,       & 
     33                            ld_extvar ) 
    3334      !!--------------------------------------------------------------------- 
    3435      !! 
     
    5253      USE iom 
    5354      USE netcdf 
     55 
    5456      !! * Arguments 
    55  
    5657      TYPE(obs_surf), INTENT(INOUT) :: obsdata       ! Observation data 
    5758      INTEGER, INTENT(IN) :: kvar    ! Index of obs type being bias corrected 
     
    6061      CHARACTER(LEN=128), DIMENSION(knumtypes), INTENT(IN) :: & 
    6162                          cl_bias_files !List of files to read 
    62       CHARACTER(LEN=128), INTENT(IN) :: cd_biasname  !Variable name in file 
     63      CHARACTER(LEN=128), INTENT(IN) :: cd_biasname ! Variable name in file 
     64      LOGICAL, OPTIONAL, INTENT(IN)  :: ld_extvar   ! If T correct an extra var 
     65 
    6366      !! * Local declarations 
    6467      INTEGER :: jobs         ! Obs loop variable 
     
    9598         & igrdj_tmp    
    9699      INTEGER :: numobsbias 
    97       INTEGER(KIND=NF90_INT) :: ifile_source 
    98       
     100      INTEGER(KIND=NF90_INT) :: ifile_source      
    99101      INTEGER :: incfile 
    100102      INTEGER :: jtype 
    101103      INTEGER :: iret  
    102104      INTEGER :: inumtype 
     105      LOGICAL :: ll_extvar 
     106 
     107      IF ( PRESENT(ld_extvar) ) THEN 
     108         ll_extvar = ld_extvar 
     109      ELSE 
     110         ll_extvar = .FALSE. 
     111      ENDIF 
     112      IF ( ll_extvar .AND. ( knumtypes /= 1 ) ) THEN 
     113         CALL ctl_stop( 'obs_app_bias: If correcting an extra variable', & 
     114            &           '              knumtypes must be 1' ) 
     115      ENDIF 
     116       
    103117      IF(lwp)WRITE(numout,*)  
    104118      IF(lwp)WRITE(numout,*) 'obs_app_bias : ' 
    105119      IF(lwp)WRITE(numout,*) '----------------- ' 
    106       IF(lwp)WRITE(numout,*) 'Read observation bias for ', TRIM(obsdata%cvars(kvar)) 
     120      IF ( ll_extvar ) THEN 
     121         IF(lwp)WRITE(numout,*) 'Read observation bias for ', TRIM(obsdata%cextvars(kvar)) 
     122      ELSE 
     123         IF(lwp)WRITE(numout,*) 'Read observation bias for ', TRIM(obsdata%cvars(kvar)) 
     124      ENDIF 
     125 
    107126      ! Open and read the files 
    108127      z_obsbias(:,:,:)=0.0_wp 
     
    114133         IF (numobsbias > 0) THEN 
    115134      
    116             !Read the bias type from the file 
    117             !No IOM get attribute command at time of writing,  
    118             !so have to use NETCDF 
    119             !routines directly - should be upgraded in the future 
    120             iret=NF90_OPEN(TRIM(cl_bias_files(jtype)), NF90_NOWRITE, incfile) 
    121             iret=NF90_GET_ATT( incfile, NF90_GLOBAL, TRIM(obsdata%cvars(kvar))//"_source", & 
    122                               ifile_source ) 
    123             ibiastypes(jtype) = ifile_source                  
    124             iret=NF90_CLOSE(incfile)        
    125             
    126             IF ( iret /= 0  ) CALL ctl_stop( & 
    127                'obs_app_bias : Cannot read bias type from file '// & 
    128                cl_bias_files(jtype) ) 
     135            IF ( .NOT. ll_extvar ) THEN 
     136               !Read the bias type from the file 
     137               !No IOM get attribute command at time of writing,  
     138               !so have to use NETCDF 
     139               !routines directly - should be upgraded in the future 
     140               iret=NF90_OPEN(TRIM(cl_bias_files(jtype)), NF90_NOWRITE, incfile) 
     141               IF ( .NOT. ll_extvar ) THEN 
     142                  iret=NF90_GET_ATT( incfile, NF90_GLOBAL, TRIM(obsdata%cvars(kvar))//"_source", & 
     143                                    ifile_source ) 
     144                  ibiastypes(jtype) = ifile_source 
     145               ENDIF 
     146               iret=NF90_CLOSE(incfile) 
     147               IF ( iret /= 0  ) CALL ctl_stop( & 
     148                  'obs_app_bias : Cannot read bias type from file '// & 
     149                  cl_bias_files(jtype) ) 
     150            ENDIF 
     151 
    129152            ! Get the bias data 
    130153            CALL iom_get( numobsbias, jpdom_data, TRIM(cd_biasname), z_obsbias_2d(:,:), 1 ) 
     
    190213         jt=1 
    191214         DO jobs = 1, obsdata%nsurf 
    192             IF ( obsdata%ntyp(jobs) == ibiastypes(jtype) ) THEN 
     215            IF ( ( obsdata%ntyp(jobs) == ibiastypes(jtype) ) .OR. & 
     216               &  ll_extvar ) THEN 
    193217               zlam = obsdata%rlam(jobs) 
    194218               zphi = obsdata%rphi(jobs) 
     
    201225               CALL obs_int_h2d( 1, 1, zweig, zbias(:,:,jt),  zext ) 
    202226               ! adjust observations with bias field 
    203                obsdata%robs(jobs,kvar) = obsdata%robs(jobs,kvar) - zext(1) 
     227               IF ( ll_extvar ) THEN 
     228                  obsdata%rext(jobs,kvar) = obsdata%rext(jobs,kvar) - zext(1) 
     229               ELSE 
     230                  obsdata%robs(jobs,kvar) = obsdata%robs(jobs,kvar) - zext(1) 
     231               ENDIF 
    204232               jt=jt+1 
    205233            ENDIF 
     
    225253         WRITE(numout,*) " " 
    226254         WRITE(numout,*) "Bias correction applied successfully" 
    227          WRITE(numout,*) "Obs types: ",ibiastypes(:), & 
    228                               " Have all been bias corrected\n" 
     255         IF ( .NOT. ll_extvar ) THEN 
     256            WRITE(numout,*) "Obs types: ",ibiastypes(:), & 
     257                                 " Have all been bias corrected" 
     258         ENDIF 
    229259      ENDIF 
    230260   END SUBROUTINE obs_app_bias 
  • NEMO/branches/UKMO/NEMO_4.0.4_generic_obs/src/OCE/OBS/obs_field.F90

    r15180 r15187  
    6363      INTEGER  :: n2dint             !: Type of horizontal interpolation method 
    6464      INTEGER  :: nmsshc             !: MSSH correction scheme 
     65      INTEGER  :: nadd_ssh           !: Index of additional variable representing SSH 
     66      INTEGER  :: next_mdt           !: Index of extra variable representing MDT 
    6567      ! 
    6668      LOGICAL  :: lenabled           !: Logical switch for group being processed and not ignored 
     
    231233         sdobsgroup%lvel3d        = .false. 
    232234         sdobsgroup%lsla          = .false. 
     235         sdobsgroup%nadd_ssh      = 0 
     236         sdobsgroup%next_mdt      = 0 
    233237 
    234238         DO jtype = 1, jpmaxntypes 
     
    265269               ELSEIF ( TRIM(sdobsgroup%cobstypes(itype)) == cobsname_sla ) THEN 
    266270                  sdobsgroup%lsla = .true. 
    267 ! THESE WILL EACH NEED TO BE 1 (ADD=SSH, EXT=MDT) 
    268                   sdobsgroup%naddvars = 0 
    269                   sdobsgroup%nextvars = 0 
     271                  ! SSH=additional, MDT=extra 
     272                  sdobsgroup%naddvars = sdobsgroup%naddvars + 1 
     273                  sdobsgroup%nextvars = sdobsgroup%nextvars + 1 
     274                  sdobsgroup%nadd_ssh = sdobsgroup%naddvars 
     275                  sdobsgroup%next_mdt = sdobsgroup%nextvars 
    270276! DO THIS FOR FBD TOO 
    271277               ENDIF 
  • NEMO/branches/UKMO/NEMO_4.0.4_generic_obs/src/OCE/OBS/obs_grid.F90

    r14075 r15187  
    687687      IF (ln_grid_search_lookup) THEN 
    688688          
    689          WRITE(numout,*) 'Calling obs_grid_setup' 
     689         IF(lwp) WRITE(numout,*) 'Calling obs_grid_setup' 
    690690          
    691691         IF(lwp) WRITE(numout,*) 
     
    724724            ! initially assume size is as defined (to be fixed) 
    725725             
    726             WRITE(numout,*) 'Reading: ',cfname 
     726            IF(lwp) WRITE(numout,*) 'Reading: ',cfname 
    727727             
    728728            CALL chkerr( nf90_open( TRIM( cfname ), nf90_nowrite, idfile ), & 
  • NEMO/branches/UKMO/NEMO_4.0.4_generic_obs/src/OCE/OBS/obs_oper.F90

    r15180 r15187  
    451451 
    452452   SUBROUTINE obs_surf_opt( surfdataqc, kt, kpi, kpj,                & 
    453       &                     kit000, kvar, kdaystp, psurf, psurfmask, & 
     453      &                     kit000, kdaystp, kvar, psurf, psurfmask, & 
    454454      &                     k2dint, ldnightav, plamscl, pphiscl,     & 
    455       &                     lindegrees ) 
     455      &                     lindegrees, kssh, kmdt ) 
    456456 
    457457      !!----------------------------------------------------------------------- 
     
    510510      LOGICAL, INTENT(IN) :: & 
    511511         & lindegrees                  ! T=> plamscl and pphiscl are specified in degrees, F=> in metres 
     512      INTEGER, OPTIONAL, INTENT(IN) :: & 
     513         & kssh                        ! Index of additional variable representing SSH 
     514      INTEGER, OPTIONAL, INTENT(IN) :: & 
     515         & kmdt                        ! Index of extra variable representing MDT 
    512516 
    513517      !! * Local declarations 
     
    739743         ELSE 
    740744 
    741             ! Get weights to average the model SLA to the observation footprint 
     745            ! Get weights to average the model field to the observation footprint 
    742746            CALL obs_avg_h2d_init( 1, 1, imaxifp, imaxjfp, k2dint, zlam,  zphi, & 
    743747               &                   zglam(:,:,iobs), zgphi(:,:,iobs), & 
     
    746750               &                   lindegrees, zweig, zobsmask ) 
    747751 
    748             ! Average the model SST to the observation footprint 
     752            ! Average the model field to the observation footprint 
    749753            CALL obs_avg_h2d( 1, 1, imaxifp, imaxjfp, & 
    750754               &              zweig, zsurftmp(:,:,iobs),  zext ) 
    751755 
    752756         ENDIF 
    753 ! WHERE BEST TO DO THIS? 
    754          IF ( TRIM(surfdataqc%cvars(1)) == 'SLA' .AND. surfdataqc%nextra == 2 ) THEN 
     757 
     758         IF ( TRIM(surfdataqc%cvars(kvar)) == 'SLA' .AND. PRESENT(kssh) .AND. PRESENT(kmdt) ) THEN 
    755759            ! ... Remove the MDT from the SSH at the observation point to get the SLA 
    756             surfdataqc%rext(jobs,1) = zext(1) 
    757             surfdataqc%rmod(jobs,1) = surfdataqc%rext(jobs,1) - surfdataqc%rext(jobs,2) 
     760            surfdataqc%radd(jobs,kssh,kvar) = zext(1) 
     761            surfdataqc%rmod(jobs,kvar) = surfdataqc%radd(jobs,kssh,kvar) - surfdataqc%rext(jobs,kmdt) 
    758762         ELSE 
    759             surfdataqc%rmod(jobs,1) = zext(1) 
    760          ENDIF 
    761           
     763            surfdataqc%rmod(jobs,kvar) = zext(1) 
     764         ENDIF 
     765! DO THIS FOR FBD TOO 
     766 
    762767         IF ( zext(1) == obfillflt ) THEN 
    763768            ! If the observation value is a fill value, set QC flag to bad 
  • NEMO/branches/UKMO/NEMO_4.0.4_generic_obs/src/OCE/OBS/obs_read_prof.F90

    r15180 r15187  
    235235            IF ( inpfiles(jj)%nvar /= kvars ) THEN 
    236236               CALL ctl_stop( 'Feedback format error: ', & 
    237                   &           ' unexpected number of vars in profile file' ) 
     237                  &           ' unexpected number of vars in feedback file', & 
     238                  &           TRIM(cdfilenames(jj)) ) 
    238239            ENDIF 
    239240 
    240241            IF ( ldmod .AND. ( inpfiles(jj)%nadd == 0 ) ) THEN 
    241                CALL ctl_stop( 'Model not in input data' ) 
     242               CALL ctl_stop( 'Model not in input data in', & 
     243                  &           TRIM(cdfilenames(jj)) ) 
     244               RETURN 
    242245            ENDIF 
    243246 
    244247            IF ( (iextr > 0) .AND. (inpfiles(jj)%next /= iextr) ) THEN 
    245248               CALL ctl_stop( 'Number of extra variables not consistent', & 
    246                   &           ' with previous files for this type' ) 
     249                  &           ' with previous files for this type in', & 
     250                  &           TRIM(cdfilenames(jj)) ) 
    247251            ELSE 
    248252               iextr = inpfiles(jj)%next 
     
    258262            END DO 
    259263            IF ( ldmod .AND. ( inpfiles(jj)%nadd == iaddin ) ) THEN 
    260                CALL ctl_stop( 'Model not in input data' ) 
     264               CALL ctl_stop( 'Model not in input data', & 
     265                  &           TRIM(cdfilenames(jj)) ) 
    261266            ENDIF 
    262267 
    263268            IF ( (iadd > 0) .AND. (iaddin /= iadd) ) THEN 
    264269               CALL ctl_stop( 'Number of additional variables not consistent', & 
    265                   &           ' with previous files for this type' ) 
     270                  &           ' with previous files for this type in', & 
     271                  &           TRIM(cdfilenames(jj)) ) 
    266272            ELSE 
    267273               iadd = iaddin 
     
    313319                  IF ( inpfiles(jj)%cname(ji) /= clvarsin(ji) ) THEN 
    314320                     CALL ctl_stop( 'Feedback file variables not consistent', & 
    315                         &           ' with previous files for this type' ) 
     321                        &           ' with previous files for this type in', & 
     322                        &           TRIM(cdfilenames(jj)) ) 
    316323                  ENDIF 
    317324               END DO 
     
    323330                        IF ( inpfiles(jj)%caddname(ji) /= claddvarsin(jadd) ) THEN 
    324331                           CALL ctl_stop( 'Feedback file additional variables not consistent', & 
    325                               &           ' with previous files for this type' ) 
     332                              &           ' with previous files for this type in', & 
     333                              &           TRIM(cdfilenames(jj)) ) 
    326334                        ENDIF 
    327335                     ENDIF 
     
    332340                     IF ( inpfiles(jj)%cextname(ji) /= clextvarsin(ji) ) THEN 
    333341                        CALL ctl_stop( 'Feedback file extra variables not consistent', & 
    334                            &           ' with previous files for this type' ) 
     342                           &           ' with previous files for this type in', & 
     343                           &           TRIM(cdfilenames(jj)) ) 
    335344                     ENDIF 
    336345                  END DO 
  • NEMO/branches/UKMO/NEMO_4.0.4_generic_obs/src/OCE/OBS/obs_read_surf.F90

    r15180 r15187  
    204204            IF ( inpfiles(jj)%nvar /= kvars ) THEN 
    205205               CALL ctl_stop( 'Feedback format error: ', & 
    206                   &           ' unexpected number of vars in feedback file' ) 
     206                  &           ' unexpected number of vars in feedback file', & 
     207                  &           TRIM(cdfilenames(jj)) ) 
    207208            ENDIF 
    208209 
    209210            IF ( ldmod .AND. ( inpfiles(jj)%nadd == 0 ) ) THEN 
    210                CALL ctl_stop( 'Model not in input data' ) 
     211               CALL ctl_stop( 'Model not in input data in', & 
     212                  &           TRIM(cdfilenames(jj)) ) 
    211213               RETURN 
    212214            ENDIF 
     
    214216            IF ( (iextr > 0) .AND. (inpfiles(jj)%next /= iextr) ) THEN 
    215217               CALL ctl_stop( 'Number of extra variables not consistent', & 
    216                   &           ' with previous files for this type' ) 
     218                  &           ' with previous files for this type in', & 
     219                  &           TRIM(cdfilenames(jj)) ) 
    217220            ELSE 
    218221               iextr = inpfiles(jj)%next 
     
    228231            END DO 
    229232            IF ( ldmod .AND. ( inpfiles(jj)%nadd == iaddin ) ) THEN 
    230                CALL ctl_stop( 'Model not in input data' ) 
     233               CALL ctl_stop( 'Model not in input data', & 
     234                  &           TRIM(cdfilenames(jj)) ) 
    231235            ENDIF 
    232236 
    233237            IF ( (iadd > 0) .AND. (iaddin /= iadd) ) THEN 
    234238               CALL ctl_stop( 'Number of additional variables not consistent', & 
    235                   &           ' with previous files for this type' ) 
     239                  &           ' with previous files for this type in', & 
     240                  &           TRIM(cdfilenames(jj)) ) 
    236241            ELSE 
    237242               iadd = iaddin 
     
    283288                  IF ( inpfiles(jj)%cname(ji) /= clvarsin(ji) ) THEN 
    284289                     CALL ctl_stop( 'Feedback file variables not consistent', & 
    285                         &           ' with previous files for this type' ) 
     290                        &           ' with previous files for this type in', & 
     291                        &           TRIM(cdfilenames(jj)) ) 
    286292                  ENDIF 
    287293               END DO 
     
    293299                        IF ( inpfiles(jj)%caddname(ji) /= claddvarsin(jadd) ) THEN 
    294300                           CALL ctl_stop( 'Feedback file additional variables not consistent', & 
    295                               &           ' with previous files for this type' ) 
     301                              &           ' with previous files for this type in', & 
     302                              &           TRIM(cdfilenames(jj)) ) 
    296303                        ENDIF 
    297304                     ENDIF 
     
    302309                     IF ( inpfiles(jj)%cextname(ji) /= clextvarsin(ji) ) THEN 
    303310                        CALL ctl_stop( 'Feedback file extra variables not consistent', & 
    304                            &           ' with previous files for this type' ) 
     311                           &           ' with previous files for this type in', & 
     312                           &           TRIM(cdfilenames(jj)) ) 
    305313                     ENDIF 
    306314                  END DO 
  • NEMO/branches/UKMO/NEMO_4.0.4_generic_obs/src/OCE/OBS/obs_readmdt.F90

    r15180 r15187  
    4444CONTAINS 
    4545 
    46    SUBROUTINE obs_rea_mdt( sladata, k2dint ) 
     46   SUBROUTINE obs_rea_mdt( sladata, k2dint, kmdt ) 
    4747      !!--------------------------------------------------------------------- 
    4848      !! 
     
    5959      TYPE(obs_surf), INTENT(inout) ::   sladata   ! SLA data 
    6060      INTEGER       , INTENT(in)    ::   k2dint    ! ? 
     61      INTEGER       , INTENT(in)    ::   kmdt      ! Index of MDT extra var 
    6162      ! 
    6263      CHARACTER(LEN=12), PARAMETER ::   cpname  = 'obs_rea_mdt' 
     
    148149         CALL obs_int_h2d( 1, 1, zweig, zmdtl(:,:,jobs),  zext ) 
    149150 
    150 ! FIGURE OUT THIS ASSIGNMENT  
    151 !         sladata%rext(jobs,2) = zext(1) 
     151         sladata%rext(jobs,kmdt) = zext(1) 
    152152 
    153153! mark any masked data with a QC flag 
     
    247247         WRITE(numout,*) '               zcorr         = ', zcorr 
    248248         WRITE(numout,*) '               nn_msshc        = ', nn_msshc 
     249 
     250         IF ( nn_msshc == 0 ) WRITE(numout,*) '           MSSH correction is not applied' 
     251         IF ( nn_msshc == 1 ) WRITE(numout,*) '           MSSH correction is applied' 
     252         IF ( nn_msshc == 2 ) WRITE(numout,*) '           User defined MSSH correction'  
    249253      ENDIF 
    250  
    251       IF ( nn_msshc == 0 ) WRITE(numout,*) '           MSSH correction is not applied' 
    252       IF ( nn_msshc == 1 ) WRITE(numout,*) '           MSSH correction is applied' 
    253       IF ( nn_msshc == 2 ) WRITE(numout,*) '           User defined MSSH correction'  
    254254 
    255255      ! 
  • NEMO/branches/UKMO/NEMO_4.0.4_generic_obs/src/OCE/OBS/obs_write.F90

    r15180 r15187  
    9191      INTEGER :: ilevel 
    9292      INTEGER :: jvar 
     93      INTEGER :: jvar2 
     94      INTEGER :: jsal 
    9395      INTEGER :: jo 
    9496      INTEGER :: jk 
     
    145147         ENDIF 
    146148      END DO 
    147 !fbdata%cextname(1)   = 'TEMP' 
    148 !fbdata%cextlong(1)   = 'Insitu temperature' 
    149 !fbdata%cextunit(1)   = 'Degrees centigrade' 
    150149 
    151150      WRITE(clfname, FMT="(A,'fb_fdbk_',I4.4,'.nc')") TRIM(clfiletype), nproc 
     
    230229                  END DO 
    231230               ENDIF 
    232 !IF ( ( jvar == 1 ) .AND. & 
    233 !   & ( TRIM(profdata%cvars(1)) == 'POTM' ) ) THEN 
    234 !   fbdata%pext(ik,jo,1) = profdata%var(jvar)%vext(jk,1) 
    235 !ENDIF  
    236231            END DO 
    237232         END DO 
    238233      END DO 
    239234 
    240 !IF ( TRIM(profdata%cvars(1)) == 'POTM' ) THEN 
    241 !   ! Convert insitu temperature to potential temperature using the model 
    242 !   ! salinity if no potential temperature 
    243 !   DO jo = 1, fbdata%nobs 
    244 !      IF ( fbdata%pphi(jo) < 9999.0 ) THEN 
    245 !         DO jk = 1, fbdata%nlev 
    246 !            IF ( ( fbdata%pob(jk,jo,1) >= 9999.0 ) .AND. & 
    247 !               & ( fbdata%pdep(jk,jo) < 9999.0 ) .AND. & 
    248 !               & ( fbdata%padd(jk,jo,1,2) < 9999.0 ) .AND. & 
    249 !               & ( fbdata%pext(jk,jo,1) < 9999.0 ) ) THEN 
    250 !               zpres = dep_to_p( REAL(fbdata%pdep(jk,jo),wp), & 
    251 !                  &              REAL(fbdata%pphi(jo),wp) ) 
    252 !               fbdata%pob(jk,jo,1) = potemp( & 
    253 !                  &                     REAL(fbdata%padd(jk,jo,1,2), wp),  & 
    254 !                  &                     REAL(fbdata%pext(jk,jo,1), wp), & 
    255 !                  &                     zpres, 0.0_wp ) 
    256 !            ENDIF 
    257 !         END DO 
    258 !      ENDIF 
    259 !   END DO 
    260 !ENDIF 
     235      ! Convert insitu temperature to potential temperature using the model 
     236      ! salinity if no potential temperature 
     237      IF (iext > 0) THEN 
     238         DO jvar = 1, profdata%nvar 
     239            IF ( TRIM(profdata%cvars(jvar)) == 'POTM' ) THEN 
     240               jsal = 0 
     241               DO jvar2 = 1, profdata%nvar 
     242                  IF ( TRIM(profdata%cvars(jvar2)) == 'PSAL' ) THEN 
     243                     jsal = jvar2 
     244                     EXIT 
     245                  ENDIF 
     246               END DO 
     247               IF (jsal > 0) THEN 
     248                  DO je = 1, iext 
     249                     IF ( TRIM(fbdata%cextname(je)) == 'TEMP' ) THEN 
     250                        DO jo = 1, fbdata%nobs 
     251                           IF ( fbdata%pphi(jo) < 9999.0 ) THEN 
     252                              DO jk = 1, fbdata%nlev 
     253                                 IF ( ( fbdata%pob(jk,jo,jvar)   >= 9999.0 ) .AND. & 
     254                                    & ( fbdata%pdep(jk,jo)        < 9999.0 ) .AND. & 
     255                                    & ( fbdata%padd(jk,jo,1,jsal) < 9999.0 ) .AND. & 
     256                                    & ( fbdata%pext(jk,jo,je)     < 9999.0 ) ) THEN 
     257                                    zpres = dep_to_p( REAL(fbdata%pdep(jk,jo),wp), & 
     258                                       &              REAL(fbdata%pphi(jo),wp) ) 
     259                                    fbdata%pob(jk,jo,jvar) = potemp( & 
     260                                       &                     REAL(fbdata%padd(jk,jo,1,jsal), wp), & 
     261                                       &                     REAL(fbdata%pext(jk,jo,je), wp),     & 
     262                                       &                     zpres, 0.0_wp ) 
     263                                 ENDIF 
     264                              END DO 
     265                           ENDIF 
     266                        END DO 
     267                        EXIT 
     268                     ENDIF 
     269                  END DO 
     270               ENDIF 
     271               EXIT 
     272            ENDIF 
     273         END DO 
     274      ENDIF 
    261275 
    262276      ! Write the obfbdata structure 
Note: See TracChangeset for help on using the changeset viewer.