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.
obs_write.F90 in branches/devukmo2010/NEMO/OPA_SRC/OBS – NEMO

source: branches/devukmo2010/NEMO/OPA_SRC/OBS/obs_write.F90 @ 2128

Last change on this file since 2128 was 2128, checked in by rfurner, 14 years ago

merged branches OBS, ASM, Rivers, BDY & mixed_dynldf ready for vn3.3 merge

File size: 32.3 KB
Line 
1MODULE obs_write
2   !!======================================================================
3   !!                       ***  MODULE obs_write   ***
4   !! Observation diagnosticss: Write observation related diagnostics
5   !!=====================================================================
6
7   !!----------------------------------------------------------------------
8   !!   cdf_wri_p3d   : Write profile observation diagnostics in NetCDF format
9   !!   obs_wri_sla   : Write SLA observation related diagnostics
10   !!   obs_wri_sst   : Write SST observation related diagnostics
11   !!   obs_wri_seaice: Write seaice observation related diagnostics
12   !!   cdf_wri_vel   : Write velocity observation diagnostics in NetCDF format
13   !!----------------------------------------------------------------------
14
15   !! * Modules used
16   USE par_kind, ONLY : &   ! Precision variables
17      & wp
18   USE in_out_manager       ! I/O manager
19   USE dom_oce              ! Ocean space and time domain variables
20   USE obs_types            ! Observation type integer to character translation
21   USE julian, ONLY : &         ! Julian date routines
22      & greg2jul
23   USE obs_utils, ONLY : &  ! Observation operator utility functions
24      & chkerr
25   USE obs_profiles_def     ! Type definitions for profiles
26   USE obs_surf_def         ! Type defintions for surface observations
27   USE obs_fbm              ! Observation feedback I/O
28   USE obs_grid             ! Grid tools
29   USE obs_conv             ! Conversion between units
30   USE obs_const
31   USE obs_sla_types
32   USE obs_rot_vel          ! Rotation of velocities
33
34   IMPLICIT NONE
35
36   !! * Routine accessibility
37   PRIVATE
38   PUBLIC obs_wri_p3d, &    ! Write profile observation related diagnostics
39      &   obs_wri_sla, &    ! Write SLA observation related diagnostics
40      &   obs_wri_sst, &    ! Write SST observation related diagnostics
41      &   obs_wri_sss, &    ! Write SSS observation related diagnostics
42      &   obs_wri_seaice, & ! Write seaice observation related diagnostics
43      &   obs_wri_vel, &    ! Write velocity observation related diagnostics
44      &   obswriinfo
45   
46   TYPE obswriinfo
47      INTEGER :: inum
48      INTEGER, POINTER, DIMENSION(:) :: ipoint
49      CHARACTER(len=ilenname), POINTER, DIMENSION(:) :: cdname
50      CHARACTER(len=ilenlong), POINTER, DIMENSION(:,:) :: cdlong
51      CHARACTER(len=ilenunit), POINTER, DIMENSION(:,:) :: cdunit
52   END TYPE obswriinfo
53
54CONTAINS
55
56   SUBROUTINE obs_wri_p3d( cprefix, profdata, padd, pext )
57      !!-----------------------------------------------------------------------
58      !!
59      !!                     *** ROUTINE obs_wri_p3d  ***
60      !!
61      !! ** Purpose : Write temperature and salinity (profile) observation
62      !!              related diagnostics
63      !!
64      !! ** Method  : NetCDF
65      !!
66      !! ** Action  :
67      !!
68      !! History :
69      !!      ! 06-04  (A. Vidard) Original
70      !!      ! 06-04  (A. Vidard) Reformatted
71      !!      ! 06-10  (A. Weaver) Cleanup
72      !!      ! 07-01  (K. Mogensen) Use profile data types
73      !!      ! 07-03  (K. Mogensen) General handling of profiles
74      !!      ! 09-01  (K. Mogensen) New feedback format
75      !!-----------------------------------------------------------------------
76
77      !! * Modules used
78
79      !! * Arguments
80      CHARACTER(LEN=*), INTENT(IN) :: cprefix        ! Prefix for output files
81      TYPE(obs_prof), INTENT(INOUT) :: profdata      ! Full set of profile data
82      TYPE(obswriinfo), OPTIONAL :: padd             ! Additional info for each variable
83      TYPE(obswriinfo), OPTIONAL :: pext             ! Extra info
84     
85      !! * Local declarations
86      TYPE(obfbdata) :: fbdata
87      CHARACTER(LEN=40) :: cfname
88      INTEGER :: ilevel
89      INTEGER :: jvar
90      INTEGER :: jo
91      INTEGER :: jk
92      INTEGER :: ik
93      INTEGER :: ja
94      INTEGER :: je
95      REAL(wp) :: zpres
96      INTEGER :: nadd
97      INTEGER :: next
98
99      IF ( PRESENT( padd ) ) THEN
100         nadd = padd%inum
101      ELSE
102         nadd = 0
103      ENDIF
104
105      IF ( PRESENT( pext ) ) THEN
106         next = pext%inum
107      ELSE
108         next = 0
109      ENDIF
110     
111      CALL init_obfbdata( fbdata )
112
113      ! Find maximum level
114      ilevel = 0
115      DO jvar = 1, 2
116         ilevel = MAX( ilevel, MAXVAL( profdata%var(jvar)%nvlidx(:) ) )
117      END DO
118      CALL alloc_obfbdata( fbdata, 2, profdata%nprof, ilevel, &
119         &                 1 + nadd, 1 + next, .TRUE. )
120
121      fbdata%cname(1)      = 'POTM'
122      fbdata%cname(2)      = 'PSAL'
123      fbdata%coblong(1)    = 'Potential temperature'
124      fbdata%coblong(2)    = 'Practical salinity'
125      fbdata%cobunit(1)    = 'Degrees centigrade'
126      fbdata%cobunit(2)    = 'PSU'
127      fbdata%cextname(1)   = 'TEMP'
128      fbdata%cextlong(1)   = 'Insitu temperature'
129      fbdata%cextunit(1)   = 'Degrees centigrade'
130      DO je = 1, next
131         fbdata%cextname(1+je) = pext%cdname(je)
132         fbdata%cextlong(1+je) = pext%cdlong(je,1)
133         fbdata%cextunit(1+je) = pext%cdunit(je,1)
134      END DO
135      fbdata%caddname(1)   = 'Hx'
136      fbdata%caddlong(1,1) = 'Model interpolated potential temperature'
137      fbdata%caddlong(1,2) = 'Model interpolated practical salinity'
138      fbdata%caddunit(1,1) = 'Degrees centigrade'
139      fbdata%caddunit(1,2) = 'PSU'
140      fbdata%cgrid(:)      = 'T'
141      DO ja = 1, nadd
142         fbdata%caddname(1+ja) = padd%cdname(ja)
143         DO jvar = 1, 2
144            fbdata%caddlong(1+ja,jvar) = padd%cdlong(ja,jvar)
145            fbdata%caddunit(1+ja,jvar) = padd%cdunit(ja,jvar)
146         END DO
147      END DO
148         
149      WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc
150
151      IF(lwp) THEN
152         WRITE(numout,*)
153         WRITE(numout,*)'obs_wri_p3d :'
154         WRITE(numout,*)'~~~~~~~~~~~~~'
155         WRITE(numout,*)'Writing profile feedback file : ',TRIM(cfname)
156      ENDIF
157
158      ! Transform obs_prof data structure into obfbdata structure
159      fbdata%cdjuldref = '19500101000000'
160      DO jo = 1, profdata%nprof
161         fbdata%plam(jo)      = profdata%rlam(jo)
162         fbdata%pphi(jo)      = profdata%rphi(jo)
163         WRITE(fbdata%cdtyp(jo),'(I4)') profdata%ntyp(jo)
164         fbdata%ivqc(jo,:)    = profdata%ivqc(jo,:)
165         fbdata%ivqcf(:,jo,:) = profdata%ivqcf(:,jo,:)
166         IF ( profdata%nqc(jo) > 10 ) THEN
167            fbdata%ioqc(jo)    = 4
168            fbdata%ioqcf(1,jo) = profdata%nqcf(1,jo)
169            fbdata%ioqcf(2,jo) = profdata%nqc(jo) - 10
170         ELSE
171            fbdata%ioqc(jo)    = profdata%nqc(jo)
172            fbdata%ioqcf(:,jo) = profdata%nqcf(:,jo)
173         ENDIF
174         fbdata%ipqc(jo)      = profdata%ipqc(jo)
175         fbdata%ipqcf(:,jo)   = profdata%ipqcf(:,jo)
176         fbdata%itqc(jo)      = profdata%itqc(jo)
177         fbdata%itqcf(:,jo)   = profdata%itqcf(:,jo)
178         fbdata%cdwmo(jo)     = profdata%cwmo(jo)
179         fbdata%kindex(jo)    = profdata%npfil(jo)
180         DO jvar = 1, profdata%nvar
181            IF (ln_grid_global) THEN
182               fbdata%iobsi(jo,jvar) = profdata%mi(jo,jvar)
183               fbdata%iobsj(jo,jvar) = profdata%mj(jo,jvar)
184            ELSE
185               fbdata%iobsi(jo,jvar) = mig(profdata%mi(jo,jvar))
186               fbdata%iobsj(jo,jvar) = mjg(profdata%mj(jo,jvar))
187            ENDIF
188         END DO
189         CALL greg2jul( 0, &
190            &           profdata%nmin(jo), &
191            &           profdata%nhou(jo), &
192            &           profdata%nday(jo), &
193            &           profdata%nmon(jo), &
194            &           profdata%nyea(jo), &
195            &           fbdata%ptim(jo),   &
196            &           krefdate = 19500101 )
197         ! Reform the profiles arrays for output
198         DO jvar = 1, 2
199            DO jk = profdata%npvsta(jo,jvar), profdata%npvend(jo,jvar)
200               ik = profdata%var(jvar)%nvlidx(jk)
201               fbdata%padd(ik,jo,1,jvar) = profdata%var(jvar)%vmod(jk)
202               fbdata%pob(ik,jo,jvar)    = profdata%var(jvar)%vobs(jk)
203               fbdata%pdep(ik,jo)        = profdata%var(jvar)%vdep(jk)
204               fbdata%idqc(ik,jo)        = profdata%var(jvar)%idqc(jk)
205               fbdata%idqcf(:,ik,jo)     = profdata%var(jvar)%idqcf(:,jk)
206               IF ( profdata%var(jvar)%nvqc(jk) > 10 ) THEN
207                  fbdata%ivlqc(ik,jo,jvar) = 4
208                  fbdata%ivlqcf(1,ik,jo,jvar) = profdata%var(jvar)%nvqcf(1,jk)
209                  fbdata%ivlqcf(2,ik,jo,jvar) = profdata%var(jvar)%nvqc(jk) - 10
210               ELSE
211                  fbdata%ivlqc(ik,jo,jvar) = profdata%var(jvar)%nvqc(jk)
212                  fbdata%ivlqcf(:,ik,jo,jvar) = profdata%var(jvar)%nvqcf(:,jk)
213               ENDIF
214               fbdata%iobsk(ik,jo,jvar)  = profdata%var(jvar)%mvk(jk)
215               DO ja = 1, nadd
216                  fbdata%padd(ik,jo,1+ja,jvar) = &
217                     & profdata%var(jvar)%vext(jk,padd%ipoint(ja))
218               END DO
219               DO je = 1, next
220                  fbdata%pext(ik,jo,1+je) = &
221                     & profdata%var(jvar)%vext(jk,pext%ipoint(je))
222               END DO
223               IF ( jvar == 1 ) THEN
224                  fbdata%pext(ik,jo,1) = profdata%var(jvar)%vext(jk,1)
225               ENDIF
226            END DO
227         END DO
228      END DO
229
230      ! Convert insitu temperature to potential temperature using the model
231      ! salinity if no potential temperature
232      DO jo = 1, fbdata%nobs
233         IF ( fbdata%pphi(jo) < 9999.0 ) THEN
234            DO jk = 1, fbdata%nlev
235               IF ( ( fbdata%pob(jk,jo,1) >= 9999.0 ) .AND. &
236                  & ( fbdata%pdep(jk,jo) < 9999.0 ) .AND. &
237                  & ( fbdata%padd(jk,jo,1,2) < 9999.0 ) .AND. &
238                  & ( fbdata%pext(jk,jo,1) < 9999.0 ) ) THEN
239                  zpres = dep_to_p( REAL(fbdata%pdep(jk,jo),wp), &
240                     &              REAL(fbdata%pphi(jo),wp) )
241                  fbdata%pob(jk,jo,1) = potemp( &
242                     &                     REAL(fbdata%padd(jk,jo,1,2), wp),  &
243                     &                     REAL(fbdata%pext(jk,jo,1), wp), &
244                     &                     zpres, 0.0_wp )
245               ENDIF
246            END DO
247         ENDIF
248      END DO
249     
250      ! Write the obfbdata structure
251      CALL write_obfbdata( cfname, fbdata )
252     
253      CALL dealloc_obfbdata( fbdata )
254     
255   END SUBROUTINE obs_wri_p3d
256
257   SUBROUTINE obs_wri_sla( cprefix, sladata, padd, pext )
258      !!-----------------------------------------------------------------------
259      !!
260      !!                     *** ROUTINE obs_wri_sla  ***
261      !!
262      !! ** Purpose : Write SLA observation diagnostics
263      !!              related
264      !!
265      !! ** Method  : NetCDF
266      !!
267      !! ** Action  :
268      !!
269      !!      ! 07-03  (K. Mogensen) Original
270      !!      ! 09-01  (K. Mogensen) New feedback format.
271      !!-----------------------------------------------------------------------
272
273      !! * Modules used
274      IMPLICIT NONE
275
276      !! * Arguments
277      CHARACTER(LEN=*), INTENT(IN) :: cprefix          ! Prefix for output files
278      TYPE(obs_surf), INTENT(INOUT) :: sladata         ! Full set of SLAa
279      TYPE(obswriinfo), OPTIONAL :: padd               ! Additional info for each variable
280      TYPE(obswriinfo), OPTIONAL :: pext               ! Extra info
281
282      !! * Local declarations
283      TYPE(obfbdata) :: fbdata
284      CHARACTER(LEN=40) :: cfname         ! netCDF filename
285      CHARACTER(LEN=12), PARAMETER :: cpname = 'obs_wri_sla'
286      INTEGER :: jo
287      INTEGER :: ja
288      INTEGER :: je
289      INTEGER :: nadd
290      INTEGER :: next
291
292      IF ( PRESENT( padd ) ) THEN
293         nadd = padd%inum
294      ELSE
295         nadd = 0
296      ENDIF
297
298      IF ( PRESENT( pext ) ) THEN
299         next = pext%inum
300      ELSE
301         next = 0
302      ENDIF
303
304      CALL init_obfbdata( fbdata )
305
306      CALL alloc_obfbdata( fbdata, 1, sladata%nsurf, 1, &
307         &                 2 + nadd, 1 + next, .TRUE. )
308
309      fbdata%cname(1)      = 'SLA'
310      fbdata%coblong(1)    = 'Sea level anomaly'
311      fbdata%cobunit(1)    = 'Metres'
312      fbdata%cextname(1)   = 'MDT'
313      fbdata%cextlong(1)   = 'Mean dynamic topography'
314      fbdata%cextunit(1)   = 'Metres'
315      DO je = 1, next
316         fbdata%cextname(1+je) = pext%cdname(je)
317         fbdata%cextlong(1+je) = pext%cdlong(je,1)
318         fbdata%cextunit(1+je) = pext%cdunit(je,1)
319      END DO
320      fbdata%caddname(1)   = 'Hx'
321      fbdata%caddlong(1,1) = 'Model interpolated SSH - MDT'
322      fbdata%caddunit(1,1) = 'Metres' 
323      fbdata%caddname(2)   = 'SSH'
324      fbdata%caddlong(2,1) = 'Model Sea surface height'
325      fbdata%caddunit(2,1) = 'Metres'
326      fbdata%cgrid(1)      = 'T'
327      DO ja = 1, nadd
328         fbdata%caddname(2+ja) = padd%cdname(ja)
329         fbdata%caddlong(2+ja,1) = padd%cdlong(ja,1)
330         fbdata%caddunit(2+ja,1) = padd%cdunit(ja,1)
331      END DO
332
333      WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc
334
335      IF(lwp) THEN
336         WRITE(numout,*)
337         WRITE(numout,*)'obs_wri_sla :'
338         WRITE(numout,*)'~~~~~~~~~~~~~'
339         WRITE(numout,*)'Writing SLA feedback file : ',TRIM(cfname)
340      ENDIF
341
342      ! Transform obs_prof data structure into obfbdata structure
343      fbdata%cdjuldref = '19500101000000'
344      DO jo = 1, sladata%nsurf
345         fbdata%plam(jo)      = sladata%rlam(jo)
346         fbdata%pphi(jo)      = sladata%rphi(jo)
347         WRITE(fbdata%cdtyp(jo),'(I4)') sladata%ntyp(jo)
348         fbdata%ivqc(jo,:)    = 0
349         fbdata%ivqcf(:,jo,:) = 0
350         IF ( sladata%nqc(jo) > 10 ) THEN
351            fbdata%ioqc(jo)    = 4
352            fbdata%ioqcf(1,jo) = 0
353            fbdata%ioqcf(2,jo) = sladata%nqc(jo) - 10
354         ELSE
355            fbdata%ioqc(jo)    = sladata%nqc(jo)
356            fbdata%ioqcf(:,jo) = 0
357         ENDIF
358         fbdata%ipqc(jo)      = 0
359         fbdata%ipqcf(:,jo)   = 0
360         fbdata%itqc(jo)      = 0
361         fbdata%itqcf(:,jo)   = 0
362         fbdata%cdwmo(jo)     = sladata%cwmo(jo)
363         fbdata%kindex(jo)    = sladata%nsfil(jo)
364         IF (ln_grid_global) THEN
365            fbdata%iobsi(jo,1) = sladata%mi(jo)
366            fbdata%iobsj(jo,1) = sladata%mj(jo)
367         ELSE
368            fbdata%iobsi(jo,1) = mig(sladata%mi(jo))
369            fbdata%iobsj(jo,1) = mjg(sladata%mj(jo))
370         ENDIF
371         CALL greg2jul( 0, &
372            &           sladata%nmin(jo), &
373            &           sladata%nhou(jo), &
374            &           sladata%nday(jo), &
375            &           sladata%nmon(jo), &
376            &           sladata%nyea(jo), &
377            &           fbdata%ptim(jo),   &
378            &           krefdate = 19500101 )
379         fbdata%padd(1,jo,1,1) = sladata%rmod(jo,1)
380         fbdata%padd(1,jo,2,1) = sladata%rext(jo,1)
381         fbdata%pob(1,jo,1)    = sladata%robs(jo,1) 
382         fbdata%pdep(1,jo)     = 0.0
383         fbdata%idqc(1,jo)     = 0
384         fbdata%idqcf(:,1,jo)  = 0
385         IF ( sladata%nqc(jo) > 10 ) THEN
386            fbdata%ivqc(jo,1)       = 4
387            fbdata%ivlqc(1,jo,1)    = 4
388            fbdata%ivlqcf(1,1,jo,1) = 0
389            fbdata%ivlqcf(2,1,jo,1) = sladata%nqc(jo) - 10
390         ELSE
391            fbdata%ivqc(jo,1)       = sladata%nqc(jo)
392            fbdata%ivlqc(1,jo,1)    = sladata%nqc(jo)
393            fbdata%ivlqcf(:,1,jo,1) = 0
394         ENDIF
395         fbdata%iobsk(1,jo,1)  = 0
396         fbdata%pext(1,jo,1) = sladata%rext(jo,2)
397         DO ja = 1, nadd
398            fbdata%padd(1,jo,2+ja,1) = &
399               & sladata%rext(jo,padd%ipoint(ja))
400         END DO
401         DO je = 1, next
402            fbdata%pext(1,jo,1+je) = &
403               & sladata%rext(jo,pext%ipoint(je))
404         END DO
405      END DO
406
407      ! Write the obfbdata structure
408      CALL write_obfbdata( cfname, fbdata )
409
410      CALL dealloc_obfbdata( fbdata )
411
412   END SUBROUTINE obs_wri_sla
413
414   SUBROUTINE obs_wri_sst( cprefix, sstdata, padd, pext )
415      !!-----------------------------------------------------------------------
416      !!
417      !!                     *** ROUTINE obs_wri_sst  ***
418      !!
419      !! ** Purpose : Write SST observation diagnostics
420      !!              related
421      !!
422      !! ** Method  : NetCDF
423      !!
424      !! ** Action  :
425      !!
426      !!      ! 07-07  (S. Ricci) Original
427      !!      ! 09-01  (K. Mogensen) New feedback format.
428      !!-----------------------------------------------------------------------
429
430      !! * Modules used
431      IMPLICIT NONE
432
433      !! * Arguments
434      CHARACTER(LEN=*), INTENT(IN) :: cprefix       ! Prefix for output files
435      TYPE(obs_surf), INTENT(INOUT) :: sstdata      ! Full set of SST
436      TYPE(obswriinfo), OPTIONAL :: padd            ! Additional info for each variable
437      TYPE(obswriinfo), OPTIONAL :: pext            ! Extra info
438
439      !! * Local declarations
440      TYPE(obfbdata) :: fbdata
441      CHARACTER(LEN=40) ::  cfname             ! netCDF filename
442      CHARACTER(LEN=12), PARAMETER :: cpname = 'obs_wri_sst'
443      INTEGER :: jo
444      INTEGER :: ja
445      INTEGER :: je
446      INTEGER :: nadd
447      INTEGER :: next
448
449      IF ( PRESENT( padd ) ) THEN
450         nadd = padd%inum
451      ELSE
452         nadd = 0
453      ENDIF
454
455      IF ( PRESENT( pext ) ) THEN
456         next = pext%inum
457      ELSE
458         next = 0
459      ENDIF
460
461      CALL init_obfbdata( fbdata )
462
463      CALL alloc_obfbdata( fbdata, 1, sstdata%nsurf, 1, &
464         &                 1 + nadd, next, .TRUE. )
465
466      fbdata%cname(1)      = 'SST'
467      fbdata%coblong(1)    = 'Sea surface temperature'
468      fbdata%cobunit(1)    = 'Degree centigrade'
469      DO je = 1, next
470         fbdata%cextname(je) = pext%cdname(je)
471         fbdata%cextlong(je) = pext%cdlong(je,1)
472         fbdata%cextunit(je) = pext%cdunit(je,1)
473      END DO
474      fbdata%caddname(1)   = 'Hx'
475      fbdata%caddlong(1,1) = 'Model interpolated SST'
476      fbdata%caddunit(1,1) = 'Degree centigrade'
477      fbdata%cgrid(1)      = 'T'
478      DO ja = 1, nadd
479         fbdata%caddname(1+ja) = padd%cdname(ja)
480         fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1)
481         fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1)
482      END DO
483
484      WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc
485
486      IF(lwp) THEN
487         WRITE(numout,*)
488         WRITE(numout,*)'obs_wri_sst :'
489         WRITE(numout,*)'~~~~~~~~~~~~~'
490         WRITE(numout,*)'Writing SST feedback file : ',TRIM(cfname)
491      ENDIF
492
493      ! Transform obs_prof data structure into obfbdata structure
494      fbdata%cdjuldref = '19500101000000'
495      DO jo = 1, sstdata%nsurf
496         fbdata%plam(jo)      = sstdata%rlam(jo)
497         fbdata%pphi(jo)      = sstdata%rphi(jo)
498         WRITE(fbdata%cdtyp(jo),'(I4)') sstdata%ntyp(jo)
499         fbdata%ivqc(jo,:)    = 0
500         fbdata%ivqcf(:,jo,:) = 0
501         IF ( sstdata%nqc(jo) > 10 ) THEN
502            fbdata%ioqc(jo)    = 4
503            fbdata%ioqcf(1,jo) = 0
504            fbdata%ioqcf(2,jo) = sstdata%nqc(jo) - 10
505         ELSE
506            fbdata%ioqc(jo)    = MAX(sstdata%nqc(jo),1)
507            fbdata%ioqcf(:,jo) = 0
508         ENDIF
509         fbdata%ipqc(jo)      = 0
510         fbdata%ipqcf(:,jo)   = 0
511         fbdata%itqc(jo)      = 0
512         fbdata%itqcf(:,jo)   = 0
513         fbdata%cdwmo(jo)     = ''
514         fbdata%kindex(jo)    = sstdata%nsfil(jo)
515         IF (ln_grid_global) THEN
516            fbdata%iobsi(jo,1) = sstdata%mi(jo)
517            fbdata%iobsj(jo,1) = sstdata%mj(jo)
518         ELSE
519            fbdata%iobsi(jo,1) = mig(sstdata%mi(jo))
520            fbdata%iobsj(jo,1) = mjg(sstdata%mj(jo))
521         ENDIF
522         CALL greg2jul( 0, &
523            &           sstdata%nmin(jo), &
524            &           sstdata%nhou(jo), &
525            &           sstdata%nday(jo), &
526            &           sstdata%nmon(jo), &
527            &           sstdata%nyea(jo), &
528            &           fbdata%ptim(jo),   &
529            &           krefdate = 19500101 )
530         fbdata%padd(1,jo,1,1) = sstdata%rmod(jo,1)
531         fbdata%pob(1,jo,1)    = sstdata%robs(jo,1)
532         fbdata%pdep(1,jo)     = 0.0
533         fbdata%idqc(1,jo)     = 0
534         fbdata%idqcf(:,1,jo)  = 0
535         IF ( sstdata%nqc(jo) > 10 ) THEN
536            fbdata%ivqc(jo,1)       = 4
537            fbdata%ivlqc(1,jo,1)    = 4
538            fbdata%ivlqcf(1,1,jo,1) = 0
539            fbdata%ivlqcf(2,1,jo,1) = sstdata%nqc(jo) - 10
540         ELSE
541            fbdata%ivqc(jo,1)       = MAX(sstdata%nqc(jo),1)
542            fbdata%ivlqc(1,jo,1)    = MAX(sstdata%nqc(jo),1)
543            fbdata%ivlqcf(:,1,jo,1) = 0
544         ENDIF
545         fbdata%iobsk(1,jo,1)  = 0
546         DO ja = 1, nadd
547            fbdata%padd(1,jo,1+ja,1) = &
548               & sstdata%rext(jo,padd%ipoint(ja))
549         END DO
550         DO je = 1, next
551            fbdata%pext(1,jo,je) = &
552               & sstdata%rext(jo,pext%ipoint(je))
553         END DO
554
555      END DO
556
557      ! Write the obfbdata structure
558
559      CALL write_obfbdata( cfname, fbdata )
560
561      CALL dealloc_obfbdata( fbdata )
562
563   END SUBROUTINE obs_wri_sst
564
565   SUBROUTINE obs_wri_sss
566   END SUBROUTINE obs_wri_sss
567
568   SUBROUTINE obs_wri_seaice( cprefix, seaicedata, padd, pext )
569      !!-----------------------------------------------------------------------
570      !!
571      !!                     *** ROUTINE obs_wri_seaice  ***
572      !!
573      !! ** Purpose : Write sea ice observation diagnostics
574      !!              related
575      !!
576      !! ** Method  : NetCDF
577      !!
578      !! ** Action  :
579      !!
580      !!      ! 07-07  (S. Ricci) Original
581      !!      ! 09-01  (K. Mogensen) New feedback format.
582      !!-----------------------------------------------------------------------
583
584      !! * Modules used
585      IMPLICIT NONE
586
587      !! * Arguments
588      CHARACTER(LEN=*), INTENT(IN) :: cprefix       ! Prefix for output files
589      TYPE(obs_surf), INTENT(INOUT) :: seaicedata   ! Full set of sea ice
590      TYPE(obswriinfo), OPTIONAL :: padd            ! Additional info for each variable
591      TYPE(obswriinfo), OPTIONAL :: pext            ! Extra info
592
593      !! * Local declarations
594      TYPE(obfbdata) :: fbdata
595      CHARACTER(LEN=40) :: cfname             ! netCDF filename
596      CHARACTER(LEN=12), PARAMETER :: cpname = 'obs_wri_seaice'
597      INTEGER :: jo
598      INTEGER :: ja
599      INTEGER :: je
600      INTEGER :: nadd
601      INTEGER :: next
602
603      IF ( PRESENT( padd ) ) THEN
604         nadd = padd%inum
605      ELSE
606         nadd = 0
607      ENDIF
608
609      IF ( PRESENT( pext ) ) THEN
610         next = pext%inum
611      ELSE
612         next = 0
613      ENDIF
614
615      CALL init_obfbdata( fbdata )
616
617      CALL alloc_obfbdata( fbdata, 1, seaicedata%nsurf, 1, 1, 0, .TRUE. )
618
619      fbdata%cname(1)      = 'SEAICE'
620      fbdata%coblong(1)    = 'Sea ice'
621      fbdata%cobunit(1)    = 'Fraction'
622      DO je = 1, next
623         fbdata%cextname(je) = pext%cdname(je)
624         fbdata%cextlong(je) = pext%cdlong(je,1)
625         fbdata%cextunit(je) = pext%cdunit(je,1)
626      END DO
627      fbdata%caddname(1)   = 'Hx'
628      fbdata%caddlong(1,1) = 'Model interpolated ICE'
629      fbdata%caddunit(1,1) = 'Fraction'
630      fbdata%cgrid(1)      = 'T'
631      DO ja = 1, nadd
632         fbdata%caddname(1+ja) = padd%cdname(ja)
633         fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1)
634         fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1)
635      END DO
636
637      WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc
638
639      IF(lwp) THEN
640         WRITE(numout,*)
641         WRITE(numout,*)'obs_wri_seaice :'
642         WRITE(numout,*)'~~~~~~~~~~~~~~~~'
643         WRITE(numout,*)'Writing SEAICE feedback file : ',TRIM(cfname)
644      ENDIF
645
646      ! Transform obs_prof data structure into obfbdata structure
647      fbdata%cdjuldref = '19500101000000'
648      DO jo = 1, seaicedata%nsurf
649         fbdata%plam(jo)      = seaicedata%rlam(jo)
650         fbdata%pphi(jo)      = seaicedata%rphi(jo)
651         WRITE(fbdata%cdtyp(jo),'(I4)') seaicedata%ntyp(jo)
652         fbdata%ivqc(jo,:)    = 0
653         fbdata%ivqcf(:,jo,:) = 0
654         IF ( seaicedata%nqc(jo) > 10 ) THEN
655            fbdata%ioqc(jo)    = 4
656            fbdata%ioqcf(1,jo) = 0
657            fbdata%ioqcf(2,jo) = seaicedata%nqc(jo) - 10
658         ELSE
659            fbdata%ioqc(jo)    = MAX(seaicedata%nqc(jo),1)
660            fbdata%ioqcf(:,jo) = 0
661         ENDIF
662         fbdata%ipqc(jo)      = 0
663         fbdata%ipqcf(:,jo)   = 0
664         fbdata%itqc(jo)      = 0
665         fbdata%itqcf(:,jo)   = 0
666         fbdata%cdwmo(jo)     = ''
667         fbdata%kindex(jo)    = seaicedata%nsfil(jo)
668         IF (ln_grid_global) THEN
669            fbdata%iobsi(jo,1) = seaicedata%mi(jo)
670            fbdata%iobsj(jo,1) = seaicedata%mj(jo)
671         ELSE
672            fbdata%iobsi(jo,1) = mig(seaicedata%mi(jo))
673            fbdata%iobsj(jo,1) = mjg(seaicedata%mj(jo))
674         ENDIF
675         CALL greg2jul( 0, &
676            &           seaicedata%nmin(jo), &
677            &           seaicedata%nhou(jo), &
678            &           seaicedata%nday(jo), &
679            &           seaicedata%nmon(jo), &
680            &           seaicedata%nyea(jo), &
681            &           fbdata%ptim(jo),   &
682            &           krefdate = 19500101 )
683         fbdata%padd(1,jo,1,1) = seaicedata%rmod(jo,1)
684         fbdata%pob(1,jo,1)    = seaicedata%robs(jo,1)
685         fbdata%pdep(1,jo)     = 0.0
686         fbdata%idqc(1,jo)     = 0
687         fbdata%idqcf(:,1,jo)  = 0
688         IF ( seaicedata%nqc(jo) > 10 ) THEN
689            fbdata%ivlqc(1,jo,1) = 4
690            fbdata%ivlqcf(1,1,jo,1) = 0
691            fbdata%ivlqcf(2,1,jo,1) = seaicedata%nqc(jo) - 10
692         ELSE
693            fbdata%ivlqc(1,jo,1) = MAX(seaicedata%nqc(jo),1)
694            fbdata%ivlqcf(:,1,jo,1) = 0
695         ENDIF
696         fbdata%iobsk(1,jo,1)  = 0
697         DO ja = 1, nadd
698            fbdata%padd(1,jo,1+ja,1) = &
699               & seaicedata%rext(jo,padd%ipoint(ja))
700         END DO
701         DO je = 1, next
702            fbdata%pext(1,jo,je) = &
703               & seaicedata%rext(jo,pext%ipoint(je))
704         END DO
705
706      END DO
707
708      ! Write the obfbdata structure
709      CALL write_obfbdata( cfname, fbdata )
710
711      CALL dealloc_obfbdata( fbdata )
712
713   END SUBROUTINE obs_wri_seaice
714
715   SUBROUTINE obs_wri_vel( cprefix, profdata, k2dint, padd, pext )
716      !!-----------------------------------------------------------------------
717      !!
718      !!                     *** ROUTINE obs_wri_p3d  ***
719      !!
720      !! ** Purpose : Write current (profile) observation
721      !!              related diagnostics
722      !!
723      !! ** Method  : NetCDF
724      !!
725      !! ** Action  :
726      !!
727      !! History :
728      !!      ! 09-01  (K. Mogensen) New feedback format routine
729      !!-----------------------------------------------------------------------
730
731      !! * Modules used
732
733      !! * Arguments
734      CHARACTER(LEN=*), INTENT(IN) :: cprefix       ! Prefix for output files
735      TYPE(obs_prof), INTENT(INOUT) :: profdata     ! Full set of profile data
736      INTEGER, INTENT(IN) :: k2dint                 ! Horizontal interpolation method
737      TYPE(obswriinfo), OPTIONAL :: padd            ! Additional info for each variable
738      TYPE(obswriinfo), OPTIONAL :: pext            ! Extra info
739
740      !! * Local declarations
741      TYPE(obfbdata) :: fbdata
742      CHARACTER(LEN=40) :: cfname
743      INTEGER :: ilevel
744      INTEGER :: jvar
745      INTEGER :: jk
746      INTEGER :: ik
747      INTEGER :: jo
748      INTEGER :: ja
749      INTEGER :: je
750      INTEGER :: nadd
751      INTEGER :: next
752      REAL(wp) :: zpres
753      REAL(wp), DIMENSION(:), ALLOCATABLE :: &
754         & zu, &
755         & zv
756
757      IF ( PRESENT( padd ) ) THEN
758         nadd = padd%inum
759      ELSE
760         nadd = 0
761      ENDIF
762
763      IF ( PRESENT( pext ) ) THEN
764         next = pext%inum
765      ELSE
766         next = 0
767      ENDIF
768
769      CALL init_obfbdata( fbdata )
770
771      ! Find maximum level
772      ilevel = 0
773      DO jvar = 1, 2
774         ilevel = MAX( ilevel, MAXVAL( profdata%var(jvar)%nvlidx(:) ) )
775      END DO
776      CALL alloc_obfbdata( fbdata, 2, profdata%nprof, ilevel, 2, 0, .TRUE. )
777
778      fbdata%cname(1)      = 'UVEL'
779      fbdata%cname(2)      = 'VVEL'
780      fbdata%coblong(1)    = 'Zonal velocity'
781      fbdata%coblong(2)    = 'Meridional velocity'
782      fbdata%cobunit(1)    = 'm/s'
783      fbdata%cobunit(2)    = 'm/s'
784      DO je = 1, next
785         fbdata%cextname(je) = pext%cdname(je)
786         fbdata%cextlong(je) = pext%cdlong(je,1)
787         fbdata%cextunit(je) = pext%cdunit(je,1)
788      END DO
789      fbdata%caddname(1)   = 'Hx'
790      fbdata%caddlong(1,1) = 'Model interpolated zonal velocity'
791      fbdata%caddlong(1,2) = 'Model interpolated meridional velocity'
792      fbdata%caddunit(1,1) = 'm/s'
793      fbdata%caddunit(1,2) = 'm/s'
794      fbdata%caddname(2)   = 'HxG'
795      fbdata%caddlong(2,1) = 'Model interpolated zonal velocity (model grid)'
796      fbdata%caddlong(2,2) = 'Model interpolated meridional velocity (model grid)'
797      fbdata%caddunit(2,1) = 'm/s'
798      fbdata%caddunit(2,2) = 'm/s' 
799      fbdata%cgrid(1)      = 'U' 
800      fbdata%cgrid(2)      = 'V'
801      DO ja = 1, nadd
802         fbdata%caddname(2+ja) = padd%cdname(ja)
803         fbdata%caddlong(2+ja,1) = padd%cdlong(ja,1)
804         fbdata%caddunit(2+ja,1) = padd%cdunit(ja,1)
805      END DO
806
807      WRITE(cfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(cprefix), nproc
808
809      IF(lwp) THEN
810         WRITE(numout,*)
811         WRITE(numout,*)'obs_wri_vel :'
812         WRITE(numout,*)'~~~~~~~~~~~~~'
813         WRITE(numout,*)'Writing velocuty feedback file : ',TRIM(cfname)
814      ENDIF
815
816      ALLOCATE( &
817         & zu(profdata%nvprot(1)), &
818         & zv(profdata%nvprot(2))  &
819         & )
820      CALL obs_rotvel( profdata, k2dint, zu, zv )
821
822      ! Transform obs_prof data structure into obfbdata structure
823      fbdata%cdjuldref = '19500101000000'
824      DO jo = 1, profdata%nprof
825         fbdata%plam(jo)      = profdata%rlam(jo)
826         fbdata%pphi(jo)      = profdata%rphi(jo)
827         WRITE(fbdata%cdtyp(jo),'(I4)') profdata%ntyp(jo)
828         fbdata%ivqc(jo,:)    = profdata%ivqc(jo,:)
829         fbdata%ivqcf(:,jo,:) = profdata%ivqcf(:,jo,:)
830         IF ( profdata%nqc(jo) > 10 ) THEN
831            fbdata%ioqc(jo)    = 4
832            fbdata%ioqcf(1,jo) = profdata%nqcf(1,jo)
833            fbdata%ioqcf(2,jo) = profdata%nqc(jo) - 10
834         ELSE
835            fbdata%ioqc(jo)    = profdata%nqc(jo)
836            fbdata%ioqcf(:,jo) = profdata%nqcf(:,jo)
837         ENDIF
838         fbdata%ipqc(jo)      = profdata%ipqc(jo)
839         fbdata%ipqcf(:,jo)   = profdata%ipqcf(:,jo)
840         fbdata%itqc(jo)      = profdata%itqc(jo)
841         fbdata%itqcf(:,jo)   = profdata%itqcf(:,jo)
842         fbdata%cdwmo(jo)     = profdata%cwmo(jo)
843         fbdata%kindex(jo)    = profdata%npfil(jo)
844         DO jvar = 1, profdata%nvar
845            IF (ln_grid_global) THEN
846               fbdata%iobsi(jo,jvar) = profdata%mi(jo,jvar)
847               fbdata%iobsj(jo,jvar) = profdata%mj(jo,jvar)
848            ELSE
849               fbdata%iobsi(jo,jvar) = mig(profdata%mi(jo,jvar))
850               fbdata%iobsj(jo,jvar) = mjg(profdata%mj(jo,jvar))
851            ENDIF
852         END DO
853         CALL greg2jul( 0, &
854            &           profdata%nmin(jo), &
855            &           profdata%nhou(jo), &
856            &           profdata%nday(jo), &
857            &           profdata%nmon(jo), &
858            &           profdata%nyea(jo), &
859            &           fbdata%ptim(jo),   &
860            &           krefdate = 19500101 )
861         ! Reform the profiles arrays for output
862         DO jvar = 1, 2
863            DO jk = profdata%npvsta(jo,jvar), profdata%npvend(jo,jvar)
864               ik = profdata%var(jvar)%nvlidx(jk)
865               IF ( jvar == 1 ) THEN
866                  fbdata%padd(ik,jo,1,jvar) = zu(jk)
867               ELSE
868                  fbdata%padd(ik,jo,1,jvar) = zv(jk)
869               ENDIF
870               fbdata%padd(ik,jo,2,jvar) = profdata%var(jvar)%vmod(jk)
871               fbdata%pob(ik,jo,jvar)    = profdata%var(jvar)%vobs(jk)
872               fbdata%pdep(ik,jo)        = profdata%var(jvar)%vdep(jk)
873               fbdata%idqc(ik,jo)        = profdata%var(jvar)%idqc(jk)
874               fbdata%idqcf(:,ik,jo)     = profdata%var(jvar)%idqcf(:,jk)
875               IF ( profdata%var(jvar)%nvqc(jk) > 10 ) THEN
876                  fbdata%ivlqc(ik,jo,jvar) = 4
877                  fbdata%ivlqcf(1,ik,jo,jvar) = profdata%var(jvar)%nvqcf(1,jk)
878                  fbdata%ivlqcf(2,ik,jo,jvar) = profdata%var(jvar)%nvqc(jk) - 10
879               ELSE
880                  fbdata%ivlqc(ik,jo,jvar) = profdata%var(jvar)%nvqc(jk)
881                  fbdata%ivlqcf(:,ik,jo,jvar) = profdata%var(jvar)%nvqcf(:,jk)
882               ENDIF
883               fbdata%iobsk(ik,jo,jvar)  = profdata%var(jvar)%mvk(jk)
884               DO ja = 1, nadd
885                  fbdata%padd(ik,jo,2+ja,jvar) = &
886                     & profdata%var(jvar)%vext(jk,padd%ipoint(ja))
887               END DO
888               DO je = 1, next
889                  fbdata%pext(ik,jo,je) = &
890                     & profdata%var(jvar)%vext(jk,pext%ipoint(je))
891               END DO
892            END DO
893         END DO
894      END DO
895
896      ! Write the obfbdata structure
897      CALL write_obfbdata( cfname, fbdata )
898     
899      CALL dealloc_obfbdata( fbdata )
900     
901      DEALLOCATE( &
902         & zu, &
903         & zv  &
904         & )
905
906   END SUBROUTINE obs_wri_vel
907
908END MODULE obs_write
Note: See TracBrowser for help on using the repository browser.