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 NEMO/branches/UKMO/NEMO_4.0_mirror_text_diagnostics/src/OCE/OBS – NEMO

source: NEMO/branches/UKMO/NEMO_4.0_mirror_text_diagnostics/src/OCE/OBS/obs_write.F90 @ 10986

Last change on this file since 10986 was 10986, checked in by andmirek, 5 years ago

GMED 462 add flush

File size: 22.3 KB
Line 
1MODULE obs_write
2   !!======================================================================
3   !!                       ***  MODULE obs_write   ***
4   !! Observation diagnosticss: Write observation related diagnostics
5   !!=====================================================================
6
7   !!----------------------------------------------------------------------
8   !!   obs_wri_prof   : Write profile observations in feedback format
9   !!   obs_wri_surf   : Write surface observations in feedback format
10   !!   obs_wri_stats  : Print basic statistics on the data being written out
11   !!----------------------------------------------------------------------
12
13   !! * Modules used
14   USE par_kind, ONLY : &   ! Precision variables
15      & wp
16   USE in_out_manager       ! I/O manager
17   USE dom_oce              ! Ocean space and time domain variables
18   USE obs_types            ! Observation type integer to character translation
19   USE julian, ONLY : &         ! Julian date routines
20      & greg2jul
21   USE obs_utils, ONLY : &  ! Observation operator utility functions
22      & chkerr
23   USE obs_profiles_def     ! Type definitions for profiles
24   USE obs_surf_def         ! Type defintions for surface observations
25   USE obs_fbm              ! Observation feedback I/O
26   USE obs_grid             ! Grid tools
27   USE obs_conv             ! Conversion between units
28   USE obs_const
29   USE obs_mpp              ! MPP support routines for observation diagnostics
30   USE lib_mpp        ! MPP routines
31
32   IMPLICIT NONE
33
34   !! * Routine accessibility
35   PRIVATE
36   PUBLIC obs_wri_prof, &    ! Write profile observation files
37      &   obs_wri_surf, &    ! Write surface observation files
38      &   obswriinfo
39   
40   TYPE obswriinfo
41      INTEGER :: inum
42      INTEGER, POINTER, DIMENSION(:) :: ipoint
43      CHARACTER(len=ilenname), POINTER, DIMENSION(:) :: cdname
44      CHARACTER(len=ilenlong), POINTER, DIMENSION(:,:) :: cdlong
45      CHARACTER(len=ilenunit), POINTER, DIMENSION(:,:) :: cdunit
46   END TYPE obswriinfo
47
48   !!----------------------------------------------------------------------
49   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
50   !! $Id$
51   !! Software governed by the CeCILL license (see ./LICENSE)
52   !!----------------------------------------------------------------------
53
54CONTAINS
55
56   SUBROUTINE obs_wri_prof( profdata, padd, pext )
57      !!-----------------------------------------------------------------------
58      !!
59      !!                     *** ROUTINE obs_wri_prof  ***
60      !!
61      !! ** Purpose : Write profile feedback files
62      !!
63      !! ** Method  : NetCDF
64      !!
65      !! ** Action  :
66      !!
67      !! History :
68      !!      ! 06-04  (A. Vidard) Original
69      !!      ! 06-04  (A. Vidard) Reformatted
70      !!      ! 06-10  (A. Weaver) Cleanup
71      !!      ! 07-01  (K. Mogensen) Use profile data types
72      !!      ! 07-03  (K. Mogensen) General handling of profiles
73      !!      ! 09-01  (K. Mogensen) New feedback format
74      !!      ! 15-02  (M. Martin) Combined routine for writing profiles
75      !!-----------------------------------------------------------------------
76
77      !! * Arguments
78      TYPE(obs_prof), INTENT(INOUT) :: profdata      ! Full set of profile data
79      TYPE(obswriinfo), OPTIONAL :: padd             ! Additional info for each variable
80      TYPE(obswriinfo), OPTIONAL :: pext             ! Extra info
81
82      !! * Local declarations
83      TYPE(obfbdata) :: fbdata
84      CHARACTER(LEN=40) :: clfname
85      CHARACTER(LEN=10) :: clfiletype
86      INTEGER :: ilevel
87      INTEGER :: jvar
88      INTEGER :: jo
89      INTEGER :: jk
90      INTEGER :: ik
91      INTEGER :: ja
92      INTEGER :: je
93      INTEGER :: iadd
94      INTEGER :: iext
95      REAL(wp) :: zpres
96
97      IF ( PRESENT( padd ) ) THEN
98         iadd = padd%inum
99      ELSE
100         iadd = 0
101      ENDIF
102
103      IF ( PRESENT( pext ) ) THEN
104         iext = pext%inum
105      ELSE
106         iext = 0
107      ENDIF
108
109      CALL init_obfbdata( fbdata )
110
111      ! Find maximum level
112      ilevel = 0
113      DO jvar = 1, 2
114         ilevel = MAX( ilevel, MAXVAL( profdata%var(jvar)%nvlidx(:) ) )
115      END DO
116
117      SELECT CASE ( TRIM(profdata%cvars(1)) )
118      CASE('POTM')
119
120         clfiletype='profb'
121         CALL alloc_obfbdata( fbdata, 2, profdata%nprof, ilevel, &
122            &                 1 + iadd, 1 + iext, .TRUE. )
123         fbdata%cname(1)      = profdata%cvars(1)
124         fbdata%cname(2)      = profdata%cvars(2)
125         fbdata%coblong(1)    = 'Potential temperature'
126         fbdata%coblong(2)    = 'Practical salinity'
127         fbdata%cobunit(1)    = 'Degrees centigrade'
128         fbdata%cobunit(2)    = 'PSU'
129         fbdata%cextname(1)   = 'TEMP'
130         fbdata%cextlong(1)   = 'Insitu temperature'
131         fbdata%cextunit(1)   = 'Degrees centigrade'
132         fbdata%caddlong(1,1) = 'Model interpolated potential temperature'
133         fbdata%caddlong(1,2) = 'Model interpolated practical salinity'
134         fbdata%caddunit(1,1) = 'Degrees centigrade'
135         fbdata%caddunit(1,2) = 'PSU'
136         fbdata%cgrid(:)      = 'T'
137         DO je = 1, iext
138            fbdata%cextname(1+je) = pext%cdname(je)
139            fbdata%cextlong(1+je) = pext%cdlong(je,1)
140            fbdata%cextunit(1+je) = pext%cdunit(je,1)
141         END DO
142         DO ja = 1, iadd
143            fbdata%caddname(1+ja) = padd%cdname(ja)
144            DO jvar = 1, 2
145               fbdata%caddlong(1+ja,jvar) = padd%cdlong(ja,jvar)
146               fbdata%caddunit(1+ja,jvar) = padd%cdunit(ja,jvar)
147            END DO
148         END DO
149
150      CASE('UVEL')
151
152         clfiletype='velfb'
153         CALL alloc_obfbdata( fbdata, 2, profdata%nprof, ilevel, 1, 0, .TRUE. )
154         fbdata%cname(1)      = profdata%cvars(1)
155         fbdata%cname(2)      = profdata%cvars(2)
156         fbdata%coblong(1)    = 'Zonal velocity'
157         fbdata%coblong(2)    = 'Meridional velocity'
158         fbdata%cobunit(1)    = 'm/s'
159         fbdata%cobunit(2)    = 'm/s'
160         DO je = 1, iext
161            fbdata%cextname(je) = pext%cdname(je)
162            fbdata%cextlong(je) = pext%cdlong(je,1)
163            fbdata%cextunit(je) = pext%cdunit(je,1)
164         END DO
165         fbdata%caddlong(1,1) = 'Model interpolated zonal velocity'
166         fbdata%caddlong(1,2) = 'Model interpolated meridional velocity'
167         fbdata%caddunit(1,1) = 'm/s'
168         fbdata%caddunit(1,2) = 'm/s'
169         fbdata%cgrid(1)      = 'U' 
170         fbdata%cgrid(2)      = 'V'
171         DO ja = 1, iadd
172            fbdata%caddname(1+ja) = padd%cdname(ja)
173            fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1)
174            fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1)
175         END DO
176
177      END SELECT
178
179      fbdata%caddname(1)   = 'Hx'
180
181      WRITE(clfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(clfiletype), nproc
182
183      IF(lwp) THEN
184         WRITE(numout,*)
185         WRITE(numout,*)'obs_wri_prof :'
186         WRITE(numout,*)'~~~~~~~~~~~~~'
187         WRITE(numout,*)'Writing '//TRIM(clfiletype)//' feedback file : ',TRIM(clfname)
188         IF(lflush) CALL FLUSH(numout)
189      ENDIF
190
191      ! Transform obs_prof data structure into obfb data structure
192      fbdata%cdjuldref = '19500101000000'
193      DO jo = 1, profdata%nprof
194         fbdata%plam(jo)      = profdata%rlam(jo)
195         fbdata%pphi(jo)      = profdata%rphi(jo)
196         WRITE(fbdata%cdtyp(jo),'(I4)') profdata%ntyp(jo)
197         fbdata%ivqc(jo,:)    = profdata%ivqc(jo,:)
198         fbdata%ivqcf(:,jo,:) = profdata%ivqcf(:,jo,:)
199         IF ( profdata%nqc(jo) > 255 ) THEN
200            fbdata%ioqc(jo)    = IBSET(profdata%nqc(jo),2)
201            fbdata%ioqcf(1,jo) = profdata%nqcf(1,jo)
202            fbdata%ioqcf(2,jo) = profdata%nqc(jo)
203         ELSE
204            fbdata%ioqc(jo)    = profdata%nqc(jo)
205            fbdata%ioqcf(:,jo) = profdata%nqcf(:,jo)
206         ENDIF
207         fbdata%ipqc(jo)      = profdata%ipqc(jo)
208         fbdata%ipqcf(:,jo)   = profdata%ipqcf(:,jo)
209         fbdata%itqc(jo)      = profdata%itqc(jo)
210         fbdata%itqcf(:,jo)   = profdata%itqcf(:,jo)
211         fbdata%cdwmo(jo)     = profdata%cwmo(jo)
212         fbdata%kindex(jo)    = profdata%npfil(jo)
213         DO jvar = 1, profdata%nvar
214            IF (ln_grid_global) THEN
215               fbdata%iobsi(jo,jvar) = profdata%mi(jo,jvar)
216               fbdata%iobsj(jo,jvar) = profdata%mj(jo,jvar)
217            ELSE
218               fbdata%iobsi(jo,jvar) = mig(profdata%mi(jo,jvar))
219               fbdata%iobsj(jo,jvar) = mjg(profdata%mj(jo,jvar))
220            ENDIF
221         END DO
222         CALL greg2jul( 0, &
223            &           profdata%nmin(jo), &
224            &           profdata%nhou(jo), &
225            &           profdata%nday(jo), &
226            &           profdata%nmon(jo), &
227            &           profdata%nyea(jo), &
228            &           fbdata%ptim(jo),   &
229            &           krefdate = 19500101 )
230         ! Reform the profiles arrays for output
231         DO jvar = 1, 2
232            DO jk = profdata%npvsta(jo,jvar), profdata%npvend(jo,jvar)
233               ik = profdata%var(jvar)%nvlidx(jk)
234               fbdata%padd(ik,jo,1,jvar) = profdata%var(jvar)%vmod(jk)
235               fbdata%pob(ik,jo,jvar)    = profdata%var(jvar)%vobs(jk)
236               fbdata%pdep(ik,jo)        = profdata%var(jvar)%vdep(jk)
237               fbdata%idqc(ik,jo)        = profdata%var(jvar)%idqc(jk)
238               fbdata%idqcf(:,ik,jo)     = profdata%var(jvar)%idqcf(:,jk)
239               IF ( profdata%var(jvar)%nvqc(jk) > 255 ) THEN
240                  fbdata%ivlqc(ik,jo,jvar) = IBSET(profdata%var(jvar)%nvqc(jk),2)
241                  fbdata%ivlqcf(1,ik,jo,jvar) = profdata%var(jvar)%nvqcf(1,jk)
242!$AGRIF_DO_NOT_TREAT
243                  fbdata%ivlqcf(2,ik,jo,jvar) = IAND(profdata%var(jvar)%nvqc(jk),b'0000000011111111')
244!$AGRIF_END_DO_NOT_TREAT
245               ELSE
246                  fbdata%ivlqc(ik,jo,jvar) = profdata%var(jvar)%nvqc(jk)
247                  fbdata%ivlqcf(:,ik,jo,jvar) = profdata%var(jvar)%nvqcf(:,jk)
248               ENDIF
249               fbdata%iobsk(ik,jo,jvar)  = profdata%var(jvar)%mvk(jk)
250               DO ja = 1, iadd
251                  fbdata%padd(ik,jo,1+ja,jvar) = &
252                     & profdata%var(jvar)%vext(jk,padd%ipoint(ja))
253               END DO
254               DO je = 1, iext
255                  fbdata%pext(ik,jo,1+je) = &
256                     & profdata%var(jvar)%vext(jk,pext%ipoint(je))
257               END DO
258               IF ( ( jvar == 1 ) .AND. &
259                  & ( TRIM(profdata%cvars(1)) == 'POTM' ) ) THEN
260                  fbdata%pext(ik,jo,1) = profdata%var(jvar)%vext(jk,1)
261               ENDIF
262            END DO
263         END DO
264      END DO
265
266      IF ( TRIM(profdata%cvars(1)) == 'POTM' ) THEN
267         ! Convert insitu temperature to potential temperature using the model
268         ! salinity if no potential temperature
269         DO jo = 1, fbdata%nobs
270            IF ( fbdata%pphi(jo) < 9999.0 ) THEN
271               DO jk = 1, fbdata%nlev
272                  IF ( ( fbdata%pob(jk,jo,1) >= 9999.0 ) .AND. &
273                     & ( fbdata%pdep(jk,jo) < 9999.0 ) .AND. &
274                     & ( fbdata%padd(jk,jo,1,2) < 9999.0 ) .AND. &
275                     & ( fbdata%pext(jk,jo,1) < 9999.0 ) ) THEN
276                     zpres = dep_to_p( REAL(fbdata%pdep(jk,jo),wp), &
277                        &              REAL(fbdata%pphi(jo),wp) )
278                     fbdata%pob(jk,jo,1) = potemp( &
279                        &                     REAL(fbdata%padd(jk,jo,1,2), wp),  &
280                        &                     REAL(fbdata%pext(jk,jo,1), wp), &
281                        &                     zpres, 0.0_wp )
282                  ENDIF
283               END DO
284            ENDIF
285         END DO
286      ENDIF
287
288      ! Write the obfbdata structure
289      CALL write_obfbdata( clfname, fbdata )
290
291      ! Output some basic statistics
292      CALL obs_wri_stats( fbdata )
293
294      CALL dealloc_obfbdata( fbdata )
295
296   END SUBROUTINE obs_wri_prof
297
298   SUBROUTINE obs_wri_surf( surfdata, padd, pext )
299      !!-----------------------------------------------------------------------
300      !!
301      !!                     *** ROUTINE obs_wri_surf  ***
302      !!
303      !! ** Purpose : Write surface observation files
304      !!
305      !! ** Method  : NetCDF
306      !!
307      !! ** Action  :
308      !!
309      !!      ! 07-03  (K. Mogensen) Original
310      !!      ! 09-01  (K. Mogensen) New feedback format.
311      !!      ! 15-02  (M. Martin) Combined surface writing routine.
312      !!-----------------------------------------------------------------------
313
314      !! * Modules used
315      IMPLICIT NONE
316
317      !! * Arguments
318      TYPE(obs_surf), INTENT(INOUT) :: surfdata         ! Full set of surface data
319      TYPE(obswriinfo), OPTIONAL :: padd               ! Additional info for each variable
320      TYPE(obswriinfo), OPTIONAL :: pext               ! Extra info
321
322      !! * Local declarations
323      TYPE(obfbdata) :: fbdata
324      CHARACTER(LEN=40) :: clfname         ! netCDF filename
325      CHARACTER(LEN=10) :: clfiletype
326      CHARACTER(LEN=12), PARAMETER :: cpname = 'obs_wri_surf'
327      INTEGER :: jo
328      INTEGER :: ja
329      INTEGER :: je
330      INTEGER :: iadd
331      INTEGER :: iext
332
333      IF ( PRESENT( padd ) ) THEN
334         iadd = padd%inum
335      ELSE
336         iadd = 0
337      ENDIF
338
339      IF ( PRESENT( pext ) ) THEN
340         iext = pext%inum
341      ELSE
342         iext = 0
343      ENDIF
344
345      CALL init_obfbdata( fbdata )
346
347      SELECT CASE ( TRIM(surfdata%cvars(1)) )
348      CASE('SLA')
349
350         CALL alloc_obfbdata( fbdata, 1, surfdata%nsurf, 1, &
351            &                 2 + iadd, 1 + iext, .TRUE. )
352
353         clfiletype = 'slafb'
354         fbdata%cname(1)      = surfdata%cvars(1)
355         fbdata%coblong(1)    = 'Sea level anomaly'
356         fbdata%cobunit(1)    = 'Metres'
357         fbdata%cextname(1)   = 'MDT'
358         fbdata%cextlong(1)   = 'Mean dynamic topography'
359         fbdata%cextunit(1)   = 'Metres'
360         DO je = 1, iext
361            fbdata%cextname(je) = pext%cdname(je)
362            fbdata%cextlong(je) = pext%cdlong(je,1)
363            fbdata%cextunit(je) = pext%cdunit(je,1)
364         END DO
365         fbdata%caddlong(1,1) = 'Model interpolated SSH - MDT'
366         fbdata%caddunit(1,1) = 'Metres' 
367         fbdata%caddname(2)   = 'SSH'
368         fbdata%caddlong(2,1) = 'Model Sea surface height'
369         fbdata%caddunit(2,1) = 'Metres'
370         fbdata%cgrid(1)      = 'T'
371         DO ja = 1, iadd
372            fbdata%caddname(2+ja) = padd%cdname(ja)
373            fbdata%caddlong(2+ja,1) = padd%cdlong(ja,1)
374            fbdata%caddunit(2+ja,1) = padd%cdunit(ja,1)
375         END DO
376
377      CASE('SST')
378
379         CALL alloc_obfbdata( fbdata, 1, surfdata%nsurf, 1, &
380            &                 1 + iadd, iext, .TRUE. )
381
382         clfiletype = 'sstfb'
383         fbdata%cname(1)      = surfdata%cvars(1)
384         fbdata%coblong(1)    = 'Sea surface temperature'
385         fbdata%cobunit(1)    = 'Degree centigrade'
386         DO je = 1, iext
387            fbdata%cextname(je) = pext%cdname(je)
388            fbdata%cextlong(je) = pext%cdlong(je,1)
389            fbdata%cextunit(je) = pext%cdunit(je,1)
390         END DO
391         fbdata%caddlong(1,1) = 'Model interpolated SST'
392         fbdata%caddunit(1,1) = 'Degree centigrade'
393         fbdata%cgrid(1)      = 'T'
394         DO ja = 1, iadd
395            fbdata%caddname(1+ja) = padd%cdname(ja)
396            fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1)
397            fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1)
398         END DO
399
400      CASE('ICECONC')
401
402         CALL alloc_obfbdata( fbdata, 1, surfdata%nsurf, 1, &
403            &                 1 + iadd, iext, .TRUE. )
404
405         clfiletype = 'sicfb'
406         fbdata%cname(1)      = surfdata%cvars(1)
407         fbdata%coblong(1)    = 'Sea ice'
408         fbdata%cobunit(1)    = 'Fraction'
409         DO je = 1, iext
410            fbdata%cextname(je) = pext%cdname(je)
411            fbdata%cextlong(je) = pext%cdlong(je,1)
412            fbdata%cextunit(je) = pext%cdunit(je,1)
413         END DO
414         fbdata%caddlong(1,1) = 'Model interpolated ICE'
415         fbdata%caddunit(1,1) = 'Fraction'
416         fbdata%cgrid(1)      = 'T'
417         DO ja = 1, iadd
418            fbdata%caddname(1+ja) = padd%cdname(ja)
419            fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1)
420            fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1)
421         END DO
422
423      CASE('SSS')
424
425         CALL alloc_obfbdata( fbdata, 1, surfdata%nsurf, 1, &
426            &                 1 + iadd, iext, .TRUE. )
427
428         clfiletype = 'sssfb'
429         fbdata%cname(1)      = surfdata%cvars(1)
430         fbdata%coblong(1)    = 'Sea surface salinity'
431         fbdata%cobunit(1)    = 'psu'
432         DO je = 1, iext
433            fbdata%cextname(je) = pext%cdname(je)
434            fbdata%cextlong(je) = pext%cdlong(je,1)
435            fbdata%cextunit(je) = pext%cdunit(je,1)
436         END DO
437         fbdata%caddlong(1,1) = 'Model interpolated SSS'
438         fbdata%caddunit(1,1) = 'psu'
439         fbdata%cgrid(1)      = 'T'
440         DO ja = 1, iadd
441            fbdata%caddname(1+ja) = padd%cdname(ja)
442            fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1)
443            fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1)
444         END DO
445
446      CASE DEFAULT
447
448         CALL ctl_stop( 'Unknown observation type '//TRIM(surfdata%cvars(1))//' in obs_wri_surf' )
449
450      END SELECT
451
452      fbdata%caddname(1)   = 'Hx'
453
454      WRITE(clfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(clfiletype), nproc
455
456      IF(lwp) THEN
457         WRITE(numout,*)
458         WRITE(numout,*)'obs_wri_surf :'
459         WRITE(numout,*)'~~~~~~~~~~~~~'
460         WRITE(numout,*)'Writing '//TRIM(surfdata%cvars(1))//' feedback file : ',TRIM(clfname)
461         IF(lflush) CALL FLUSH(numout)
462      ENDIF
463
464      ! Transform surf data structure into obfbdata structure
465      fbdata%cdjuldref = '19500101000000'
466      DO jo = 1, surfdata%nsurf
467         fbdata%plam(jo)      = surfdata%rlam(jo)
468         fbdata%pphi(jo)      = surfdata%rphi(jo)
469         WRITE(fbdata%cdtyp(jo),'(I4)') surfdata%ntyp(jo)
470         fbdata%ivqc(jo,:)    = 0
471         fbdata%ivqcf(:,jo,:) = 0
472         IF ( surfdata%nqc(jo) > 255 ) THEN
473            fbdata%ioqc(jo)    = 4
474            fbdata%ioqcf(1,jo) = 0
475!$AGRIF_DO_NOT_TREAT
476            fbdata%ioqcf(2,jo) = IAND(surfdata%nqc(jo),b'0000000011111111')
477!$AGRIF_END_DO_NOT_TREAT
478         ELSE
479            fbdata%ioqc(jo)    = surfdata%nqc(jo)
480            fbdata%ioqcf(:,jo) = 0
481         ENDIF
482         fbdata%ipqc(jo)      = 0
483         fbdata%ipqcf(:,jo)   = 0
484         fbdata%itqc(jo)      = 0
485         fbdata%itqcf(:,jo)   = 0
486         fbdata%cdwmo(jo)     = surfdata%cwmo(jo)
487         fbdata%kindex(jo)    = surfdata%nsfil(jo)
488         IF (ln_grid_global) THEN
489            fbdata%iobsi(jo,1) = surfdata%mi(jo)
490            fbdata%iobsj(jo,1) = surfdata%mj(jo)
491         ELSE
492            fbdata%iobsi(jo,1) = mig(surfdata%mi(jo))
493            fbdata%iobsj(jo,1) = mjg(surfdata%mj(jo))
494         ENDIF
495         CALL greg2jul( 0, &
496            &           surfdata%nmin(jo), &
497            &           surfdata%nhou(jo), &
498            &           surfdata%nday(jo), &
499            &           surfdata%nmon(jo), &
500            &           surfdata%nyea(jo), &
501            &           fbdata%ptim(jo),   &
502            &           krefdate = 19500101 )
503         fbdata%padd(1,jo,1,1) = surfdata%rmod(jo,1)
504         IF ( TRIM(surfdata%cvars(1)) == 'SLA' ) fbdata%padd(1,jo,2,1) = surfdata%rext(jo,1)
505         fbdata%pob(1,jo,1)    = surfdata%robs(jo,1) 
506         fbdata%pdep(1,jo)     = 0.0
507         fbdata%idqc(1,jo)     = 0
508         fbdata%idqcf(:,1,jo)  = 0
509         IF ( surfdata%nqc(jo) > 255 ) THEN
510            fbdata%ivqc(jo,1)       = 4
511            fbdata%ivlqc(1,jo,1)    = 4
512            fbdata%ivlqcf(1,1,jo,1) = 0
513!$AGRIF_DO_NOT_TREAT
514            fbdata%ivlqcf(2,1,jo,1) = IAND(surfdata%nqc(jo),b'0000000011111111')
515!$AGRIF_END_DO_NOT_TREAT
516         ELSE
517            fbdata%ivqc(jo,1)       = surfdata%nqc(jo)
518            fbdata%ivlqc(1,jo,1)    = surfdata%nqc(jo)
519            fbdata%ivlqcf(:,1,jo,1) = 0
520         ENDIF
521         fbdata%iobsk(1,jo,1)  = 0
522         IF ( TRIM(surfdata%cvars(1)) == 'SLA' ) fbdata%pext(1,jo,1) = surfdata%rext(jo,2)
523         DO ja = 1, iadd
524            fbdata%padd(1,jo,2+ja,1) = &
525               & surfdata%rext(jo,padd%ipoint(ja))
526         END DO
527         DO je = 1, iext
528            fbdata%pext(1,jo,1+je) = &
529               & surfdata%rext(jo,pext%ipoint(je))
530         END DO
531      END DO
532
533      ! Write the obfbdata structure
534      CALL write_obfbdata( clfname, fbdata )
535
536      ! Output some basic statistics
537      CALL obs_wri_stats( fbdata )
538
539      CALL dealloc_obfbdata( fbdata )
540
541   END SUBROUTINE obs_wri_surf
542
543   SUBROUTINE obs_wri_stats( fbdata )
544      !!-----------------------------------------------------------------------
545      !!
546      !!                     *** ROUTINE obs_wri_stats  ***
547      !!
548      !! ** Purpose : Output some basic statistics of the data being written out
549      !!
550      !! ** Method  :
551      !!
552      !! ** Action  :
553      !!
554      !!      ! 2014-08  (D. Lea) Initial version
555      !!-----------------------------------------------------------------------
556
557      !! * Arguments
558      TYPE(obfbdata) :: fbdata
559
560      !! * Local declarations
561      INTEGER :: jvar
562      INTEGER :: jo
563      INTEGER :: jk
564      INTEGER :: inumgoodobs
565      INTEGER :: inumgoodobsmpp
566      REAL(wp) :: zsumx
567      REAL(wp) :: zsumx2
568      REAL(wp) :: zomb
569     
570
571      IF (lwp) THEN
572         WRITE(numout,*) ''
573         WRITE(numout,*) 'obs_wri_stats :'
574         WRITE(numout,*) '~~~~~~~~~~~~~~~'
575         IF(lflush) CALL FLUSH(numout)
576      ENDIF
577
578      DO jvar = 1, fbdata%nvar
579         zsumx=0.0_wp
580         zsumx2=0.0_wp
581         inumgoodobs=0
582         DO jo = 1, fbdata%nobs
583            DO jk = 1, fbdata%nlev
584               IF ( ( fbdata%pob(jk,jo,jvar) < 9999.0 ) .AND. &
585                  & ( fbdata%pdep(jk,jo) < 9999.0 ) .AND. &
586                  & ( fbdata%padd(jk,jo,1,jvar) < 9999.0 ) ) THEN
587
588                  zomb=fbdata%pob(jk, jo, jvar)-fbdata%padd(jk, jo, 1, jvar)
589                  zsumx=zsumx+zomb
590                  zsumx2=zsumx2+zomb**2
591                  inumgoodobs=inumgoodobs+1
592               ENDIF
593            ENDDO
594         ENDDO
595
596         CALL obs_mpp_sum_integer( inumgoodobs, inumgoodobsmpp )
597         CALL mpp_sum('obs_write', zsumx)
598         CALL mpp_sum('obs_write', zsumx2)
599
600         IF (lwp) THEN
601            WRITE(numout,*) 'Type: ',fbdata%cname(jvar),'  Total number of good observations: ',inumgoodobsmpp 
602            WRITE(numout,*) 'Overall mean obs minus model of the good observations: ',zsumx/inumgoodobsmpp
603            WRITE(numout,*) 'Overall RMS obs minus model of the good observations: ',sqrt( zsumx2/inumgoodobsmpp )
604            WRITE(numout,*) ''
605            IF(lflush) CALL FLUSH(numout)
606         ENDIF
607
608      ENDDO
609
610   END SUBROUTINE obs_wri_stats
611
612END MODULE obs_write
Note: See TracBrowser for help on using the repository browser.