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

Ignore:
Timestamp:
2016-08-08T12:26:45+02:00 (8 years ago)
Author:
dford
Message:

Initial implementation of observation operator for LogChl?.

File:
1 edited

Legend:

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

    r5838 r6854  
    1111   !!   obs_wri_seaice: Write seaice observation related diagnostics 
    1212   !!   obs_wri_vel   : Write velocity observation diagnostics in NetCDF format 
     13   !!   obs_wri_logchl: Write logchl observation related diagnostics 
    1314   !!   obs_wri_stats : Print basic statistics on the data being written out 
    1415   !!---------------------------------------------------------------------- 
     
    4546      &   obs_wri_seaice, & ! Write seaice observation related diagnostics 
    4647      &   obs_wri_vel, &    ! Write velocity observation related diagnostics 
     48      &   obs_wri_logchl, & ! Write logchl observation related diagnostics 
    4749      &   obswriinfo 
    4850    
     
    930932   END SUBROUTINE obs_wri_vel 
    931933 
     934   SUBROUTINE obs_wri_logchl( cprefix, logchldata, padd, pext ) 
     935      !!----------------------------------------------------------------------- 
     936      !! 
     937      !!                     *** ROUTINE obs_wri_logchl  *** 
     938      !! 
     939      !! ** Purpose : Write logchl observation diagnostics 
     940      !!              related  
     941      !! 
     942      !! ** Method  : NetCDF 
     943      !!  
     944      !! ** Action  : 
     945      !! 
     946      !!----------------------------------------------------------------------- 
     947 
     948      !! * Modules used 
     949      IMPLICIT NONE 
     950 
     951      !! * Arguments 
     952      CHARACTER(LEN=*), INTENT(IN) :: cprefix       ! Prefix for output files 
     953      TYPE(obs_surf), INTENT(INOUT) :: logchldata   ! Full set of logchl 
     954      TYPE(obswriinfo), OPTIONAL :: padd            ! Additional info for each variable 
     955      TYPE(obswriinfo), OPTIONAL :: pext            ! Extra info 
     956 
     957      !! * Local declarations  
     958      TYPE(obfbdata) :: fbdata 
     959      CHARACTER(LEN=40) :: cfname             ! netCDF filename 
     960      CHARACTER(LEN=12), PARAMETER :: cpname = 'obs_wri_logchl' 
     961      INTEGER :: jo 
     962      INTEGER :: ja 
     963      INTEGER :: je 
     964      INTEGER :: nadd 
     965      INTEGER :: next 
     966 
     967      IF ( PRESENT( padd ) ) THEN 
     968         nadd = padd%inum 
     969      ELSE 
     970         nadd = 0 
     971      ENDIF 
     972 
     973      IF ( PRESENT( pext ) ) THEN 
     974         next = pext%inum 
     975      ELSE 
     976         next = 0 
     977      ENDIF 
     978 
     979      CALL init_obfbdata( fbdata ) 
     980 
     981      CALL alloc_obfbdata( fbdata, 1, logchldata%nsurf, 1, & 
     982         &                 1 + nadd, next, .TRUE. ) 
     983 
     984      fbdata%cname(1)      = 'LOGCHL' 
     985      fbdata%coblong(1)    = 'logchl concentration' 
     986      fbdata%cobunit(1)    = 'mg/m3' 
     987      DO je = 1, next 
     988         fbdata%cextname(je) = pext%cdname(je) 
     989         fbdata%cextlong(je) = pext%cdlong(je,1) 
     990         fbdata%cextunit(je) = pext%cdunit(je,1) 
     991      END DO 
     992      fbdata%caddname(1)   = 'Hx' 
     993      fbdata%caddlong(1,1) = 'Model interpolated LOGCHL' 
     994      fbdata%caddunit(1,1) = 'mg/m3' 
     995      fbdata%cgrid(1)      = 'T' 
     996      DO ja = 1, nadd 
     997         fbdata%caddname(1+ja) = padd%cdname(ja) 
     998         fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1) 
     999         fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1) 
     1000      END DO 
     1001 
     1002      WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc 
     1003 
     1004      IF(lwp) THEN 
     1005         WRITE(numout,*) 
     1006         WRITE(numout,*)'obs_wri_logchl :' 
     1007         WRITE(numout,*)'~~~~~~~~~~~~~~~~' 
     1008         WRITE(numout,*)'Writing logchl feedback file : ',TRIM(cfname) 
     1009      ENDIF 
     1010 
     1011      ! Transform obs_prof data structure into obfbdata structure 
     1012      fbdata%cdjuldref = '19500101000000' 
     1013      DO jo = 1, logchldata%nsurf 
     1014         fbdata%plam(jo)      = logchldata%rlam(jo) 
     1015         fbdata%pphi(jo)      = logchldata%rphi(jo) 
     1016         WRITE(fbdata%cdtyp(jo),'(I4)') logchldata%ntyp(jo) 
     1017         fbdata%ivqc(jo,:)    = 0 
     1018         fbdata%ivqcf(:,jo,:) = 0 
     1019         IF ( logchldata%nqc(jo) > 10 ) THEN 
     1020            fbdata%ioqc(jo)    = 4 
     1021            fbdata%ioqcf(1,jo) = 0 
     1022            fbdata%ioqcf(2,jo) = logchldata%nqc(jo) - 10 
     1023         ELSE 
     1024            fbdata%ioqc(jo)    = MAX(logchldata%nqc(jo),1) 
     1025            fbdata%ioqcf(:,jo) = 0 
     1026         ENDIF 
     1027         fbdata%ipqc(jo)      = 0 
     1028         fbdata%ipqcf(:,jo)   = 0 
     1029         fbdata%itqc(jo)      = 0 
     1030         fbdata%itqcf(:,jo)   = 0 
     1031         fbdata%cdwmo(jo)     = '' 
     1032         fbdata%kindex(jo)    = logchldata%nsfil(jo) 
     1033         IF (ln_grid_global) THEN 
     1034            fbdata%iobsi(jo,1) = logchldata%mi(jo) 
     1035            fbdata%iobsj(jo,1) = logchldata%mj(jo) 
     1036         ELSE 
     1037            fbdata%iobsi(jo,1) = mig(logchldata%mi(jo)) 
     1038            fbdata%iobsj(jo,1) = mjg(logchldata%mj(jo)) 
     1039         ENDIF 
     1040         CALL greg2jul( 0, & 
     1041            &           logchldata%nmin(jo), & 
     1042            &           logchldata%nhou(jo), & 
     1043            &           logchldata%nday(jo), & 
     1044            &           logchldata%nmon(jo), & 
     1045            &           logchldata%nyea(jo), & 
     1046            &           fbdata%ptim(jo),   & 
     1047            &           krefdate = 19500101 ) 
     1048         fbdata%padd(1,jo,1,1) = logchldata%rmod(jo,1) 
     1049         fbdata%pob(1,jo,1)    = logchldata%robs(jo,1) 
     1050         fbdata%pdep(1,jo)     = 0.0 
     1051         fbdata%idqc(1,jo)     = 0 
     1052         fbdata%idqcf(:,1,jo)  = 0 
     1053         IF ( logchldata%nqc(jo) > 10 ) THEN 
     1054            fbdata%ivlqc(1,jo,1) = 4 
     1055            fbdata%ivlqcf(1,1,jo,1) = 0 
     1056            fbdata%ivlqcf(2,1,jo,1) = logchldata%nqc(jo) - 10 
     1057         ELSE 
     1058            fbdata%ivlqc(1,jo,1) = MAX(logchldata%nqc(jo),1) 
     1059            fbdata%ivlqcf(:,1,jo,1) = 0 
     1060         ENDIF 
     1061         fbdata%iobsk(1,jo,1)  = 0 
     1062         DO ja = 1, nadd 
     1063            fbdata%padd(1,jo,1+ja,1) = & 
     1064               & logchldata%rext(jo,padd%ipoint(ja)) 
     1065         END DO 
     1066         DO je = 1, next 
     1067            fbdata%pext(1,jo,je) = & 
     1068               & logchldata%rext(jo,pext%ipoint(je)) 
     1069         END DO 
     1070 
     1071      END DO 
     1072 
     1073      ! Write the obfbdata structure 
     1074      CALL write_obfbdata( cfname, fbdata ) 
     1075       
     1076      ! Output some basic statistics 
     1077      CALL obs_wri_stats( fbdata ) 
     1078 
     1079      CALL dealloc_obfbdata( fbdata ) 
     1080 
     1081   END SUBROUTINE obs_wri_logchl 
     1082 
    9321083   SUBROUTINE obs_wri_stats( fbdata ) 
    9331084      !!----------------------------------------------------------------------- 
Note: See TracChangeset for help on using the changeset viewer.