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

Ignore:
Timestamp:
2017-02-22T12:40:19+01:00 (7 years ago)
Author:
dford
Message:

Add observation operator code for surface log10(chlorophyll), SPM, pCO2 and fCO2, for use with FABM-ERSEM, HadOCC and MEDUSA. See internal Met Office NEMO ticket 660.

File:
1 edited

Legend:

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

    r5838 r7713  
    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 
     14   !!   obs_wri_spm   : Write spm observation related diagnostics 
     15   !!   obs_wri_fco2  : Write fco2 observation related diagnostics 
     16   !!   obs_wri_pco2  : Write fco2 observation related diagnostics 
    1317   !!   obs_wri_stats : Print basic statistics on the data being written out 
    1418   !!---------------------------------------------------------------------- 
     
    4549      &   obs_wri_seaice, & ! Write seaice observation related diagnostics 
    4650      &   obs_wri_vel, &    ! Write velocity observation related diagnostics 
     51      &   obs_wri_logchl, & ! Write logchl observation related diagnostics 
     52      &   obs_wri_spm, &    ! Write spm observation related diagnostics 
     53      &   obs_wri_fco2, &   ! Write fco2 observation related diagnostics 
     54      &   obs_wri_pco2, &   ! Write pco2 observation related diagnostics 
    4755      &   obswriinfo 
    4856    
     
    930938   END SUBROUTINE obs_wri_vel 
    931939 
     940   SUBROUTINE obs_wri_logchl( cprefix, logchldata, padd, pext ) 
     941      !!----------------------------------------------------------------------- 
     942      !! 
     943      !!                     *** ROUTINE obs_wri_logchl  *** 
     944      !! 
     945      !! ** Purpose : Write logchl observation diagnostics 
     946      !!              related  
     947      !! 
     948      !! ** Method  : NetCDF 
     949      !!  
     950      !! ** Action  : 
     951      !! 
     952      !!----------------------------------------------------------------------- 
     953 
     954      !! * Modules used 
     955      IMPLICIT NONE 
     956 
     957      !! * Arguments 
     958      CHARACTER(LEN=*), INTENT(IN) :: cprefix       ! Prefix for output files 
     959      TYPE(obs_surf), INTENT(INOUT) :: logchldata   ! Full set of logchl 
     960      TYPE(obswriinfo), OPTIONAL :: padd            ! Additional info for each variable 
     961      TYPE(obswriinfo), OPTIONAL :: pext            ! Extra info 
     962 
     963      !! * Local declarations  
     964      TYPE(obfbdata) :: fbdata 
     965      CHARACTER(LEN=40) :: cfname             ! netCDF filename 
     966      CHARACTER(LEN=12), PARAMETER :: cpname = 'obs_wri_logchl' 
     967      INTEGER :: jo 
     968      INTEGER :: ja 
     969      INTEGER :: je 
     970      INTEGER :: nadd 
     971      INTEGER :: next 
     972 
     973      IF ( PRESENT( padd ) ) THEN 
     974         nadd = padd%inum 
     975      ELSE 
     976         nadd = 0 
     977      ENDIF 
     978 
     979      IF ( PRESENT( pext ) ) THEN 
     980         next = pext%inum 
     981      ELSE 
     982         next = 0 
     983      ENDIF 
     984 
     985      CALL init_obfbdata( fbdata ) 
     986 
     987      CALL alloc_obfbdata( fbdata, 1, logchldata%nsurf, 1, & 
     988         &                 1 + nadd, next, .TRUE. ) 
     989 
     990      fbdata%cname(1)      = 'LOGCHL' 
     991      fbdata%coblong(1)    = 'logchl concentration' 
     992      fbdata%cobunit(1)    = 'mg/m3' 
     993      DO je = 1, next 
     994         fbdata%cextname(je) = pext%cdname(je) 
     995         fbdata%cextlong(je) = pext%cdlong(je,1) 
     996         fbdata%cextunit(je) = pext%cdunit(je,1) 
     997      END DO 
     998      fbdata%caddname(1)   = 'Hx' 
     999      fbdata%caddlong(1,1) = 'Model interpolated LOGCHL' 
     1000      fbdata%caddunit(1,1) = 'mg/m3' 
     1001      fbdata%cgrid(1)      = 'T' 
     1002      DO ja = 1, nadd 
     1003         fbdata%caddname(1+ja) = padd%cdname(ja) 
     1004         fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1) 
     1005         fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1) 
     1006      END DO 
     1007 
     1008      WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc 
     1009 
     1010      IF(lwp) THEN 
     1011         WRITE(numout,*) 
     1012         WRITE(numout,*)'obs_wri_logchl :' 
     1013         WRITE(numout,*)'~~~~~~~~~~~~~~~~' 
     1014         WRITE(numout,*)'Writing logchl feedback file : ',TRIM(cfname) 
     1015      ENDIF 
     1016 
     1017      ! Transform obs_prof data structure into obfbdata structure 
     1018      fbdata%cdjuldref = '19500101000000' 
     1019      DO jo = 1, logchldata%nsurf 
     1020         fbdata%plam(jo)      = logchldata%rlam(jo) 
     1021         fbdata%pphi(jo)      = logchldata%rphi(jo) 
     1022         WRITE(fbdata%cdtyp(jo),'(I4)') logchldata%ntyp(jo) 
     1023         fbdata%ivqc(jo,:)    = 0 
     1024         fbdata%ivqcf(:,jo,:) = 0 
     1025         IF ( logchldata%nqc(jo) > 10 ) THEN 
     1026            fbdata%ioqc(jo)    = 4 
     1027            fbdata%ioqcf(1,jo) = 0 
     1028            fbdata%ioqcf(2,jo) = logchldata%nqc(jo) - 10 
     1029         ELSE 
     1030            fbdata%ioqc(jo)    = MAX(logchldata%nqc(jo),1) 
     1031            fbdata%ioqcf(:,jo) = 0 
     1032         ENDIF 
     1033         fbdata%ipqc(jo)      = 0 
     1034         fbdata%ipqcf(:,jo)   = 0 
     1035         fbdata%itqc(jo)      = 0 
     1036         fbdata%itqcf(:,jo)   = 0 
     1037         fbdata%cdwmo(jo)     = '' 
     1038         fbdata%kindex(jo)    = logchldata%nsfil(jo) 
     1039         IF (ln_grid_global) THEN 
     1040            fbdata%iobsi(jo,1) = logchldata%mi(jo) 
     1041            fbdata%iobsj(jo,1) = logchldata%mj(jo) 
     1042         ELSE 
     1043            fbdata%iobsi(jo,1) = mig(logchldata%mi(jo)) 
     1044            fbdata%iobsj(jo,1) = mjg(logchldata%mj(jo)) 
     1045         ENDIF 
     1046         CALL greg2jul( 0, & 
     1047            &           logchldata%nmin(jo), & 
     1048            &           logchldata%nhou(jo), & 
     1049            &           logchldata%nday(jo), & 
     1050            &           logchldata%nmon(jo), & 
     1051            &           logchldata%nyea(jo), & 
     1052            &           fbdata%ptim(jo),   & 
     1053            &           krefdate = 19500101 ) 
     1054         fbdata%padd(1,jo,1,1) = logchldata%rmod(jo,1) 
     1055         fbdata%pob(1,jo,1)    = logchldata%robs(jo,1) 
     1056         fbdata%pdep(1,jo)     = 0.0 
     1057         fbdata%idqc(1,jo)     = 0 
     1058         fbdata%idqcf(:,1,jo)  = 0 
     1059         IF ( logchldata%nqc(jo) > 10 ) THEN 
     1060            fbdata%ivlqc(1,jo,1) = 4 
     1061            fbdata%ivlqcf(1,1,jo,1) = 0 
     1062            fbdata%ivlqcf(2,1,jo,1) = logchldata%nqc(jo) - 10 
     1063         ELSE 
     1064            fbdata%ivlqc(1,jo,1) = MAX(logchldata%nqc(jo),1) 
     1065            fbdata%ivlqcf(:,1,jo,1) = 0 
     1066         ENDIF 
     1067         fbdata%iobsk(1,jo,1)  = 0 
     1068         DO ja = 1, nadd 
     1069            fbdata%padd(1,jo,1+ja,1) = & 
     1070               & logchldata%rext(jo,padd%ipoint(ja)) 
     1071         END DO 
     1072         DO je = 1, next 
     1073            fbdata%pext(1,jo,je) = & 
     1074               & logchldata%rext(jo,pext%ipoint(je)) 
     1075         END DO 
     1076 
     1077      END DO 
     1078 
     1079      ! Write the obfbdata structure 
     1080      CALL write_obfbdata( cfname, fbdata ) 
     1081       
     1082      ! Output some basic statistics 
     1083      CALL obs_wri_stats( fbdata ) 
     1084 
     1085      CALL dealloc_obfbdata( fbdata ) 
     1086 
     1087   END SUBROUTINE obs_wri_logchl 
     1088 
     1089   SUBROUTINE obs_wri_spm( cprefix, spmdata, padd, pext ) 
     1090      !!----------------------------------------------------------------------- 
     1091      !! 
     1092      !!                     *** ROUTINE obs_wri_spm  *** 
     1093      !! 
     1094      !! ** Purpose : Write spm observation diagnostics 
     1095      !!              related  
     1096      !! 
     1097      !! ** Method  : NetCDF 
     1098      !!  
     1099      !! ** Action  : 
     1100      !! 
     1101      !!----------------------------------------------------------------------- 
     1102 
     1103      !! * Modules used 
     1104      IMPLICIT NONE 
     1105 
     1106      !! * Arguments 
     1107      CHARACTER(LEN=*), INTENT(IN) :: cprefix       ! Prefix for output files 
     1108      TYPE(obs_surf), INTENT(INOUT) :: spmdata   ! Full set of spm 
     1109      TYPE(obswriinfo), OPTIONAL :: padd            ! Additional info for each variable 
     1110      TYPE(obswriinfo), OPTIONAL :: pext            ! Extra info 
     1111 
     1112      !! * Local declarations  
     1113      TYPE(obfbdata) :: fbdata 
     1114      CHARACTER(LEN=40) :: cfname             ! netCDF filename 
     1115      CHARACTER(LEN=12), PARAMETER :: cpname = 'obs_wri_spm' 
     1116      INTEGER :: jo 
     1117      INTEGER :: ja 
     1118      INTEGER :: je 
     1119      INTEGER :: nadd 
     1120      INTEGER :: next 
     1121 
     1122      IF ( PRESENT( padd ) ) THEN 
     1123         nadd = padd%inum 
     1124      ELSE 
     1125         nadd = 0 
     1126      ENDIF 
     1127 
     1128      IF ( PRESENT( pext ) ) THEN 
     1129         next = pext%inum 
     1130      ELSE 
     1131         next = 0 
     1132      ENDIF 
     1133 
     1134      CALL init_obfbdata( fbdata ) 
     1135 
     1136      CALL alloc_obfbdata( fbdata, 1, spmdata%nsurf, 1, & 
     1137         &                 1 + nadd, next, .TRUE. ) 
     1138 
     1139      fbdata%cname(1)      = 'spm' 
     1140      fbdata%coblong(1)    = 'spm' 
     1141      fbdata%cobunit(1)    = 'g/m3' 
     1142      DO je = 1, next 
     1143         fbdata%cextname(je) = pext%cdname(je) 
     1144         fbdata%cextlong(je) = pext%cdlong(je,1) 
     1145         fbdata%cextunit(je) = pext%cdunit(je,1) 
     1146      END DO 
     1147      fbdata%caddname(1)   = 'Hx' 
     1148      fbdata%caddlong(1,1) = 'Model interpolated spm' 
     1149      fbdata%caddunit(1,1) = 'g/m3' 
     1150      fbdata%cgrid(1)      = 'T' 
     1151      DO ja = 1, nadd 
     1152         fbdata%caddname(1+ja) = padd%cdname(ja) 
     1153         fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1) 
     1154         fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1) 
     1155      END DO 
     1156 
     1157      WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc 
     1158 
     1159      IF(lwp) THEN 
     1160         WRITE(numout,*) 
     1161         WRITE(numout,*)'obs_wri_spm :' 
     1162         WRITE(numout,*)'~~~~~~~~~~~~~~~~' 
     1163         WRITE(numout,*)'Writing spm feedback file : ',TRIM(cfname) 
     1164      ENDIF 
     1165 
     1166      ! Transform obs_prof data structure into obfbdata structure 
     1167      fbdata%cdjuldref = '19500101000000' 
     1168      DO jo = 1, spmdata%nsurf 
     1169         fbdata%plam(jo)      = spmdata%rlam(jo) 
     1170         fbdata%pphi(jo)      = spmdata%rphi(jo) 
     1171         WRITE(fbdata%cdtyp(jo),'(I4)') spmdata%ntyp(jo) 
     1172         fbdata%ivqc(jo,:)    = 0 
     1173         fbdata%ivqcf(:,jo,:) = 0 
     1174         IF ( spmdata%nqc(jo) > 10 ) THEN 
     1175            fbdata%ioqc(jo)    = 4 
     1176            fbdata%ioqcf(1,jo) = 0 
     1177            fbdata%ioqcf(2,jo) = spmdata%nqc(jo) - 10 
     1178         ELSE 
     1179            fbdata%ioqc(jo)    = MAX(spmdata%nqc(jo),1) 
     1180            fbdata%ioqcf(:,jo) = 0 
     1181         ENDIF 
     1182         fbdata%ipqc(jo)      = 0 
     1183         fbdata%ipqcf(:,jo)   = 0 
     1184         fbdata%itqc(jo)      = 0 
     1185         fbdata%itqcf(:,jo)   = 0 
     1186         fbdata%cdwmo(jo)     = '' 
     1187         fbdata%kindex(jo)    = spmdata%nsfil(jo) 
     1188         IF (ln_grid_global) THEN 
     1189            fbdata%iobsi(jo,1) = spmdata%mi(jo) 
     1190            fbdata%iobsj(jo,1) = spmdata%mj(jo) 
     1191         ELSE 
     1192            fbdata%iobsi(jo,1) = mig(spmdata%mi(jo)) 
     1193            fbdata%iobsj(jo,1) = mjg(spmdata%mj(jo)) 
     1194         ENDIF 
     1195         CALL greg2jul( 0, & 
     1196            &           spmdata%nmin(jo), & 
     1197            &           spmdata%nhou(jo), & 
     1198            &           spmdata%nday(jo), & 
     1199            &           spmdata%nmon(jo), & 
     1200            &           spmdata%nyea(jo), & 
     1201            &           fbdata%ptim(jo),   & 
     1202            &           krefdate = 19500101 ) 
     1203         fbdata%padd(1,jo,1,1) = spmdata%rmod(jo,1) 
     1204         fbdata%pob(1,jo,1)    = spmdata%robs(jo,1) 
     1205         fbdata%pdep(1,jo)     = 0.0 
     1206         fbdata%idqc(1,jo)     = 0 
     1207         fbdata%idqcf(:,1,jo)  = 0 
     1208         IF ( spmdata%nqc(jo) > 10 ) THEN 
     1209            fbdata%ivlqc(1,jo,1) = 4 
     1210            fbdata%ivlqcf(1,1,jo,1) = 0 
     1211            fbdata%ivlqcf(2,1,jo,1) = spmdata%nqc(jo) - 10 
     1212         ELSE 
     1213            fbdata%ivlqc(1,jo,1) = MAX(spmdata%nqc(jo),1) 
     1214            fbdata%ivlqcf(:,1,jo,1) = 0 
     1215         ENDIF 
     1216         fbdata%iobsk(1,jo,1)  = 0 
     1217         DO ja = 1, nadd 
     1218            fbdata%padd(1,jo,1+ja,1) = & 
     1219               & spmdata%rext(jo,padd%ipoint(ja)) 
     1220         END DO 
     1221         DO je = 1, next 
     1222            fbdata%pext(1,jo,je) = & 
     1223               & spmdata%rext(jo,pext%ipoint(je)) 
     1224         END DO 
     1225 
     1226      END DO 
     1227 
     1228      ! Write the obfbdata structure 
     1229      CALL write_obfbdata( cfname, fbdata ) 
     1230       
     1231      ! Output some basic statistics 
     1232      CALL obs_wri_stats( fbdata ) 
     1233 
     1234      CALL dealloc_obfbdata( fbdata ) 
     1235 
     1236   END SUBROUTINE obs_wri_spm 
     1237 
     1238   SUBROUTINE obs_wri_fco2( cprefix, fco2data, padd, pext ) 
     1239      !!----------------------------------------------------------------------- 
     1240      !! 
     1241      !!                     *** ROUTINE obs_wri_fco2  *** 
     1242      !! 
     1243      !! ** Purpose : Write fco2 observation diagnostics 
     1244      !!              related  
     1245      !! 
     1246      !! ** Method  : NetCDF 
     1247      !!  
     1248      !! ** Action  : 
     1249      !! 
     1250      !!----------------------------------------------------------------------- 
     1251 
     1252      !! * Modules used 
     1253      IMPLICIT NONE 
     1254 
     1255      !! * Arguments 
     1256      CHARACTER(LEN=*), INTENT(IN) :: cprefix       ! Prefix for output files 
     1257      TYPE(obs_surf), INTENT(INOUT) :: fco2data   ! Full set of fco2 
     1258      TYPE(obswriinfo), OPTIONAL :: padd            ! Additional info for each variable 
     1259      TYPE(obswriinfo), OPTIONAL :: pext            ! Extra info 
     1260 
     1261      !! * Local declarations  
     1262      TYPE(obfbdata) :: fbdata 
     1263      CHARACTER(LEN=40) :: cfname             ! netCDF filename 
     1264      CHARACTER(LEN=12), PARAMETER :: cpname = 'obs_wri_fco2' 
     1265      INTEGER :: jo 
     1266      INTEGER :: ja 
     1267      INTEGER :: je 
     1268      INTEGER :: nadd 
     1269      INTEGER :: next 
     1270 
     1271      IF ( PRESENT( padd ) ) THEN 
     1272         nadd = padd%inum 
     1273      ELSE 
     1274         nadd = 0 
     1275      ENDIF 
     1276 
     1277      IF ( PRESENT( pext ) ) THEN 
     1278         next = pext%inum 
     1279      ELSE 
     1280         next = 0 
     1281      ENDIF 
     1282 
     1283      CALL init_obfbdata( fbdata ) 
     1284 
     1285      CALL alloc_obfbdata( fbdata, 1, fco2data%nsurf, 1, & 
     1286         &                 1 + nadd, next, .TRUE. ) 
     1287 
     1288      fbdata%cname(1)      = 'fco2' 
     1289      fbdata%coblong(1)    = 'fco2' 
     1290      fbdata%cobunit(1)    = 'uatm' 
     1291      DO je = 1, next 
     1292         fbdata%cextname(je) = pext%cdname(je) 
     1293         fbdata%cextlong(je) = pext%cdlong(je,1) 
     1294         fbdata%cextunit(je) = pext%cdunit(je,1) 
     1295      END DO 
     1296      fbdata%caddname(1)   = 'Hx' 
     1297      fbdata%caddlong(1,1) = 'Model interpolated fco2' 
     1298      fbdata%caddunit(1,1) = 'uatm' 
     1299      fbdata%cgrid(1)      = 'T' 
     1300      DO ja = 1, nadd 
     1301         fbdata%caddname(1+ja) = padd%cdname(ja) 
     1302         fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1) 
     1303         fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1) 
     1304      END DO 
     1305 
     1306      WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc 
     1307 
     1308      IF(lwp) THEN 
     1309         WRITE(numout,*) 
     1310         WRITE(numout,*)'obs_wri_fco2 :' 
     1311         WRITE(numout,*)'~~~~~~~~~~~~~~~~' 
     1312         WRITE(numout,*)'Writing fco2 feedback file : ',TRIM(cfname) 
     1313      ENDIF 
     1314 
     1315      ! Transform obs_prof data structure into obfbdata structure 
     1316      fbdata%cdjuldref = '19500101000000' 
     1317      DO jo = 1, fco2data%nsurf 
     1318         fbdata%plam(jo)      = fco2data%rlam(jo) 
     1319         fbdata%pphi(jo)      = fco2data%rphi(jo) 
     1320         WRITE(fbdata%cdtyp(jo),'(I4)') fco2data%ntyp(jo) 
     1321         fbdata%ivqc(jo,:)    = 0 
     1322         fbdata%ivqcf(:,jo,:) = 0 
     1323         IF ( fco2data%nqc(jo) > 10 ) THEN 
     1324            fbdata%ioqc(jo)    = 4 
     1325            fbdata%ioqcf(1,jo) = 0 
     1326            fbdata%ioqcf(2,jo) = fco2data%nqc(jo) - 10 
     1327         ELSE 
     1328            fbdata%ioqc(jo)    = MAX(fco2data%nqc(jo),1) 
     1329            fbdata%ioqcf(:,jo) = 0 
     1330         ENDIF 
     1331         fbdata%ipqc(jo)      = 0 
     1332         fbdata%ipqcf(:,jo)   = 0 
     1333         fbdata%itqc(jo)      = 0 
     1334         fbdata%itqcf(:,jo)   = 0 
     1335         fbdata%cdwmo(jo)     = '' 
     1336         fbdata%kindex(jo)    = fco2data%nsfil(jo) 
     1337         IF (ln_grid_global) THEN 
     1338            fbdata%iobsi(jo,1) = fco2data%mi(jo) 
     1339            fbdata%iobsj(jo,1) = fco2data%mj(jo) 
     1340         ELSE 
     1341            fbdata%iobsi(jo,1) = mig(fco2data%mi(jo)) 
     1342            fbdata%iobsj(jo,1) = mjg(fco2data%mj(jo)) 
     1343         ENDIF 
     1344         CALL greg2jul( 0, & 
     1345            &           fco2data%nmin(jo), & 
     1346            &           fco2data%nhou(jo), & 
     1347            &           fco2data%nday(jo), & 
     1348            &           fco2data%nmon(jo), & 
     1349            &           fco2data%nyea(jo), & 
     1350            &           fbdata%ptim(jo),   & 
     1351            &           krefdate = 19500101 ) 
     1352         fbdata%padd(1,jo,1,1) = fco2data%rmod(jo,1) 
     1353         fbdata%pob(1,jo,1)    = fco2data%robs(jo,1) 
     1354         fbdata%pdep(1,jo)     = 0.0 
     1355         fbdata%idqc(1,jo)     = 0 
     1356         fbdata%idqcf(:,1,jo)  = 0 
     1357         IF ( fco2data%nqc(jo) > 10 ) THEN 
     1358            fbdata%ivlqc(1,jo,1) = 4 
     1359            fbdata%ivlqcf(1,1,jo,1) = 0 
     1360            fbdata%ivlqcf(2,1,jo,1) = fco2data%nqc(jo) - 10 
     1361         ELSE 
     1362            fbdata%ivlqc(1,jo,1) = MAX(fco2data%nqc(jo),1) 
     1363            fbdata%ivlqcf(:,1,jo,1) = 0 
     1364         ENDIF 
     1365         fbdata%iobsk(1,jo,1)  = 0 
     1366         DO ja = 1, nadd 
     1367            fbdata%padd(1,jo,1+ja,1) = & 
     1368               & fco2data%rext(jo,padd%ipoint(ja)) 
     1369         END DO 
     1370         DO je = 1, next 
     1371            fbdata%pext(1,jo,je) = & 
     1372               & fco2data%rext(jo,pext%ipoint(je)) 
     1373         END DO 
     1374 
     1375      END DO 
     1376 
     1377      ! Write the obfbdata structure 
     1378      CALL write_obfbdata( cfname, fbdata ) 
     1379       
     1380      ! Output some basic statistics 
     1381      CALL obs_wri_stats( fbdata ) 
     1382 
     1383      CALL dealloc_obfbdata( fbdata ) 
     1384 
     1385   END SUBROUTINE obs_wri_fco2 
     1386 
     1387   SUBROUTINE obs_wri_pco2( cprefix, pco2data, padd, pext ) 
     1388      !!----------------------------------------------------------------------- 
     1389      !! 
     1390      !!                     *** ROUTINE obs_wri_pco2  *** 
     1391      !! 
     1392      !! ** Purpose : Write pco2 observation diagnostics 
     1393      !!              related  
     1394      !! 
     1395      !! ** Method  : NetCDF 
     1396      !!  
     1397      !! ** Action  : 
     1398      !! 
     1399      !!----------------------------------------------------------------------- 
     1400 
     1401      !! * Modules used 
     1402      IMPLICIT NONE 
     1403 
     1404      !! * Arguments 
     1405      CHARACTER(LEN=*), INTENT(IN) :: cprefix       ! Prefix for output files 
     1406      TYPE(obs_surf), INTENT(INOUT) :: pco2data   ! Full set of pco2 
     1407      TYPE(obswriinfo), OPTIONAL :: padd            ! Additional info for each variable 
     1408      TYPE(obswriinfo), OPTIONAL :: pext            ! Extra info 
     1409 
     1410      !! * Local declarations  
     1411      TYPE(obfbdata) :: fbdata 
     1412      CHARACTER(LEN=40) :: cfname             ! netCDF filename 
     1413      CHARACTER(LEN=12), PARAMETER :: cpname = 'obs_wri_pco2' 
     1414      INTEGER :: jo 
     1415      INTEGER :: ja 
     1416      INTEGER :: je 
     1417      INTEGER :: nadd 
     1418      INTEGER :: next 
     1419 
     1420      IF ( PRESENT( padd ) ) THEN 
     1421         nadd = padd%inum 
     1422      ELSE 
     1423         nadd = 0 
     1424      ENDIF 
     1425 
     1426      IF ( PRESENT( pext ) ) THEN 
     1427         next = pext%inum 
     1428      ELSE 
     1429         next = 0 
     1430      ENDIF 
     1431 
     1432      CALL init_obfbdata( fbdata ) 
     1433 
     1434      CALL alloc_obfbdata( fbdata, 1, pco2data%nsurf, 1, & 
     1435         &                 1 + nadd, next, .TRUE. ) 
     1436 
     1437      fbdata%cname(1)      = 'pco2' 
     1438      fbdata%coblong(1)    = 'pco2' 
     1439      fbdata%cobunit(1)    = 'uatm' 
     1440      DO je = 1, next 
     1441         fbdata%cextname(je) = pext%cdname(je) 
     1442         fbdata%cextlong(je) = pext%cdlong(je,1) 
     1443         fbdata%cextunit(je) = pext%cdunit(je,1) 
     1444      END DO 
     1445      fbdata%caddname(1)   = 'Hx' 
     1446      fbdata%caddlong(1,1) = 'Model interpolated pco2' 
     1447      fbdata%caddunit(1,1) = 'uatm' 
     1448      fbdata%cgrid(1)      = 'T' 
     1449      DO ja = 1, nadd 
     1450         fbdata%caddname(1+ja) = padd%cdname(ja) 
     1451         fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1) 
     1452         fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1) 
     1453      END DO 
     1454 
     1455      WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc 
     1456 
     1457      IF(lwp) THEN 
     1458         WRITE(numout,*) 
     1459         WRITE(numout,*)'obs_wri_pco2 :' 
     1460         WRITE(numout,*)'~~~~~~~~~~~~~~~~' 
     1461         WRITE(numout,*)'Writing pco2 feedback file : ',TRIM(cfname) 
     1462      ENDIF 
     1463 
     1464      ! Transform obs_prof data structure into obfbdata structure 
     1465      fbdata%cdjuldref = '19500101000000' 
     1466      DO jo = 1, pco2data%nsurf 
     1467         fbdata%plam(jo)      = pco2data%rlam(jo) 
     1468         fbdata%pphi(jo)      = pco2data%rphi(jo) 
     1469         WRITE(fbdata%cdtyp(jo),'(I4)') pco2data%ntyp(jo) 
     1470         fbdata%ivqc(jo,:)    = 0 
     1471         fbdata%ivqcf(:,jo,:) = 0 
     1472         IF ( pco2data%nqc(jo) > 10 ) THEN 
     1473            fbdata%ioqc(jo)    = 4 
     1474            fbdata%ioqcf(1,jo) = 0 
     1475            fbdata%ioqcf(2,jo) = pco2data%nqc(jo) - 10 
     1476         ELSE 
     1477            fbdata%ioqc(jo)    = MAX(pco2data%nqc(jo),1) 
     1478            fbdata%ioqcf(:,jo) = 0 
     1479         ENDIF 
     1480         fbdata%ipqc(jo)      = 0 
     1481         fbdata%ipqcf(:,jo)   = 0 
     1482         fbdata%itqc(jo)      = 0 
     1483         fbdata%itqcf(:,jo)   = 0 
     1484         fbdata%cdwmo(jo)     = '' 
     1485         fbdata%kindex(jo)    = pco2data%nsfil(jo) 
     1486         IF (ln_grid_global) THEN 
     1487            fbdata%iobsi(jo,1) = pco2data%mi(jo) 
     1488            fbdata%iobsj(jo,1) = pco2data%mj(jo) 
     1489         ELSE 
     1490            fbdata%iobsi(jo,1) = mig(pco2data%mi(jo)) 
     1491            fbdata%iobsj(jo,1) = mjg(pco2data%mj(jo)) 
     1492         ENDIF 
     1493         CALL greg2jul( 0, & 
     1494            &           pco2data%nmin(jo), & 
     1495            &           pco2data%nhou(jo), & 
     1496            &           pco2data%nday(jo), & 
     1497            &           pco2data%nmon(jo), & 
     1498            &           pco2data%nyea(jo), & 
     1499            &           fbdata%ptim(jo),   & 
     1500            &           krefdate = 19500101 ) 
     1501         fbdata%padd(1,jo,1,1) = pco2data%rmod(jo,1) 
     1502         fbdata%pob(1,jo,1)    = pco2data%robs(jo,1) 
     1503         fbdata%pdep(1,jo)     = 0.0 
     1504         fbdata%idqc(1,jo)     = 0 
     1505         fbdata%idqcf(:,1,jo)  = 0 
     1506         IF ( pco2data%nqc(jo) > 10 ) THEN 
     1507            fbdata%ivlqc(1,jo,1) = 4 
     1508            fbdata%ivlqcf(1,1,jo,1) = 0 
     1509            fbdata%ivlqcf(2,1,jo,1) = pco2data%nqc(jo) - 10 
     1510         ELSE 
     1511            fbdata%ivlqc(1,jo,1) = MAX(pco2data%nqc(jo),1) 
     1512            fbdata%ivlqcf(:,1,jo,1) = 0 
     1513         ENDIF 
     1514         fbdata%iobsk(1,jo,1)  = 0 
     1515         DO ja = 1, nadd 
     1516            fbdata%padd(1,jo,1+ja,1) = & 
     1517               & pco2data%rext(jo,padd%ipoint(ja)) 
     1518         END DO 
     1519         DO je = 1, next 
     1520            fbdata%pext(1,jo,je) = & 
     1521               & pco2data%rext(jo,pext%ipoint(je)) 
     1522         END DO 
     1523 
     1524      END DO 
     1525 
     1526      ! Write the obfbdata structure 
     1527      CALL write_obfbdata( cfname, fbdata ) 
     1528       
     1529      ! Output some basic statistics 
     1530      CALL obs_wri_stats( fbdata ) 
     1531 
     1532      CALL dealloc_obfbdata( fbdata ) 
     1533 
     1534   END SUBROUTINE obs_wri_pco2 
     1535 
    9321536   SUBROUTINE obs_wri_stats( fbdata ) 
    9331537      !!----------------------------------------------------------------------- 
Note: See TracChangeset for help on using the changeset viewer.