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.
c4comb.F90 in branches/2014/dev_r4650_UKMO12_CFL_diags_take2/NEMOGCM/TOOLS/OBSTOOLS/src – NEMO

source: branches/2014/dev_r4650_UKMO12_CFL_diags_take2/NEMOGCM/TOOLS/OBSTOOLS/src/c4comb.F90 @ 5947

Last change on this file since 5947 was 5947, checked in by timgraham, 8 years ago

Reinstate svn Id keywords before merge

  • Property svn:keywords set to Id
File size: 37.0 KB
Line 
1PROGRAM c4comb
2   !!---------------------------------------------------------------------
3   !!
4   !!                     ** PROGRAM c4comb **
5   !!
6   !!  ** Purpose : Combine MPI decomposed class4 files into one file
7   !!
8   !!  ** Method  : Use of utilities from obs_utils, ooo_utils.
9   !!
10   !!  ** Action  :
11   !!
12   !!   Usage:
13   !!     c4comb.exe outputfile inputfile1 inputfile2 ...
14   !!
15   !!   History :
16   !!----------------------------------------------------------------------
17   USE netcdf
18   USE obs_const
19   USE obs_utils
20   USE ooo_utils, ONLY: date_format, obfilldbl
21   USE toolspar_kind
22   IMPLICIT NONE
23   !! Command line setup
24#ifndef NOIARGCPROTO
25   INTEGER,EXTERNAL :: iargc
26#endif
27   INTEGER :: nargs,    & !: number of command line arguments
28            & ia,       & !: argument loop index
29            & ninfiles    !: number of input files
30   !! Routine arguments
31   CHARACTER(len=256) :: cdoutfile
32   CHARACTER(len=256),ALLOCATABLE :: cdinfile(:)
33   !! Routine variables
34   CHARACTER(len=80) :: cpname
35   INTEGER,PARAMETER :: nstr=8, n128=128
36   INTEGER           ::    ncid,   & !: netcdf file id
37                        &  dimid,  & !: netcdf dimension id
38                        &  dpdim,  & !: netcdf dimension ids
39                        &  fcdim,  &
40                        &  vrdim,  &
41                        &  obdim,  &
42                        &  stdim,  &
43                        &  sxdim,  &
44                        &  fdvid,  & !: netcdf variable ids
45                        &  lonid,  &
46                        &  latid,  &
47                        &  depid,  &
48                        &  varid,  &
49                        &  unitid, &
50                        &  obvid,  &
51                        &  fcvid,  &
52                        &  prvid,  &
53                        &  clvid,  &
54                        &  dm2id,  &
55                        &  dm1id,  &
56                        &  mdtid,  &
57                        &  altid,  &
58                        &  qcvid,  &
59                        &  jdvid,  &
60                        &  mjdid,  &
61                        &  typid,  &
62                        &  idvid,  &
63                        &  ndeps,  & !: number depths
64                        &  nfcst,  & !: number forecast
65                        &  nvars,  & !: number variables
66                        &  nobs,   & !: number obs
67                        &  sdeps,  &
68                        &  sobs,   &
69                        &  l_dex,  &
70                        &  u_dex
71
72   INTEGER                   :: iob, idep, istat
73   INTEGER, DIMENSION(2)     :: dim2a, dim2b, dim2c, dim2d
74   INTEGER, DIMENSION(3)     :: dim3a
75   INTEGER, DIMENSION(4)     :: dim4a
76   INTEGER,  ALLOCATABLE, DIMENSION(:)           :: fcday
77   REAL(wp), ALLOCATABLE, DIMENSION(:)           :: modjd
78   !: Global Attributes
79   CHARACTER(len=40)                             :: nam_str, &
80                                                 & version,  &
81                                                 & contact,  &
82                                                 & sys_str,  &
83                                                 & cfg_str,  &
84                                                 & ins_str,  &
85                                                 & val_str,  &
86                                                 & dat_str,  &
87                                                 & obs_str
88   !: Variable Attributes
89   CHARACTER(len=100)                            :: lon_units, &
90                                                 &  lat_units, &
91                                                 &  dep_units, &
92                                                 &  jul_units, &
93                                                 &  mjd_units, &
94                                                 &  fcd_units, &
95                                                 &  lead_comment, &
96                                                 &  fcst_comment, &
97                                                 &  per_comment,  &
98                                                 &  cli_comment,  &
99                                                 &  dm2_comment,  &
100                                                 &  dm1_comment
101   CHARACTER(len=128)                            :: qc_comment,   &
102                                                 &  qc_flag_meaning
103   INTEGER, DIMENSION(2)                         :: qc_flag_value
104
105   !: Global Arrays
106   REAL(wp), ALLOCATABLE, DIMENSION(:)           :: g_lam, &
107                                                 &  g_phi, & 
108                                                 &  gjuld 
109   CHARACTER(len=n128),ALLOCATABLE,DIMENSION(:)  :: gtype
110   CHARACTER(len=nstr),ALLOCATABLE,DIMENSION(:)  ::        &
111                                                 &  g_id,  & 
112                                                 &  gvnam, & 
113                                                 &  gunit
114   REAL(wp), ALLOCATABLE, DIMENSION(:,:)         :: g_dep
115   REAL(wp), ALLOCATABLE, DIMENSION(:,:,:)       :: g3dob, &
116                                                 &  g3dcl, &
117                                                 &  g3mdt, &
118                                                 &  g3alt, &
119                                                 &  g3dm2, &
120                                                 &  g3dm1
121   INTEGER(ik),  ALLOCATABLE, DIMENSION(:,:,:)   :: g3dqc
122   REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:)     :: g3dmc, &
123                                                 &  g3dpr
124   !: Small Arrays 
125   REAL(wp), ALLOCATABLE, DIMENSION(:)       ::     s_lam, &
126                                                 &  s_phi, & 
127                                                 &  sjuld
128
129   CHARACTER(len=n128),ALLOCATABLE,DIMENSION(:) ::  stype
130   CHARACTER(len=nstr),ALLOCATABLE,DIMENSION(:) ::        &
131                                                 &  s_id 
132
133   REAL(wp), ALLOCATABLE, DIMENSION(:,:)        :: s_dep
134   REAL(wp), ALLOCATABLE, DIMENSION(:,:,:)      :: s3dob, &
135                                                &  s3dcl, &
136                                                &  s3mdt, &
137                                                &  s3alt, &
138                                                &  s3dm2, &
139                                                &  s3dm1
140   INTEGER(ik),  ALLOCATABLE, DIMENSION(:,:,:)  :: s3dqc
141   REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:)    :: s3dmc, &
142                                                &  s3dpr
143
144   !: File creation logical
145   LOGICAL                                     :: ln_cre
146
147   !: Optional variable logicals
148   LOGICAL                                     :: ln_init, &
149                                                & ln_mdt,  &
150                                                & ln_altbias, &
151                                                & ln_best
152   !! Command name
153   cpname='c4comb.exe'
154
155   !! Process command line
156   nargs = IARGC()
157   IF (nargs /= 2) THEN
158      WRITE(*, *) "Usage: c4comb.exe outputfile inputfile1 inputfile2 ..."
159      CALL abort()
160   END IF
161   CALL GETARG(1, cdoutfile)
162
163   !! Process input files
164   !! Set output file creation to off
165   ln_cre = .false.
166
167   !! Turn optional variables off
168   ln_init = .false.
169   ln_best = .false.
170   ln_altbias = .false.
171   ln_mdt = .false.
172
173   !! Compute size of output file
174   nobs = 0
175   ndeps= 0
176   ALLOCATE( cdinfile( nargs - 1 ) )
177   ninfiles = nargs - 1
178   DO ia = 1, ninfiles
179      CALL GETARG(ia+1, cdinfile(ia))
180      WRITE(*,*) "Opening : ", TRIM(cdinfile(ia))
181      !! Open Netcdf file
182      istat = nf90_open(TRIM(cdinfile(ia)),nf90_nowrite,ncid)
183      IF (istat == nf90_noerr) THEN
184         !! Turn output file creation on
185         ln_cre = .true.
186         !! Get Dimensions
187         CALL chkerr( nf90_inq_dimid(ncid, 'numobs',  dimid),              cpname, __LINE__ )
188         CALL chkerr( nf90_inquire_dimension(ncid,    dimid, len=sobs  ),  cpname, __LINE__ )
189         CALL chkerr( nf90_inq_dimid(ncid, 'numdeps', dimid),              cpname, __LINE__ )
190         CALL chkerr( nf90_inquire_dimension(ncid,    dimid, len=sdeps ),  cpname, __LINE__ )
191         CALL chkerr( nf90_inq_dimid(ncid, 'numfcsts',dimid),              cpname, __LINE__ )
192         CALL chkerr( nf90_inquire_dimension(ncid,    dimid, len=nfcst ),  cpname, __LINE__ )
193         CALL chkerr( nf90_inq_dimid(ncid, 'numvars', dimid),              cpname, __LINE__ )
194         CALL chkerr( nf90_inquire_dimension(ncid,    dimid, len=nvars ),  cpname, __LINE__ )
195         !! Close Netcdf file
196         CALL chkerr( nf90_close(ncid), cpname, __LINE__ )
197         !! Report on file contents
198         WRITE(*,'(2A)')'File = ', TRIM(cdinfile(ia))
199         WRITE(*,'(A,I9,A)')'has', sobs, ' observations'
200         !! Increment size
201         nobs  = nobs + sobs       !: Accumulate number of profiles
202         ndeps = MAX(ndeps, sdeps) !: Define maximum number of levels needed
203      END IF ! istat
204   END DO
205
206   !! Allocate global arrays
207   ALLOCATE( g_phi(nobs),                      &
208          &  g_lam(nobs),                      &
209          &  g_dep(ndeps, nobs),               &
210          &  g3dob(ndeps,        nvars, nobs), &
211          &  g3dmc(ndeps, nfcst, nvars, nobs), &
212          &  g3dpr(ndeps, nfcst, nvars, nobs), &
213          &  g3dcl(ndeps,        nvars, nobs), &
214          &  g3dm2(ndeps,        nvars, nobs), &
215          &  g3dm1(ndeps,        nvars, nobs), &
216          &  g3mdt(ndeps,        nvars, nobs), &
217          &  g3alt(ndeps,        nvars, nobs), &
218          &  g3dqc(ndeps,        nvars, nobs), &
219          &  gjuld(nobs),                      &
220          &  gtype(nobs),                      &
221          &  g_id(nobs),                       &
222          &  gvnam(nvars),                     &
223          &  gunit(nvars) )
224   ALLOCATE(fcday(nfcst), modjd(nfcst))
225
226   !! Fill with missing data value
227   g_dep(:,:)     = 99999.
228   g3dmc(:,:,:,:) = 99999.
229   g3dpr(:,:,:,:) = 99999.
230   g3dob(:,:,:)   = 99999.
231   g3dcl(:,:,:)   = 99999.
232   g3dm2(:,:,:)   = 99999.
233   g3dm1(:,:,:)   = 99999.
234   g3mdt(:,:,:)   = 99999.
235   g3alt(:,:,:)   = 99999.
236   g3dqc(:,:,:)   = NF90_FILL_SHORT
237
238   !! Read in each file
239   !  initialise global matrix indices
240   l_dex = 0
241   u_dex = 0
242
243   !! initialise Global attribute strings
244   nam_str = ''
245   version = ''
246   contact = ''
247   sys_str = ''
248   cfg_str = ''
249   ins_str = ''
250   val_str = ''
251   dat_str = ''
252   obs_str = ''
253
254   !! initialise Variable attribute strings
255   fcd_units = ''
256   lon_units = ''
257   lat_units = ''
258   dep_units = ''
259   jul_units = ''
260   mjd_units = ''
261   lead_comment = ''
262   fcst_comment = ''
263   per_comment  = ''
264   cli_comment  = ''
265   dm2_comment  = ''
266   dm1_comment  = ''
267   qc_comment  = ''
268   qc_flag_meaning  = ''
269
270   DO ia = 1, ninfiles
271      WRITE(*,*) "Opening : ", TRIM(cdinfile(ia))
272      !! Open Netcdf file
273      istat = nf90_open(TRIM(cdinfile(ia)),nf90_nowrite,ncid)
274      IF (istat == nf90_noerr) THEN
275         !! Get Global Attributes
276         CALL chkerr( nf90_get_att(ncid, nf90_global,'title',         nam_str),cpname, __LINE__)
277         CALL chkerr( nf90_get_att(ncid, nf90_global,'version',       version),cpname, __LINE__)
278         CALL chkerr( nf90_get_att(ncid, nf90_global,'contact',       contact),cpname, __LINE__)
279         CALL chkerr( nf90_get_att(ncid, nf90_global,'obs_type',      obs_str),cpname, __LINE__)
280         CALL chkerr( nf90_get_att(ncid, nf90_global,'system',        sys_str),cpname, __LINE__)
281         CALL chkerr( nf90_get_att(ncid, nf90_global,'configuration', cfg_str),cpname, __LINE__)
282         CALL chkerr( nf90_get_att(ncid, nf90_global,'institution',   ins_str),cpname, __LINE__)
283         CALL chkerr( nf90_get_att(ncid, nf90_global,'validity_time', val_str),cpname, __LINE__)
284         !! Get Dimensions of single file
285         CALL chkerr( nf90_inq_dimid(ncid, 'numdeps', dimid),              cpname, __LINE__ )
286         CALL chkerr( nf90_inquire_dimension(ncid,    dimid, len=sdeps ),  cpname, __LINE__ )
287         CALL chkerr( nf90_inq_dimid(ncid, 'numfcsts',dimid),              cpname, __LINE__ )
288         CALL chkerr( nf90_inquire_dimension(ncid,    dimid, len=nfcst ),  cpname, __LINE__ )
289         CALL chkerr( nf90_inq_dimid(ncid, 'numvars', dimid),              cpname, __LINE__ )
290         CALL chkerr( nf90_inquire_dimension(ncid,    dimid, len=nvars ),  cpname, __LINE__ )
291         CALL chkerr( nf90_inq_dimid(ncid, 'numobs',  dimid),              cpname, __LINE__ )
292         CALL chkerr( nf90_inquire_dimension(ncid,    dimid, len=sobs  ),  cpname, __LINE__ )
293         !! Check for Optional variables in first file
294         IF (ia == 1) THEN
295            !! Best estimate
296            istat = nf90_inq_varid(ncid,'best_estimate',dm2id)
297            IF (istat == nf90_noerr) THEN
298               ln_best = .TRUE.
299            ENDIF 
300            !! nrt_analysis
301            istat = nf90_inq_varid(ncid,'nrt_analysis',dm1id)
302            IF (istat == nf90_noerr) THEN
303               ln_init = .TRUE.
304            ENDIF 
305            !! Mean Dynamic Topography
306            istat = nf90_inq_varid(ncid,'mdt_reference',mdtid)
307            IF (istat == nf90_noerr) THEN
308               ln_mdt = .TRUE.
309            ENDIF 
310            !! Altimeter bias
311            istat = nf90_inq_varid(ncid,'altimeter_bias',altid)
312            IF (istat == nf90_noerr) THEN
313               ln_altbias = .TRUE.
314            ENDIF
315         END IF
316         WRITE(*,*) TRIM(cdinfile(ia)), " contains ", sobs, " observations"
317         WRITE(*,*) TRIM(cdinfile(ia)), " contains ", sdeps, " depths"
318         WRITE(*,*) TRIM(cdinfile(ia)), " contains ", nfcst, " forecasts"
319         WRITE(*,*) TRIM(cdinfile(ia)), " contains ", nvars, " vars"
320         !! Read Variables
321         IF (sobs /= 0) THEN
322            !! Get Variable ids
323            CALL chkerr(nf90_inq_varid(ncid,'leadtime', fdvid) ,cpname, __LINE__ )
324            CALL chkerr(nf90_inq_varid(ncid,'longitude',   lonid) ,cpname, __LINE__ )
325            CALL chkerr(nf90_inq_varid(ncid,'latitude',    latid) ,cpname, __LINE__ )
326            CALL chkerr(nf90_inq_varid(ncid,'depth',       depid) ,cpname, __LINE__ )
327            CALL chkerr(nf90_inq_varid(ncid,'varname',     varid) ,cpname, __LINE__ )
328            CALL chkerr(nf90_inq_varid(ncid,'unitname',   unitid) ,cpname, __LINE__ )
329            CALL chkerr(nf90_inq_varid(ncid,'observation', obvid) ,cpname, __LINE__ )
330            CALL chkerr(nf90_inq_varid(ncid,'forecast',    fcvid) ,cpname, __LINE__ )
331            CALL chkerr(nf90_inq_varid(ncid,'persistence', prvid) ,cpname, __LINE__ )
332            CALL chkerr(nf90_inq_varid(ncid,'climatology', clvid) ,cpname, __LINE__ )
333            CALL chkerr(nf90_inq_varid(ncid,'qc',          qcvid) ,cpname, __LINE__ )
334            CALL chkerr(nf90_inq_varid(ncid,'juld',        jdvid) ,cpname, __LINE__ )
335            CALL chkerr(nf90_inq_varid(ncid,'modeljuld',   mjdid) ,cpname, __LINE__ )
336            CALL chkerr(nf90_inq_varid(ncid,'type',        typid) ,cpname, __LINE__ )
337            CALL chkerr(nf90_inq_varid(ncid,'id',          idvid) ,cpname, __LINE__ )
338            !! Get variable attributes
339            CALL chkerr(nf90_get_att(ncid, fdvid, 'units', fcd_units) ,cpname, __LINE__ )
340            CALL chkerr(nf90_get_att(ncid, lonid, 'units', lon_units) ,cpname, __LINE__ )
341            CALL chkerr(nf90_get_att(ncid, latid, 'units', lat_units) ,cpname, __LINE__ )
342            CALL chkerr(nf90_get_att(ncid, depid, 'units', dep_units) ,cpname, __LINE__ )
343            CALL chkerr(nf90_get_att(ncid, jdvid, 'units', jul_units) ,cpname, __LINE__ )
344            CALL chkerr(nf90_get_att(ncid, mjdid, 'units', mjd_units) ,cpname, __LINE__ )
345            CALL chkerr(nf90_get_att(ncid, fcvid, 'comment', fcst_comment) ,cpname, __LINE__ )
346            CALL chkerr(nf90_get_att(ncid, prvid, 'comment', per_comment)  ,cpname, __LINE__ )
347            CALL chkerr(nf90_get_att(ncid, clvid, 'comment', cli_comment)  ,cpname, __LINE__ )
348            CALL chkerr(nf90_get_att(ncid, fdvid, 'comment', lead_comment) ,cpname, __LINE__ )
349            CALL chkerr(nf90_get_att(ncid, qcvid, 'comment', qc_comment) ,cpname, __LINE__ )
350            CALL chkerr(nf90_get_att(ncid, qcvid, 'flag_value', qc_flag_value) ,cpname, __LINE__ )
351            CALL chkerr(nf90_get_att(ncid, qcvid, 'flag_meaning', qc_flag_meaning) ,cpname, __LINE__ )
352            !! Optional variables
353            IF (ln_best) THEN
354               CALL chkerr(nf90_inq_varid(ncid,'best_estimate',dm2id) ,cpname, __LINE__ )
355               CALL chkerr(nf90_get_att(ncid, dm2id, 'comment', dm2_comment)  ,cpname, __LINE__ )
356            ENDIF
357            IF (ln_init) THEN
358               CALL chkerr(nf90_inq_varid(ncid,'nrt_analysis', dm1id) ,cpname, __LINE__ )
359               CALL chkerr(nf90_get_att(ncid, dm1id, 'comment', dm1_comment)  ,cpname, __LINE__ )
360            ENDIF
361            IF (ln_mdt) THEN
362               CALL chkerr(nf90_inq_varid(ncid,'mdt_reference', mdtid) ,cpname, __LINE__ )
363            ENDIF
364            IF (ln_altbias) THEN
365               CALL chkerr(nf90_inq_varid(ncid,'altimeter_bias', altid) ,cpname, __LINE__ )
366            ENDIF
367
368            !! Allocate small arrays
369            ALLOCATE( s_lam(sobs), s_phi(sobs),  s_dep(sdeps, sobs),               &
370                   &  s3dob(sdeps,        nvars, sobs), & !: observations
371                   &  s3dmc(sdeps, nfcst, nvars, sobs), & !: model data
372                   &  s3dpr(sdeps, nfcst, nvars, sobs), & !: persistence
373                   &  s3dcl(sdeps,        nvars, sobs), & !: climatology
374                   &  s3dm2(sdeps,        nvars, sobs), & !: best estimate
375                   &  s3dm1(sdeps,        nvars, sobs), & !: nrt_analysis
376                   &  s3mdt(sdeps,        nvars, sobs), & !: mdt
377                   &  s3alt(sdeps,        nvars, sobs), & !: altbias
378                   &  s3dqc(sdeps,        nvars, sobs), & !: QC
379                   &  sjuld(sobs), stype( sobs),  &
380                   &  s_id(sobs) )
381            !! Fill with missing data value
382            s3dmc(:,:,:,:) = 99999.
383            s3dpr(:,:,:,:) = 99999.
384            s3dob(:,:,:)   = 99999.
385            s3dcl(:,:,:)   = 99999.
386            s3dm2(:,:,:)   = 99999.
387            s3dm1(:,:,:)   = 99999.
388            s3mdt(:,:,:)   = 99999.
389            s3alt(:,:,:)   = 99999.
390            s3dqc(:,:,:)   = NF90_FILL_SHORT
391
392            !! Read variables into small arrays
393            CALL chkerr( nf90_get_var(ncid, fdvid, fcday), cpname, __LINE__ )       
394            CALL chkerr( nf90_get_var(ncid, lonid, s_lam), cpname, __LINE__ )       
395            CALL chkerr( nf90_get_var(ncid, latid, s_phi), cpname, __LINE__ )       
396            CALL chkerr( nf90_get_var(ncid, depid, s_dep), cpname, __LINE__ )       
397            CALL chkerr( nf90_get_var(ncid, obvid, s3dob), cpname, __LINE__ )       
398            CALL chkerr( nf90_get_var(ncid, fcvid, s3dmc), cpname, __LINE__ )       
399            CALL chkerr( nf90_get_var(ncid, prvid, s3dpr), cpname, __LINE__ )       
400            CALL chkerr( nf90_get_var(ncid, clvid, s3dcl), cpname, __LINE__ )       
401            CALL chkerr( nf90_get_var(ncid, qcvid, s3dqc), cpname, __LINE__ )       
402            CALL chkerr( nf90_get_var(ncid, jdvid, sjuld), cpname, __LINE__ )       
403            CALL chkerr( nf90_get_var(ncid, mjdid, modjd), cpname, __LINE__ )       
404            CALL chkerr( nf90_get_var(ncid, typid, stype), cpname, __LINE__ )       
405            CALL chkerr( nf90_get_var(ncid, idvid, s_id),  cpname, __LINE__ )       
406            !! Read unitname and varname into global arrays
407            CALL chkerr( nf90_get_var(ncid, varid, gvnam), cpname, __LINE__ )       
408            CALL chkerr( nf90_get_var(ncid, unitid,gunit), cpname, __LINE__ )       
409            !! Optional variables read
410            IF (ln_best) THEN
411               CALL chkerr( nf90_get_var(ncid, dm2id, s3dm2), cpname, __LINE__ )       
412            ENDIF
413            IF (ln_init) THEN
414               CALL chkerr( nf90_get_var(ncid, dm1id, s3dm1), cpname, __LINE__ )       
415            ENDIF
416            IF (ln_mdt) THEN
417               CALL chkerr( nf90_get_var(ncid, mdtid, s3mdt), cpname, __LINE__ )       
418            ENDIF
419            IF (ln_altbias) THEN
420               CALL chkerr( nf90_get_var(ncid, altid, s3alt), cpname, __LINE__ )       
421            ENDIF
422
423            !! Fill Global arrays
424            !  increment numobs indices
425            l_dex = u_dex + 1
426            u_dex = l_dex + sobs -1
427
428            g_lam(l_dex:u_dex)                          = s_lam(:)
429            g_phi(l_dex:u_dex)                          = s_phi(:)
430            g_dep(1:sdeps,  l_dex:u_dex)                = s_dep(1:sdeps,:)
431            g3dob(1:sdeps,1:nvars,l_dex:u_dex)          = s3dob(1:sdeps,1:nvars,:)
432            g3dmc(1:sdeps,1:nfcst,1:nvars,l_dex:u_dex)  = s3dmc(1:sdeps,1:nfcst,1:nvars,:)
433            g3dpr(1:sdeps,1:nfcst,1:nvars,l_dex:u_dex)  = s3dpr(1:sdeps,1:nfcst,1:nvars,:)
434            g3dcl(1:sdeps,1:nvars,l_dex:u_dex)          = s3dcl(1:sdeps,1:nvars,:)
435            g3dm2(1:sdeps,1:nvars,l_dex:u_dex)          = s3dm2(1:sdeps,1:nvars,:)
436            g3dm1(1:sdeps,1:nvars,l_dex:u_dex)          = s3dm1(1:sdeps,1:nvars,:)
437            g3mdt(1:sdeps,1:nvars,l_dex:u_dex)          = s3mdt(1:sdeps,1:nvars,:)
438            g3alt(1:sdeps,1:nvars,l_dex:u_dex)          = s3alt(1:sdeps,1:nvars,:)
439            g3dqc(1:sdeps,1:nvars,l_dex:u_dex)          = s3dqc(1:sdeps,1:nvars,:)
440            gjuld(l_dex:u_dex)                          = sjuld(:)
441            gtype(l_dex:u_dex)                          = stype(:)
442            g_id(l_dex:u_dex)                           = s_id(:)
443
444            !! Deallocate small array
445            DEALLOCATE( s_lam, s_phi, s_dep, s3dob, s3dmc, s3dpr, s3dcl, s3dqc, s3dm2, s3dm1, s3mdt, s3alt, sjuld, stype, s_id)
446         ENDIF ! sobs   
447         !! Close Netcdf file
448         CALL chkerr( nf90_close(ncid), cpname, __LINE__ )
449      END IF ! istat
450   END DO
451
452   !! Create Output file
453   IF (ln_cre) THEN
454      WRITE(*,*) 'Create the output file, ',trim(cdoutfile)
455      CALL chkerr( nf90_create(trim(cdoutfile),nf90_clobber,ncid),      cpname, __LINE__ )
456      !! Put Global Attributes
457      CALL date_format(dat_str) 
458      CALL chkerr( nf90_put_att(ncid, nf90_global,'title',          trim(nam_str)),cpname, __LINE__)
459      CALL chkerr( nf90_put_att(ncid, nf90_global,'version',        trim(version)),cpname, __LINE__)
460      CALL chkerr( nf90_put_att(ncid, nf90_global,'creation_date',  trim(dat_str)),cpname, __LINE__)
461      CALL chkerr( nf90_put_att(ncid, nf90_global,'contact',        trim(contact)),cpname, __LINE__)
462      CALL chkerr( nf90_put_att(ncid, nf90_global,'obs_type',       trim(obs_str)),cpname, __LINE__)
463      CALL chkerr( nf90_put_att(ncid, nf90_global,'system',         trim(sys_str)),cpname, __LINE__)
464      CALL chkerr( nf90_put_att(ncid, nf90_global,'configuration',  trim(cfg_str)),cpname, __LINE__)
465      CALL chkerr( nf90_put_att(ncid, nf90_global,'institution',    trim(ins_str)),cpname, __LINE__)
466      CALL chkerr( nf90_put_att(ncid, nf90_global,'validity_time',  trim(val_str)),cpname, __LINE__)
467      CALL chkerr( nf90_put_att(ncid, nf90_global,'best_estimate_description', &
468                              &       'analysis produced 2 days behind real time'),cpname, __LINE__)
469      CALL chkerr( nf90_put_att(ncid, nf90_global,'time_interp', 'daily average fields'),cpname, __LINE__)
470      WRITE(*,*) 'Succesfully put global attributes '
471
472      !! Define Dimensions
473      CALL chkerr( nf90_def_dim(ncid, 'numdeps',          ndeps, dpdim) ,cpname, __LINE__ )
474      CALL chkerr( nf90_def_dim(ncid, 'numfcsts',         nfcst, fcdim) ,cpname, __LINE__ )
475      CALL chkerr( nf90_def_dim(ncid, 'numvars',          nvars, vrdim) ,cpname, __LINE__ )
476      CALL chkerr( nf90_def_dim(ncid, 'numobs',           nobs,  obdim) ,cpname, __LINE__ )
477      CALL chkerr( nf90_def_dim(ncid, 'string_length8',   nstr,  stdim) ,cpname, __LINE__ )
478      CALL chkerr( nf90_def_dim(ncid, 'string_length128', n128,  sxdim) ,cpname, __LINE__ )
479      WRITE(*,*) 'Succesfully defined dimensions'
480
481      !! Define possible dimension permutations
482      ! 2d
483      dim2a(:) = (/ dpdim, obdim /) !: (/ ndeps, nobs  /)
484      dim2b(:) = (/ stdim, obdim /) !: (/ nstr,  nobs  /)
485      dim2c(:) = (/ stdim, vrdim /) !: (/ nstr,  nvars /)
486      dim2d(:) = (/ sxdim, obdim /) !: (/ nstr,  nobs  /)
487      ! 3d
488      dim3a(:) = (/ dpdim, vrdim, obdim/) !: (/ ndeps, nvars, nobs /)
489      ! 4d
490      dim4a(:) = (/ dpdim, fcdim, vrdim, obdim /) !: (/ ndeps, nfcst, nvars, nobs /)
491
492
493      !! Create the variables
494      !  Forecast day
495      CALL chkerr( nf90_def_var(ncid, 'leadtime',  nf90_double, fcdim, fdvid)          ,cpname, __LINE__ )
496      CALL chkerr( nf90_put_att(ncid, fdvid, 'long_name', 'Model forecast day offset') ,cpname, __LINE__ )
497      CALL chkerr( nf90_put_att(ncid, fdvid, 'units',     trim(fcd_units))             ,cpname, __LINE__ )
498      CALL chkerr( nf90_put_att(ncid, fdvid, 'comment',   trim(lead_comment))          ,cpname, __LINE__ )
499      WRITE(*,*) 'leadtime created'
500      !  longitude
501      CALL chkerr( nf90_def_var(ncid, 'longitude',    nf90_float, obdim, lonid)        ,cpname, __LINE__ )
502      CALL chkerr( nf90_put_att(ncid, lonid, 'long_name', 'Longitudes')                ,cpname, __LINE__ )
503      CALL chkerr( nf90_put_att(ncid, lonid, 'units',     trim(lon_units))             ,cpname, __LINE__ )
504      WRITE(*,*) 'lon created'
505      !  latitude
506      CALL chkerr( nf90_def_var(ncid, 'latitude',     nf90_float, obdim, latid)        ,cpname, __LINE__ )
507      CALL chkerr( nf90_put_att(ncid, latid, 'long_name', 'Latitudes')                 ,cpname, __LINE__ )
508      CALL chkerr( nf90_put_att(ncid, latid, 'units',     trim(lat_units))             ,cpname, __LINE__ )
509      WRITE(*,*) 'lat created'
510      !  depth
511      CALL chkerr( nf90_def_var(ncid, 'depth',        nf90_float, dim2a, depid)        ,cpname, __LINE__ )
512      CALL chkerr( nf90_put_att(ncid, depid, 'long_name', 'Depths')                    ,cpname, __LINE__ )
513      CALL chkerr( nf90_put_att(ncid, depid, 'units',     trim(dep_units))             ,cpname, __LINE__ )
514      CALL chkerr( nf90_put_att(ncid, depid, '_FillValue',obfillflt)                    ,cpname, __LINE__ )
515      WRITE(*,*) 'dep created'
516      !  varname
517      CALL chkerr( nf90_def_var(ncid, 'varname',      nf90_char,  dim2c, varid)        ,cpname, __LINE__ )
518      CALL chkerr( nf90_put_att(ncid, varid, 'long_name', 'Variable name')             ,cpname, __LINE__ )
519      WRITE(*,*) 'varname created'
520      !  unitname
521      CALL chkerr( nf90_def_var(ncid, 'unitname',     nf90_char,  dim2c, unitid)        ,cpname, __LINE__ )
522      CALL chkerr( nf90_put_att(ncid,  unitid, 'long_name', 'Unit name')              ,cpname, __LINE__ )
523      WRITE(*,*) 'unitname created'
524      !  obs
525      CALL chkerr( nf90_def_var(ncid, 'observation', nf90_float,  dim3a,  obvid)        ,cpname, __LINE__ )
526      CALL chkerr( nf90_put_att(ncid, obvid, '_FillValue',obfillflt)                    ,cpname, __LINE__ )
527      CALL chkerr( nf90_put_att(ncid,  obvid, 'long_name', 'Observation value')         ,cpname, __LINE__ )
528      WRITE(*,*) 'obs created'
529      !  forecast
530      CALL chkerr( nf90_def_var(ncid, 'forecast',    nf90_float,  dim4a,  fcvid)        ,cpname, __LINE__ )
531      CALL chkerr( nf90_put_att(ncid, fcvid, '_FillValue',obfillflt)                    ,cpname, __LINE__ )
532      CALL chkerr( nf90_put_att(ncid, fcvid, 'long_name', 'Model forecast counterpart of obs. value')   ,cpname, __LINE__ )
533      CALL chkerr( nf90_put_att(ncid, fcvid, 'comment',   trim(fcst_comment))           ,cpname, __LINE__ )
534      WRITE(*,*) 'forecast created'
535      !  persistence
536      CALL chkerr( nf90_def_var(ncid, 'persistence', nf90_float,  dim4a,  prvid)        ,cpname, __LINE__ )
537      CALL chkerr( nf90_put_att(ncid, prvid, '_FillValue',obfillflt)                    ,cpname, __LINE__ )
538      CALL chkerr( nf90_put_att(ncid, prvid, 'long_name', 'Model persistence counterpart of obs. value'),cpname, __LINE__ )
539      CALL chkerr( nf90_put_att(ncid, prvid, 'comment',   trim(per_comment))            ,cpname, __LINE__ )
540      WRITE(*,*) 'persistence created'
541      !  clim
542      CALL chkerr( nf90_def_var(ncid, 'climatology', nf90_float,  dim3a, clvid)        ,cpname, __LINE__ )
543      CALL chkerr( nf90_put_att(ncid, clvid, '_FillValue',obfillflt)                    ,cpname, __LINE__ )
544      CALL chkerr( nf90_put_att(ncid, clvid, 'long_name', 'Climatological value')       ,cpname, __LINE__ )
545      CALL chkerr( nf90_put_att(ncid, clvid, 'comment',   trim(cli_comment))            ,cpname, __LINE__ )
546      WRITE(*,*) 'clim created'
547      IF (ln_best) THEN
548         !  daym2
549         CALL chkerr( nf90_def_var(ncid, 'best_estimate', nf90_float,  dim3a, dm2id)    ,cpname, __LINE__ )
550         CALL chkerr( nf90_put_att(ncid, dm2id, '_FillValue',obfillflt)                 ,cpname, __LINE__ )
551         CALL chkerr( nf90_put_att(ncid, dm2id, 'long_name', 'Best estimate')           ,cpname, __LINE__ )
552         CALL chkerr( nf90_put_att(ncid, dm2id, 'comment',   trim(dm2_comment))         ,cpname, __LINE__ )
553         WRITE(*,*) 'daym2 created'
554      ENDIF
555      IF (ln_init) THEN
556         !  daym1
557         CALL chkerr( nf90_def_var(ncid, 'nrt_analysis', nf90_float,  dim3a, dm1id)        ,cpname, __LINE__ )
558         CALL chkerr( nf90_put_att(ncid, dm1id, '_FillValue',obfillflt)                    ,cpname, __LINE__ )
559         CALL chkerr( nf90_put_att(ncid, dm1id, 'long_name', 'Near real time analysis')    ,cpname, __LINE__ )
560         CALL chkerr( nf90_put_att(ncid, dm1id, 'comment',   trim(dm1_comment))            ,cpname, __LINE__ )
561         WRITE(*,*) 'daym1 created'
562      ENDIF
563      IF (ln_mdt) THEN
564         !  mdt
565         CALL chkerr( nf90_def_var(ncid, 'mdt_reference', nf90_float,  dim3a, mdtid)      ,cpname, __LINE__ )
566         CALL chkerr( nf90_put_att(ncid, mdtid, '_FillValue',obfillflt)                    ,cpname, __LINE__ )
567         CALL chkerr( nf90_put_att(ncid, mdtid, 'long_name', 'Mean dynamic topography')    ,cpname, __LINE__ )
568         WRITE(*,*) 'mdt created'
569      ENDIF
570      IF (ln_altbias) THEN
571         !  altbias
572         CALL chkerr( nf90_def_var(ncid, 'altimeter_bias', nf90_float,  dim3a, altid)     ,cpname, __LINE__ )
573         CALL chkerr( nf90_put_att(ncid, altid, '_FillValue',obfillflt)                    ,cpname, __LINE__ )
574         CALL chkerr( nf90_put_att(ncid, altid, 'long_name', 'Altimeter bias')             ,cpname, __LINE__ )
575         WRITE(*,*) 'altbias created'
576      ENDIF
577      !  qc
578      CALL chkerr( nf90_def_var(ncid, 'qc',          nf90_short,  dim3a, qcvid)         ,cpname, __LINE__ )
579      CALL chkerr( nf90_put_att(ncid, qcvid, '_FillValue', NF90_FILL_SHORT)             ,cpname, __LINE__ )
580      CALL chkerr( nf90_put_att(ncid, qcvid, 'long_name', 'Quality flags')              ,cpname, __LINE__ )
581      CALL chkerr( nf90_put_att(ncid, qcvid, 'flag_value', qc_flag_value)               ,cpname, __LINE__ )
582      CALL chkerr( nf90_put_att(ncid, qcvid, 'flag_meaning', qc_flag_meaning)           ,cpname, __LINE__ )
583      CALL chkerr( nf90_put_att(ncid, qcvid, 'comment', qc_comment)                     ,cpname, __LINE__ )
584      WRITE(*,*) 'qc created'
585      !  juld
586      CALL chkerr( nf90_def_var(ncid, 'juld',       nf90_double,  obdim, jdvid)         ,cpname, __LINE__ )
587      CALL chkerr( nf90_put_att(ncid, jdvid, '_FillValue',99999.)                       ,cpname, __LINE__ )
588      CALL chkerr( nf90_put_att(ncid, jdvid, 'long_name', 'Observation time in Julian days'),cpname, __LINE__ )
589      CALL chkerr( nf90_put_att(ncid, jdvid, 'units',     trim(jul_units))              ,cpname, __LINE__ )
590      WRITE(*,*) 'juld created'
591      !  modeljuld
592      CALL chkerr( nf90_def_var(ncid, 'modeljuld',  nf90_double,  fcdim, mjdid)         ,cpname, __LINE__ )
593      CALL chkerr( nf90_put_att(ncid, mjdid, 'long_name', 'Model field date in Julian days'),cpname, __LINE__ )
594      CALL chkerr( nf90_put_att(ncid, mjdid, 'units',     trim(mjd_units))              ,cpname, __LINE__ )
595      WRITE(*,*) 'modeljuld created'
596      !  type
597      CALL chkerr( nf90_def_var(ncid, 'type',     nf90_char,      dim2d, typid)         ,cpname, __LINE__ )
598      CALL chkerr( nf90_put_att(ncid, typid, 'long_name', 'Observation type')           ,cpname, __LINE__ )
599      WRITE(*,*) 'type created'
600      !  id
601      CALL chkerr( nf90_def_var(ncid, 'id',       nf90_char,      dim2b, idvid)         ,cpname, __LINE__ )
602      CALL chkerr( nf90_put_att(ncid, idvid, 'long_name', 'Observation id')             ,cpname, __LINE__ )
603      WRITE(*,*) 'id created'
604      ! Close Netcdf file
605      CALL chkerr( nf90_close(ncid) ,cpname, __LINE__ ) 
606      !! Fill in the variables
607      CALL chkerr( nf90_open(trim(cdoutfile),nf90_write,ncid),                   cpname, __LINE__ )
608      WRITE(*,*) 'Create the variables ',trim(cdoutfile)
609      !  Forecast day
610      CALL chkerr( nf90_inq_varid(ncid, 'leadtime', fdvid)                ,cpname, __LINE__ )
611      CALL chkerr( nf90_put_var(ncid, fdvid, fcday) ,cpname, __LINE__ )
612      WRITE(*,*) 'forecast day put'
613      !  longitude
614      CALL chkerr( nf90_inq_varid(ncid, 'longitude', lonid)                  ,cpname, __LINE__ )
615      CALL chkerr( nf90_put_var(ncid, lonid, g_lam) ,cpname, __LINE__ )
616      WRITE(*,*) 'lon put'
617      !  latitude
618      CALL chkerr( nf90_inq_varid(ncid, 'latitude', latid)                   ,cpname, __LINE__ )
619      CALL chkerr( nf90_put_var(ncid, latid, g_phi) ,cpname, __LINE__ )
620      WRITE(*,*) 'lat put'
621      !  depth
622      CALL chkerr( nf90_inq_varid(ncid, 'depth',depid)                       ,cpname, __LINE__ )
623      CALL chkerr( nf90_put_var(ncid, depid, g_dep) ,cpname, __LINE__ )
624      WRITE(*,*) 'dep put'
625      ! varname
626      CALL chkerr( nf90_inq_varid(ncid, 'varname', varid)                    ,cpname, __LINE__ )
627      CALL chkerr( nf90_put_var(ncid, varid, gvnam,(/1,1/), (/nstr,nvars/) ) ,cpname, __LINE__ )   
628      WRITE(*,*) 'var put'
629      ! unitname
630      CALL chkerr( nf90_inq_varid(ncid, 'unitname',unitid)                   ,cpname, __LINE__ )
631      CALL chkerr( nf90_put_var(ncid, unitid, gunit,(/1,1/),(/nstr,nvars/) ) ,cpname, __LINE__ )   
632      WRITE(*,*) 'unitnam put'
633      ! obs
634      CALL chkerr( nf90_inq_varid(ncid, 'observation', obvid)                ,cpname, __LINE__ )
635      CALL chkerr( nf90_put_var(ncid, obvid,g3dob ) ,cpname, __LINE__ )   
636      WRITE(*,*) 'obs put'
637      ! clim
638      CALL chkerr( nf90_inq_varid(ncid, 'climatology', clvid)                ,cpname, __LINE__ )
639      CALL chkerr( nf90_put_var(ncid, clvid,g3dcl ) ,cpname, __LINE__ )   
640      WRITE(*,*) 'cli put'
641      IF (ln_best) THEN
642         ! daym2
643         CALL chkerr( nf90_inq_varid(ncid, 'best_estimate',dm2id)            ,cpname, __LINE__ )
644         CALL chkerr( nf90_put_var(ncid, dm2id,g3dm2 ) ,cpname, __LINE__ )   
645         WRITE(*,*) 'daym2 put'
646      ENDIF
647      IF (ln_init) THEN 
648         ! daym1
649         CALL chkerr( nf90_inq_varid(ncid, 'nrt_analysis',dm1id)              ,cpname, __LINE__ )
650         CALL chkerr( nf90_put_var(ncid, dm1id,g3dm1 ) ,cpname, __LINE__ )   
651         WRITE(*,*) 'daym1 put'
652      ENDIF
653      IF (ln_mdt) THEN 
654         ! mdt
655         CALL chkerr( nf90_inq_varid(ncid, 'mdt_reference', mdtid)              ,cpname, __LINE__ )
656         CALL chkerr( nf90_put_var(ncid, mdtid, g3mdt ) ,cpname, __LINE__ )   
657         WRITE(*,*) 'mdt put'
658      ENDIF
659      IF (ln_altbias) THEN 
660         ! altbias
661         CALL chkerr( nf90_inq_varid(ncid, 'altimeter_bias', altid)             ,cpname, __LINE__ )
662         CALL chkerr( nf90_put_var(ncid, altid, g3alt ) ,cpname, __LINE__ )   
663         WRITE(*,*) 'altbias put'
664      ENDIF
665      ! persistence
666      CALL chkerr( nf90_inq_varid(ncid, 'persistence',prvid)                                 ,cpname, __LINE__ )
667      CALL chkerr( nf90_put_var(ncid, prvid, g3dpr, (/1,1,1,1/) ,(/ ndeps,nfcst,nvars,nobs/) ) ,cpname, __LINE__ )   
668      WRITE(*,*) 'per put'
669      ! forecast
670      CALL chkerr( nf90_inq_varid(ncid, 'forecast',fcvid)                                      ,cpname, __LINE__ )
671      CALL chkerr( nf90_put_var(ncid, fcvid, g3dmc, (/1,1,1,1/), (/ ndeps,nfcst,nvars,nobs/) ) ,cpname, __LINE__ )   
672      WRITE(*,*) 'fcst put'
673      ! qc
674      CALL chkerr( nf90_inq_varid(ncid, 'qc', qcvid)                         ,cpname, __LINE__ )
675      CALL chkerr( nf90_put_var(ncid, qcvid,g3dqc ) ,cpname, __LINE__ )   
676      WRITE(*,*) 'qc put'
677      ! juld
678      CALL chkerr( nf90_inq_varid(ncid, 'juld',jdvid)                        ,cpname, __LINE__ )
679      CALL chkerr( nf90_put_var(ncid, jdvid, gjuld) ,cpname, __LINE__ )   
680      WRITE(*,*) 'juld put'
681      ! modeljuld
682      CALL chkerr( nf90_inq_varid(ncid, 'modeljuld', mjdid)                  ,cpname, __LINE__ )
683      CALL chkerr( nf90_put_var(ncid, mjdid, modjd,(/1/),(/nfcst/))          ,cpname, __LINE__ )   
684      WRITE(*,*) 'modjuld put'
685      ! type
686      CALL chkerr( nf90_inq_varid(ncid, 'type', typid)                       ,cpname, __LINE__ )
687      CALL chkerr( nf90_put_var(ncid, typid, gtype,(/1,1/) , (/n128,nobs/) ) ,cpname, __LINE__ )   
688      WRITE(*,*) 'type put'
689      ! id
690      CALL chkerr( nf90_inq_varid(ncid, 'id', idvid)                         ,cpname, __LINE__ )
691      CALL chkerr( nf90_put_var(ncid, idvid, g_id,(/1,1/)  , (/nstr,nobs/) ) ,cpname, __LINE__ )   
692      WRITE(*,*) 'id put'
693      ! Close netcdf file
694      CALL chkerr( nf90_close(ncid),                                          cpname, __LINE__ )
695   END IF ! ln_cre
696   !! Deallocate Global arrays
697   DEALLOCATE( g_lam, g_phi, g_dep, g3dob, g3dmc, g3dpr, g3dcl, g3dm2, g3dm1, g3mdt, g3alt, g3dqc, gjuld, gtype, g_id, gvnam, gunit)
698   DEALLOCATE( fcday, modjd )
699
700   !! Deallocate input argument list
701   DEALLOCATE(cdinfile)
702END PROGRAM c4comb
Note: See TracBrowser for help on using the repository browser.