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.
fbprint.F90 in branches/UKMO/r5518_INGV1_WAVE-coupling/NEMOGCM/TOOLS/OBSTOOLS/src – NEMO

source: branches/UKMO/r5518_INGV1_WAVE-coupling/NEMOGCM/TOOLS/OBSTOOLS/src/fbprint.F90 @ 7152

Last change on this file since 7152 was 7152, checked in by jcastill, 7 years ago

Initial implementation of wave coupling branch - INGV wave branch + UKMO wave coupling branch

File size: 17.7 KB
Line 
1PROGRAM fbprint
2   !!---------------------------------------------------------------------
3   !!
4   !!                     ** PROGRAM fbprint **
5   !!
6   !!  ** Purpose : Print feedback file contents as text
7   !!
8   !!  ** Method  : Use of utilities from obs_fbm.
9   !!
10   !!  ** Action  :
11   !!
12   !!   Usage :
13   !!     fbprint.exe [options] inputfile
14   !!   Options :
15   !!     -b            shorter output
16   !!     -q            QC flags (nqc=1) select observations based on QC flags
17   !!     -Q            QC flags (nqc=2) select observations based on QC flags
18   !!     -B            QC flags (nqc=3) select observations based on QC flags
19   !!     -u            unsorted
20   !!     -s ID         select station ID 
21   !!     -t TYPE       select observation type
22   !!     -v NUM1-NUM2  select variable range to print by number (default all)
23   !!     -a NUM1-NUM2  select additional variable range to print by number (default all)
24   !!     -e NUM1-NUM2  select extra variable range to print by number (default all)
25   !!     -d            output date range
26   !!     -D            print depths
27   !!     -z            use zipped files
28   !!
29   !!   History :
30   !!        ! 2010 (K. Mogensen) Initial version
31   !!----------------------------------------------------------------------
32   USE toolspar_kind 
33   USE obs_fbm
34   USE index_sort
35   USE date_utils
36   USE proftools
37
38   IMPLICIT NONE
39   !
40   ! Command line arguments input file
41   !
42#ifndef NOIARGCPROTO
43   INTEGER,EXTERNAL :: iargc
44#endif
45   INTEGER :: nargs
46   CHARACTER(len=256) :: cdinfile, cdbrief
47   LOGICAL :: lbrief, lqcflags, lstat, ltyp, lsort, ldaterange, lzinp, ldepths
48   CHARACTER(len=ilenwmo) :: cdstat
49   CHARACTER(len=ilentyp) :: cdtyp
50   INTEGER :: nqc
51   INTEGER :: nvar1, nvar2, nadd1, nadd2, next1, next2
52   !
53   ! Input data
54   !
55   TYPE(obfbdata) :: obsdata
56   !
57   ! Loop variables
58   !
59   INTEGER :: ii, iarg, ip
60   !
61   ! Sorting
62   !
63   INTEGER :: iwmo
64   REAL(KIND=dp),ALLOCATABLE :: zsort(:,:)
65   INTEGER,ALLOCATABLE  :: iindex(:)
66   !
67   ! Get number of command line arguments
68   !
69   nargs  = IARGC()
70   lbrief = .FALSE.
71   lstat = .FALSE.
72   ltyp = .FALSE.
73   ldaterange = .FALSE.
74   ldepths = .FALSE.
75   cdstat = 'XXXXXXX'
76   cdtyp = 'XXXX'
77   lsort = .TRUE.
78   nqc = 0
79   nvar1 = -1
80   nvar2 = -1
81   nadd1 = -1
82   nadd2 = -1
83   next1 = -1
84   next2 = -1
85   lzinp = .FALSE.
86   IF ( nargs < 1 ) THEN
87      CALL usage()
88   ENDIF
89   iarg = 1
90   DO
91      IF ( iarg == nargs ) EXIT
92      CALL getarg( iarg, cdbrief )
93      IF ( cdbrief == '-b' ) THEN
94         lbrief = .TRUE.
95         iarg = iarg + 1
96      ELSEIF( cdbrief == '-q' ) THEN
97         lqcflags = .TRUE.
98         nqc=1
99         iarg = iarg + 1
100      ELSEIF( cdbrief == '-Q' ) THEN
101         lqcflags = .TRUE.
102         nqc=2
103         iarg = iarg + 1
104      ELSEIF( cdbrief == '-B' ) THEN
105         lqcflags = .TRUE.
106         nqc=3
107         iarg = iarg + 1
108      ELSEIF( cdbrief == '-u' ) THEN
109         lsort = .FALSE.
110         iarg = iarg + 1
111      ELSEIF ( cdbrief == '-s' ) THEN
112         lstat = .TRUE.
113         CALL getarg( iarg + 1, cdstat )
114         iarg = iarg + 2 
115      ELSEIF ( cdbrief == '-t' ) THEN
116         ltyp = .TRUE.
117         CALL getarg( iarg + 1, cdtyp )
118         iarg = iarg + 2 
119      ELSEIF ( cdbrief == '-v' ) THEN
120         CALL getarg( iarg + 1, cdbrief )
121         ip=INDEX(cdbrief,'-')
122         IF (ip==0) THEN
123            READ(cdbrief,'(I10)') nvar1
124            IF (nvar1==0) THEN
125               nvar2=-1
126            ELSE
127               nvar2 = nvar1
128            ENDIF
129         ELSEIF(ip==1) THEN
130            nvar1=1
131            READ(cdbrief(ip+1:),'(I10)') nvar2
132         ELSE
133            READ(cdbrief(1:ip-1),'(I10)') nvar1
134            READ(cdbrief(ip+1:),'(I10)') nvar2
135         ENDIF
136         iarg = iarg + 2
137      ELSEIF ( cdbrief == '-a' ) THEN
138         CALL getarg( iarg + 1, cdbrief )
139         ip=INDEX(cdbrief,'-')
140         IF (ip==0) THEN
141            READ(cdbrief,'(I10)') nadd1
142            IF (nadd1==0) THEN
143               nadd2=-1
144            ELSE
145               nadd2 = nadd1
146            ENDIF
147         ELSEIF(ip==1) THEN
148            nadd1=1
149            READ(cdbrief(ip+1:),'(I10)') nadd2
150         ELSE
151            READ(cdbrief(1:ip-1),'(I10)') nadd1
152            READ(cdbrief(ip+1:),'(I10)') nadd2
153         ENDIF
154         iarg = iarg + 2
155      ELSEIF ( cdbrief == '-e' ) THEN
156         CALL getarg( iarg + 1, cdbrief )
157         ip=INDEX(cdbrief,'-')
158         IF (ip==0) THEN
159            READ(cdbrief,'(I10)') next1
160            IF (next1==0) THEN
161               next2=-1
162            ELSE
163               next2 = next1
164            ENDIF
165         ELSEIF(ip==1) THEN
166            next1=1
167            READ(cdbrief(ip+1:),'(I10)') next2
168         ELSE
169            READ(cdbrief(1:ip-1),'(I10)') next1
170            READ(cdbrief(ip+1:),'(I10)') next2
171         ENDIF
172         iarg = iarg + 2
173      ELSEIF ( cdbrief == '-d' ) THEN
174         ldaterange=.TRUE.
175         iarg = iarg + 1
176      ELSEIF ( cdbrief == '-D' ) THEN
177         ldepths=.TRUE.
178         iarg = iarg + 1
179      ELSEIF ( cdbrief == '-z' ) THEN
180         lzinp=.TRUE.
181         iarg = iarg + 1
182      ELSE
183         CALL usage
184      ENDIF
185   ENDDO
186   CALL getarg( nargs, cdinfile )
187   !
188   ! Get input data
189   !
190   IF (lzinp) THEN
191#if defined NOSYSTEM
192      WRITE(*,*)'Compressed files need the system subroutine call'
193      CALL abort
194#else
195      CALL system( 'cp '//TRIM(cdinfile)//' fbprint_tmp.nc.gz' )
196      CALL system( 'gzip -df fbprint_tmp.nc.gz' )
197      CALL read_obfbdata( 'fbprint_tmp.nc', obsdata )
198      CALL system( 'rm -f fbprint_tmp.nc' )
199#endif
200   ELSE
201      CALL read_obfbdata( TRIM(cdinfile), obsdata )
202   ENDIF
203   CALL sealsfromargo( obsdata )
204   WRITE(*,'(2A,I9,A,I9,A)')TRIM(cdinfile), ' has ', obsdata%nobs ,&
205      & ' observations and a maximum of ', obsdata%nlev, ' levels'
206   IF (nvar1<0) THEN
207      nvar1 = 1
208      nvar2 = obsdata%nvar
209   ENDIF
210   IF (nadd1<0) THEN
211      nadd1 = 1
212      nadd2 = obsdata%nadd
213   ENDIF
214   IF (next1<0) THEN
215      next1 = 1
216      next2 = obsdata%next
217   ENDIF
218   !
219   ! Sort the data
220   !   
221   ALLOCATE(zsort(5,obsdata%nobs),iindex(obsdata%nobs))
222   IF (lsort) THEN
223      DO ii=1,obsdata%nobs
224         zsort(1,ii)=obsdata%ptim(ii)
225         zsort(2,ii)=obsdata%pphi(ii)
226         zsort(3,ii)=obsdata%plam(ii)
227         iwmo = TRANSFER( obsdata%cdwmo(ii)(1:4), iwmo )
228         zsort(4,ii) = iwmo
229         iwmo = TRANSFER( obsdata%cdwmo(ii)(5:8), iwmo )
230         zsort(5,ii) = iwmo
231      ENDDO
232      CALL index_sort_dp_n(zsort,5,iindex,obsdata%nobs)
233   ELSE
234      DO ii=1,obsdata%nobs
235         iindex(ii)=ii
236      ENDDO
237   ENDIF
238   IF (ldaterange) THEN
239      IF (obsdata%nobs>0) THEN
240         WRITE(*,'(A)')'First observation'
241         CALL print_time(obsdata%ptim(1))
242         WRITE(*,'(A)')'Last observation'
243         CALL print_time(obsdata%ptim(obsdata%nobs))
244      ENDIF
245   ELSE
246      !
247      ! Print the sorted list
248      !   
249      DO ii=1,obsdata%nobs
250         IF (lstat) THEN
251            IF (TRIM(ADJUSTL(cdstat)) /= &
252               &TRIM(ADJUSTL(obsdata%cdwmo(iindex(ii))))) CYCLE
253         ENDIF
254         IF (ltyp) THEN
255            IF (TRIM(ADJUSTL(cdtyp)) /= &
256               &TRIM(ADJUSTL(obsdata%cdtyp(iindex(ii))))) CYCLE
257         ENDIF
258         IF (ldepths) THEN
259            CALL print_depths(obsdata,iindex(ii))
260         ELSE
261            IF (lqcflags) THEN
262               CALL print_obs_qc(obsdata,iindex(ii),nqc,nvar1,nvar2)
263            ELSE
264               CALL print_obs(obsdata,iindex(ii),lbrief,&
265                  &           nvar1,nvar2,nadd1,nadd2,next1,next2)
266            ENDIF
267         ENDIF
268      ENDDO
269
270   ENDIF
271
272CONTAINS
273
274   SUBROUTINE usage
275      WRITE(*,'(A)')'Usage:'
276      WRITE(*,'(A)')'fbprint [options] inputfile'
277      CALL abort()
278   END SUBROUTINE usage
279   
280   SUBROUTINE print_depths(obsdata,iindex)
281      IMPLICIT NONE
282      TYPE(obfbdata) :: obsdata
283      INTEGER :: iindex
284      INTEGER :: kj
285      REAL :: mindep,maxdep
286
287      mindep=10000
288      maxdep=0
289      DO kj=1,obsdata%nlev
290         IF (obsdata%pdep(kj,iindex)<99999.0) THEN
291            IF (obsdata%pdep(kj,iindex)>maxdep) maxdep=obsdata%pdep(kj,iindex)
292            IF (obsdata%pdep(kj,iindex)<mindep) mindep=obsdata%pdep(kj,iindex)
293         ENDIF
294      ENDDO
295     
296      WRITE(*,*)'Fileindex           = ',obsdata%kindex(iindex)
297      WRITE(*,*)'Station identifier  = ',obsdata%cdwmo(iindex)
298      WRITE(*,*)'Station type        = ',obsdata%cdtyp(iindex)
299      WRITE(*,*)'Latitude            = ',obsdata%pphi(iindex)
300      WRITE(*,*)'Longtude            = ',obsdata%plam(iindex)
301      CALL print_time( obsdata%ptim(iindex) )
302      WRITE(*,*)'Position QC         = ',obsdata%ipqc(iindex)
303      WRITE(*,*)'Observation QC      = ',obsdata%ioqc(iindex)
304      WRITE(*,*)'Minimum obs. depth  = ',mindep
305      WRITE(*,*)'Maximum obs. depth  = ',maxdep
306      WRITE(*,*)
307
308   END SUBROUTINE print_depths
309
310   SUBROUTINE print_obs(obsdata,iindex,lshort,&
311      &                 kvar1,kvar2,kadd1,kadd2,kext1,kext2)
312      IMPLICIT NONE
313      TYPE(obfbdata) :: obsdata
314      INTEGER :: iindex
315      LOGICAL :: lshort
316      INTEGER :: kvar1,kvar2,kadd1,kadd2,kext1,kext2
317      INTEGER :: jv,ja,je,jk
318      INTEGER :: kj
319      LOGICAL :: lskip
320      CHARACTER(len=1024) :: cdfmt1,cdfmt2
321      CHARACTER(len=16) :: cdtmp
322
323      WRITE(*,*)'Fileindex           = ',obsdata%kindex(iindex)
324      WRITE(*,*)'Station identifier  = ',obsdata%cdwmo(iindex)
325      WRITE(*,*)'Station type        = ',obsdata%cdtyp(iindex)
326      WRITE(*,*)'Latitude            = ',obsdata%pphi(iindex)
327      WRITE(*,*)'Longtude            = ',obsdata%plam(iindex)
328      CALL print_time( obsdata%ptim(iindex) )
329      WRITE(*,*)'Position QC         = ',obsdata%ipqc(iindex)
330      WRITE(*,*)'Observation QC      = ',obsdata%ioqc(iindex)
331      IF (.NOT.lshort) THEN
332         DO jv = kvar1, kvar2
333            WRITE(*,*)'Variable name       = ',obsdata%cname(jv)
334            WRITE(*,*)'Variable QC         = ',obsdata%ivqc(iindex,jv)
335            IF (obsdata%lgrid) THEN
336               WRITE(*,*)'Grid I              = ',obsdata%iobsi(iindex,jv)
337               WRITE(*,*)'Grid J              = ',obsdata%iobsj(iindex,jv)
338            ENDIF
339         ENDDO
340         cdfmt1='(1X,A8,1X,A8'
341         cdfmt2='(1X,F8.2,1X,I8'
342         DO jv = kvar1, kvar2
343            cdfmt1 = TRIM(cdfmt1)//',1X,A15,1X,A8'
344            cdfmt2 = TRIM(cdfmt2)//',1X,E15.9,1X,I8'
345            IF (kadd2-kadd1+1>0) THEN
346               WRITE(cdtmp,'(I10)')kadd2-kadd1+1
347               cdfmt1 = TRIM(cdfmt1)//',1X,'//TRIM(ADJUSTL(cdtmp))//'(1X,A15)'
348               cdfmt2 = TRIM(cdfmt2)//',1X,'//TRIM(ADJUSTL(cdtmp))//'(1X,E15.9)'
349            ENDIF
350            IF (obsdata%lgrid) THEN
351               cdfmt1 = TRIM(cdfmt1)//',1X,A10'
352               cdfmt2 = TRIM(cdfmt2)//',1X,I10'
353            ENDIF
354         ENDDO
355         IF (kext2-kext1+1>0) THEN
356            WRITE(cdtmp,'(I10)')kext2-kext1+1
357            cdfmt1 = TRIM(cdfmt1)//',1X,'//TRIM(ADJUSTL(cdtmp))//'(1X,A15)'
358            cdfmt2 = TRIM(cdfmt2)//',1X,'//TRIM(ADJUSTL(cdtmp))//'(1X,E15.9)'
359         ENDIF
360         cdfmt1=TRIM(cdfmt1)//')'
361         cdfmt2=TRIM(cdfmt2)//')'
362         IF (obsdata%lgrid) THEN
363            WRITE(*,FMT=cdfmt1)&
364               & 'DEPTH', 'DEP_QC', &
365               & (TRIM(obsdata%cname(jv))//'_OBS', &
366               & TRIM(obsdata%cname(jv))//'_QC' , &
367               & (TRIM(obsdata%cname(jv))//'_'//TRIM(obsdata%caddname(ja)),&
368               & ja = kadd1, kadd2 ), &
369               & TRIM(obsdata%cname(jv))//'_K' , &
370               & jv = kvar1, kvar2 ), &
371               & ( TRIM(obsdata%cextname(ja)),&
372               & ja = kext1, kext2 )
373            DO kj=1,obsdata%nlev
374               IF (obsdata%pdep(kj,iindex)<99999.0) THEN
375                  WRITE (*,FMT=cdfmt2) &
376                     & obsdata%pdep(kj,iindex),   &
377                     & obsdata%idqc(kj,iindex),   &
378                     & ( obsdata%pob(kj,iindex,jv), obsdata%ivlqc(kj,iindex,jv), &
379                     & ( obsdata%padd(kj,iindex,ja,jv) , ja = kadd1, kadd2 ), &
380                     & obsdata%iobsk(kj,iindex,jv), &
381                     & jv = kvar1, kvar2 ), &
382                     & ( obsdata%pext(kj,iindex,ja), ja = kext1, kext2 )
383               ENDIF
384            ENDDO
385         ELSE
386            cdfmt1=TRIM(cdfmt1)//')'
387            cdfmt2=TRIM(cdfmt2)//')'
388            WRITE(*,FMT=cdfmt1)&
389               & 'DEPTH', 'DEP_QC', &
390               & (TRIM(obsdata%cname(jv))//'_OBS', &
391               & TRIM(obsdata%cname(jv))//'_QC' , &
392               & (TRIM(obsdata%cname(jv))//TRIM(obsdata%caddname(ja)),&
393               & ja = kadd1, kadd2 ), &
394               & jv = kvar1, kvar2 ), &
395               & ( TRIM(obsdata%cextname(ja)),&
396               & ja = kext1, kext2 )
397            DO kj=1,obsdata%nlev
398               IF (obsdata%pdep(kj,iindex)<99999.0) THEN
399                  WRITE (*,FMT=cdfmt2) &
400                     & obsdata%pdep(kj,iindex),   &
401                     & obsdata%idqc(kj,iindex),   &
402                     & ( obsdata%pob(kj,iindex,jv), obsdata%ivlqc(kj,iindex,jv), &
403                     & ( obsdata%padd(kj,iindex,ja,jv) , ja = kadd1, kadd2 ), &
404                     & jv = kvar1, kvar2 ), &
405                     & ( obsdata%pext(kj,iindex,ja), ja = kext1, kext2 )
406               ENDIF
407            ENDDO
408         ENDIF
409      ENDIF
410      WRITE(*,*)
411   END SUBROUTINE print_obs
412   
413   SUBROUTINE print_obs_qc(obsdata,iindex,kqc,kvar1,kvar2)
414      IMPLICIT NONE
415      TYPE(obfbdata) :: obsdata
416      INTEGER :: iindex
417      LOGICAL :: lqc
418      INTEGER :: kqc
419      INTEGER :: kvar1,kvar2
420      INTEGER :: jv,ja,je,jk
421      INTEGER :: kj
422      LOGICAL :: lskip
423      CHARACTER(len=1024) :: cdfmt1,cdfmt2
424      CHARACTER(len=16) :: cdtmp
425      INTEGER :: iqcf
426     
427      IF (kqc==2) THEN
428         lskip=.TRUE.
429         IF (obsdata%ipqc(iindex)>1) lskip=.FALSE.
430         IF (obsdata%ioqc(iindex)>1) lskip=.FALSE.
431         DO jv = kvar1, kvar2
432            IF (obsdata%ivqc(iindex,jv)>1) lskip=.FALSE.
433         ENDDO
434         DO kj=1,obsdata%nlev
435            IF (obsdata%pdep(kj,iindex)<99999.0) THEN
436               IF (obsdata%idqc(kj,iindex)>1) lskip=.FALSE.
437               DO jv = kvar1, kvar2
438                  IF (obsdata%ivlqc(kj,iindex,jv)>1) lskip=.FALSE.
439               ENDDO
440            ENDIF
441         ENDDO
442         IF (lskip) RETURN
443      ELSEIF (kqc==3) THEN
444         lskip=.TRUE.
445         DO kj=1,obsdata%nlev
446            IF (obsdata%pdep(kj,iindex)<99999.0) THEN
447               iqcf=0
448               DO jv = kvar1, kvar2
449                  IF (obsdata%ivlqc(kj,iindex,jv)>1) iqcf=iqcf+1
450                  IF (iqcf==obsdata%nvar) lskip=.FALSE.
451               ENDDO
452            ENDIF
453         ENDDO
454         IF (lskip) RETURN
455      ENDIF
456      WRITE(*,*)'Fileindex           = ',obsdata%kindex(iindex)
457      WRITE(*,*)'Station identifier  = ',obsdata%cdwmo(iindex)
458      WRITE(*,*)'Station type        = ',obsdata%cdtyp(iindex)
459      WRITE(*,*)'Latitude            = ',obsdata%pphi(iindex)
460      WRITE(*,*)'Longtude            = ',obsdata%plam(iindex)
461      CALL print_time( obsdata%ptim(iindex) )
462      WRITE(*,*)'Position QC         = ',obsdata%ipqc(iindex)
463      WRITE(*,*)'Position QC flags   = ',obsdata%ipqcf(:,iindex)
464      WRITE(*,*)'Observation QC      = ',obsdata%ioqc(iindex)
465      WRITE(*,*)'Observation QC flags= ',obsdata%ioqcf(:,iindex)
466      DO jv = kvar1, kvar2
467         WRITE(*,*)'Variable name       = ',obsdata%cname(jv)
468         WRITE(*,*)'Variable QC         = ',obsdata%ivqc(iindex,jv)
469         WRITE(*,*)'Variable QC flags   = ',obsdata%ivqcf(:,iindex,jv)
470      ENDDO
471      cdfmt1='(1X,A8,1X,A8'
472      cdfmt2='(1X,F8.2,1X,I8'
473      WRITE(cdtmp,'(I10)')obsdata%nqcf
474      cdfmt1 = TRIM(cdfmt1)//',1X,A18'
475      cdfmt2 = TRIM(cdfmt2)//',1X,'//TRIM(ADJUSTL(cdtmp))//'(I9)'
476      DO jv = kvar1, kvar2
477         cdfmt1 = TRIM(cdfmt1)//',1X,A15,1X,A8'
478         cdfmt2 = TRIM(cdfmt2)//',1X,E15.9,1X,I8'
479         WRITE(cdtmp,'(I10)')obsdata%nqcf
480         cdfmt1 = TRIM(cdfmt1)//',1X,A18'
481         cdfmt2 = TRIM(cdfmt2)//',1X,'//TRIM(ADJUSTL(cdtmp))//'(I9)'
482      ENDDO
483      IF (obsdata%next>0) THEN
484         WRITE(cdtmp,'(I10)')obsdata%next
485         cdfmt1 = TRIM(cdfmt1)//',1X,'//TRIM(ADJUSTL(cdtmp))//'(1X,A15)'
486         cdfmt2 = TRIM(cdfmt2)//',1X,'//TRIM(ADJUSTL(cdtmp))//'(1X,E15.9)'
487      ENDIF
488      cdfmt1=TRIM(cdfmt1)//')'
489      cdfmt2=TRIM(cdfmt2)//')'
490      WRITE(*,FMT=cdfmt1)&
491         & 'DEPTH', 'DEP_QC', 'DEP_QC_FLAGS', &
492         & (TRIM(obsdata%cname(jv))//'_OBS', &
493         & TRIM(obsdata%cname(jv))//'_QC' , &
494         & TRIM(obsdata%cname(jv))//'_QC_FLAGS',&
495         & jv = kvar1, kvar2 ), &
496         & ( TRIM(obsdata%cextname(ja)),&
497         & ja = 1, obsdata%next )
498      DO kj=1,obsdata%nlev
499         IF (kqc>=2)  THEN
500            lskip=.TRUE.
501            IF (obsdata%idqc(kj,iindex)>1) lskip=.FALSE.
502            DO jv = kvar1, kvar2
503               IF (obsdata%ivlqc(kj,iindex,jv)>1) lskip=.FALSE.
504            ENDDO
505            IF (lskip) CYCLE
506         ENDIF
507         IF (obsdata%pdep(kj,iindex)<99999.0) THEN
508            WRITE (*,FMT=cdfmt2) &
509               & obsdata%pdep(kj,iindex),   &
510               & obsdata%idqc(kj,iindex),   &
511               & ( obsdata%idqcf(ja,kj,iindex), ja = 1, obsdata%nqcf ), &
512               & ( obsdata%pob(kj,iindex,jv), obsdata%ivlqc(kj,iindex,jv), &
513               & ( obsdata%ivlqcf(ja,kj,iindex,jv) , ja=1, obsdata%nqcf ), &
514               & jv = kvar1, kvar2 ), &
515               & ( obsdata%pext(kj,iindex,ja), ja=1, obsdata%next )
516         ENDIF
517      ENDDO
518      WRITE(*,*)
519     
520   END SUBROUTINE print_obs_qc
521
522   SUBROUTINE print_time(ptim)
523      IMPLICIT NONE
524      REAL(fbdp) :: ptim
525      INTEGER:: iyr,imon,iday,ihou,imin,isec     
526      WRITE(*,*)'Julian date         = ',ptim
527      CALL jul2greg(isec,imin,ihou,iday,imon,iyr,ptim)
528      WRITE(*,'(1X,A,I4,2I2.2)') &
529         &      'Gregorian date      = ',iyr,imon,iday
530      WRITE(*,'(1X,A,I2.2,A1,I2.2,A1,I2.2)') &
531         &      'Time                = ',ihou,':',imin,':',isec
532   END  SUBROUTINE print_time
533
534END PROGRAM fbprint
535
536   
Note: See TracBrowser for help on using the repository browser.