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/UKMO/dev_r8183_ICEMODEL_svn_removed/NEMOGCM/NEMO/OPA_SRC/OBS – NEMO

source: branches/UKMO/dev_r8183_ICEMODEL_svn_removed/NEMOGCM/NEMO/OPA_SRC/OBS/obs_write.F90 @ 8733

Last change on this file since 8733 was 8733, checked in by dancopsey, 6 years ago

Remove svn keywords.

File size: 20.9 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/OPA 3.3 , NEMO Consortium (2010)
50   !! $Id$
51   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
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=6) :: 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      ENDIF
189
190      ! Transform obs_prof data structure into obfb data structure
191      fbdata%cdjuldref = '19500101000000'
192      DO jo = 1, profdata%nprof
193         fbdata%plam(jo)      = profdata%rlam(jo)
194         fbdata%pphi(jo)      = profdata%rphi(jo)
195         WRITE(fbdata%cdtyp(jo),'(I4)') profdata%ntyp(jo)
196         fbdata%ivqc(jo,:)    = profdata%ivqc(jo,:)
197         fbdata%ivqcf(:,jo,:) = profdata%ivqcf(:,jo,:)
198         IF ( profdata%nqc(jo) > 10 ) THEN
199            fbdata%ioqc(jo)    = 4
200            fbdata%ioqcf(1,jo) = profdata%nqcf(1,jo)
201            fbdata%ioqcf(2,jo) = profdata%nqc(jo) - 10
202         ELSE
203            fbdata%ioqc(jo)    = profdata%nqc(jo)
204            fbdata%ioqcf(:,jo) = profdata%nqcf(:,jo)
205         ENDIF
206         fbdata%ipqc(jo)      = profdata%ipqc(jo)
207         fbdata%ipqcf(:,jo)   = profdata%ipqcf(:,jo)
208         fbdata%itqc(jo)      = profdata%itqc(jo)
209         fbdata%itqcf(:,jo)   = profdata%itqcf(:,jo)
210         fbdata%cdwmo(jo)     = profdata%cwmo(jo)
211         fbdata%kindex(jo)    = profdata%npfil(jo)
212         DO jvar = 1, profdata%nvar
213            IF (ln_grid_global) THEN
214               fbdata%iobsi(jo,jvar) = profdata%mi(jo,jvar)
215               fbdata%iobsj(jo,jvar) = profdata%mj(jo,jvar)
216            ELSE
217               fbdata%iobsi(jo,jvar) = mig(profdata%mi(jo,jvar))
218               fbdata%iobsj(jo,jvar) = mjg(profdata%mj(jo,jvar))
219            ENDIF
220         END DO
221         CALL greg2jul( 0, &
222            &           profdata%nmin(jo), &
223            &           profdata%nhou(jo), &
224            &           profdata%nday(jo), &
225            &           profdata%nmon(jo), &
226            &           profdata%nyea(jo), &
227            &           fbdata%ptim(jo),   &
228            &           krefdate = 19500101 )
229         ! Reform the profiles arrays for output
230         DO jvar = 1, 2
231            DO jk = profdata%npvsta(jo,jvar), profdata%npvend(jo,jvar)
232               ik = profdata%var(jvar)%nvlidx(jk)
233               fbdata%padd(ik,jo,1,jvar) = profdata%var(jvar)%vmod(jk)
234               fbdata%pob(ik,jo,jvar)    = profdata%var(jvar)%vobs(jk)
235               fbdata%pdep(ik,jo)        = profdata%var(jvar)%vdep(jk)
236               fbdata%idqc(ik,jo)        = profdata%var(jvar)%idqc(jk)
237               fbdata%idqcf(:,ik,jo)     = profdata%var(jvar)%idqcf(:,jk)
238               IF ( profdata%var(jvar)%nvqc(jk) > 10 ) THEN
239                  fbdata%ivlqc(ik,jo,jvar) = 4
240                  fbdata%ivlqcf(1,ik,jo,jvar) = profdata%var(jvar)%nvqcf(1,jk)
241                  fbdata%ivlqcf(2,ik,jo,jvar) = profdata%var(jvar)%nvqc(jk) - 10
242               ELSE
243                  fbdata%ivlqc(ik,jo,jvar) = profdata%var(jvar)%nvqc(jk)
244                  fbdata%ivlqcf(:,ik,jo,jvar) = profdata%var(jvar)%nvqcf(:,jk)
245               ENDIF
246               fbdata%iobsk(ik,jo,jvar)  = profdata%var(jvar)%mvk(jk)
247               DO ja = 1, iadd
248                  fbdata%padd(ik,jo,1+ja,jvar) = &
249                     & profdata%var(jvar)%vext(jk,padd%ipoint(ja))
250               END DO
251               DO je = 1, iext
252                  fbdata%pext(ik,jo,1+je) = &
253                     & profdata%var(jvar)%vext(jk,pext%ipoint(je))
254               END DO
255               IF ( ( jvar == 1 ) .AND. &
256                  & ( TRIM(profdata%cvars(1)) == 'POTM' ) ) THEN
257                  fbdata%pext(ik,jo,1) = profdata%var(jvar)%vext(jk,1)
258               ENDIF
259            END DO
260         END DO
261      END DO
262
263      IF ( TRIM(profdata%cvars(1)) == 'POTM' ) THEN
264         ! Convert insitu temperature to potential temperature using the model
265         ! salinity if no potential temperature
266         DO jo = 1, fbdata%nobs
267            IF ( fbdata%pphi(jo) < 9999.0 ) THEN
268               DO jk = 1, fbdata%nlev
269                  IF ( ( fbdata%pob(jk,jo,1) >= 9999.0 ) .AND. &
270                     & ( fbdata%pdep(jk,jo) < 9999.0 ) .AND. &
271                     & ( fbdata%padd(jk,jo,1,2) < 9999.0 ) .AND. &
272                     & ( fbdata%pext(jk,jo,1) < 9999.0 ) ) THEN
273                     zpres = dep_to_p( REAL(fbdata%pdep(jk,jo),wp), &
274                        &              REAL(fbdata%pphi(jo),wp) )
275                     fbdata%pob(jk,jo,1) = potemp( &
276                        &                     REAL(fbdata%padd(jk,jo,1,2), wp),  &
277                        &                     REAL(fbdata%pext(jk,jo,1), wp), &
278                        &                     zpres, 0.0_wp )
279                  ENDIF
280               END DO
281            ENDIF
282         END DO
283      ENDIF
284
285      ! Write the obfbdata structure
286      CALL write_obfbdata( clfname, fbdata )
287
288      ! Output some basic statistics
289      CALL obs_wri_stats( fbdata )
290
291      CALL dealloc_obfbdata( fbdata )
292
293   END SUBROUTINE obs_wri_prof
294
295   SUBROUTINE obs_wri_surf( surfdata, padd, pext )
296      !!-----------------------------------------------------------------------
297      !!
298      !!                     *** ROUTINE obs_wri_surf  ***
299      !!
300      !! ** Purpose : Write surface observation files
301      !!
302      !! ** Method  : NetCDF
303      !!
304      !! ** Action  :
305      !!
306      !!      ! 07-03  (K. Mogensen) Original
307      !!      ! 09-01  (K. Mogensen) New feedback format.
308      !!      ! 15-02  (M. Martin) Combined surface writing routine.
309      !!-----------------------------------------------------------------------
310
311      !! * Modules used
312      IMPLICIT NONE
313
314      !! * Arguments
315      TYPE(obs_surf), INTENT(INOUT) :: surfdata         ! Full set of surface data
316      TYPE(obswriinfo), OPTIONAL :: padd               ! Additional info for each variable
317      TYPE(obswriinfo), OPTIONAL :: pext               ! Extra info
318
319      !! * Local declarations
320      TYPE(obfbdata) :: fbdata
321      CHARACTER(LEN=40) :: clfname         ! netCDF filename
322      CHARACTER(LEN=6)  :: clfiletype
323      CHARACTER(LEN=12), PARAMETER :: cpname = 'obs_wri_surf'
324      INTEGER :: jo
325      INTEGER :: ja
326      INTEGER :: je
327      INTEGER :: iadd
328      INTEGER :: iext
329
330      IF ( PRESENT( padd ) ) THEN
331         iadd = padd%inum
332      ELSE
333         iadd = 0
334      ENDIF
335
336      IF ( PRESENT( pext ) ) THEN
337         iext = pext%inum
338      ELSE
339         iext = 0
340      ENDIF
341
342      CALL init_obfbdata( fbdata )
343
344      SELECT CASE ( TRIM(surfdata%cvars(1)) )
345      CASE('SLA')
346
347         CALL alloc_obfbdata( fbdata, 1, surfdata%nsurf, 1, &
348            &                 2 + iadd, 1 + iext, .TRUE. )
349
350         clfiletype = 'slafb'
351         fbdata%cname(1)      = surfdata%cvars(1)
352         fbdata%coblong(1)    = 'Sea level anomaly'
353         fbdata%cobunit(1)    = 'Metres'
354         fbdata%cextname(1)   = 'MDT'
355         fbdata%cextlong(1)   = 'Mean dynamic topography'
356         fbdata%cextunit(1)   = 'Metres'
357         DO je = 1, iext
358            fbdata%cextname(je) = pext%cdname(je)
359            fbdata%cextlong(je) = pext%cdlong(je,1)
360            fbdata%cextunit(je) = pext%cdunit(je,1)
361         END DO
362         fbdata%caddlong(1,1) = 'Model interpolated SSH - MDT'
363         fbdata%caddunit(1,1) = 'Metres' 
364         fbdata%caddname(2)   = 'SSH'
365         fbdata%caddlong(2,1) = 'Model Sea surface height'
366         fbdata%caddunit(2,1) = 'Metres'
367         fbdata%cgrid(1)      = 'T'
368         DO ja = 1, iadd
369            fbdata%caddname(2+ja) = padd%cdname(ja)
370            fbdata%caddlong(2+ja,1) = padd%cdlong(ja,1)
371            fbdata%caddunit(2+ja,1) = padd%cdunit(ja,1)
372         END DO
373
374      CASE('SST')
375
376         CALL alloc_obfbdata( fbdata, 1, surfdata%nsurf, 1, &
377            &                 1 + iadd, iext, .TRUE. )
378
379         clfiletype = 'sstfb'
380         fbdata%cname(1)      = surfdata%cvars(1)
381         fbdata%coblong(1)    = 'Sea surface temperature'
382         fbdata%cobunit(1)    = 'Degree centigrade'
383         DO je = 1, iext
384            fbdata%cextname(je) = pext%cdname(je)
385            fbdata%cextlong(je) = pext%cdlong(je,1)
386            fbdata%cextunit(je) = pext%cdunit(je,1)
387         END DO
388         fbdata%caddlong(1,1) = 'Model interpolated SST'
389         fbdata%caddunit(1,1) = 'Degree centigrade'
390         fbdata%cgrid(1)      = 'T'
391         DO ja = 1, iadd
392            fbdata%caddname(1+ja) = padd%cdname(ja)
393            fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1)
394            fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1)
395         END DO
396
397      CASE('ICECON')
398
399         CALL alloc_obfbdata( fbdata, 1, surfdata%nsurf, 1, &
400            &                 1 + iadd, iext, .TRUE. )
401
402         clfiletype = 'sicfb'
403         fbdata%cname(1)      = surfdata%cvars(1)
404         fbdata%coblong(1)    = 'Sea ice'
405         fbdata%cobunit(1)    = 'Fraction'
406         DO je = 1, iext
407            fbdata%cextname(je) = pext%cdname(je)
408            fbdata%cextlong(je) = pext%cdlong(je,1)
409            fbdata%cextunit(je) = pext%cdunit(je,1)
410         END DO
411         fbdata%caddlong(1,1) = 'Model interpolated ICE'
412         fbdata%caddunit(1,1) = 'Fraction'
413         fbdata%cgrid(1)      = 'T'
414         DO ja = 1, iadd
415            fbdata%caddname(1+ja) = padd%cdname(ja)
416            fbdata%caddlong(1+ja,1) = padd%cdlong(ja,1)
417            fbdata%caddunit(1+ja,1) = padd%cdunit(ja,1)
418         END DO
419
420      END SELECT
421
422      fbdata%caddname(1)   = 'Hx'
423
424      WRITE(clfname, FMT="(A,'_fdbk_',I4.4,'.nc')") TRIM(clfiletype), nproc
425
426      IF(lwp) THEN
427         WRITE(numout,*)
428         WRITE(numout,*)'obs_wri_surf :'
429         WRITE(numout,*)'~~~~~~~~~~~~~'
430         WRITE(numout,*)'Writing '//TRIM(surfdata%cvars(1))//' feedback file : ',TRIM(clfname)
431      ENDIF
432
433      ! Transform surf data structure into obfbdata structure
434      fbdata%cdjuldref = '19500101000000'
435      DO jo = 1, surfdata%nsurf
436         fbdata%plam(jo)      = surfdata%rlam(jo)
437         fbdata%pphi(jo)      = surfdata%rphi(jo)
438         WRITE(fbdata%cdtyp(jo),'(I4)') surfdata%ntyp(jo)
439         fbdata%ivqc(jo,:)    = 0
440         fbdata%ivqcf(:,jo,:) = 0
441         IF ( surfdata%nqc(jo) > 10 ) THEN
442            fbdata%ioqc(jo)    = 4
443            fbdata%ioqcf(1,jo) = 0
444            fbdata%ioqcf(2,jo) = surfdata%nqc(jo) - 10
445         ELSE
446            fbdata%ioqc(jo)    = surfdata%nqc(jo)
447            fbdata%ioqcf(:,jo) = 0
448         ENDIF
449         fbdata%ipqc(jo)      = 0
450         fbdata%ipqcf(:,jo)   = 0
451         fbdata%itqc(jo)      = 0
452         fbdata%itqcf(:,jo)   = 0
453         fbdata%cdwmo(jo)     = surfdata%cwmo(jo)
454         fbdata%kindex(jo)    = surfdata%nsfil(jo)
455         IF (ln_grid_global) THEN
456            fbdata%iobsi(jo,1) = surfdata%mi(jo)
457            fbdata%iobsj(jo,1) = surfdata%mj(jo)
458         ELSE
459            fbdata%iobsi(jo,1) = mig(surfdata%mi(jo))
460            fbdata%iobsj(jo,1) = mjg(surfdata%mj(jo))
461         ENDIF
462         CALL greg2jul( 0, &
463            &           surfdata%nmin(jo), &
464            &           surfdata%nhou(jo), &
465            &           surfdata%nday(jo), &
466            &           surfdata%nmon(jo), &
467            &           surfdata%nyea(jo), &
468            &           fbdata%ptim(jo),   &
469            &           krefdate = 19500101 )
470         fbdata%padd(1,jo,1,1) = surfdata%rmod(jo,1)
471         IF ( TRIM(surfdata%cvars(1)) == 'SLA' ) fbdata%padd(1,jo,2,1) = surfdata%rext(jo,1)
472         fbdata%pob(1,jo,1)    = surfdata%robs(jo,1) 
473         fbdata%pdep(1,jo)     = 0.0
474         fbdata%idqc(1,jo)     = 0
475         fbdata%idqcf(:,1,jo)  = 0
476         IF ( surfdata%nqc(jo) > 10 ) THEN
477            fbdata%ivqc(jo,1)       = 4
478            fbdata%ivlqc(1,jo,1)    = 4
479            fbdata%ivlqcf(1,1,jo,1) = 0
480            fbdata%ivlqcf(2,1,jo,1) = surfdata%nqc(jo) - 10
481         ELSE
482            fbdata%ivqc(jo,1)       = surfdata%nqc(jo)
483            fbdata%ivlqc(1,jo,1)    = surfdata%nqc(jo)
484            fbdata%ivlqcf(:,1,jo,1) = 0
485         ENDIF
486         fbdata%iobsk(1,jo,1)  = 0
487         IF ( TRIM(surfdata%cvars(1)) == 'SLA' ) fbdata%pext(1,jo,1) = surfdata%rext(jo,2)
488         DO ja = 1, iadd
489            fbdata%padd(1,jo,2+ja,1) = &
490               & surfdata%rext(jo,padd%ipoint(ja))
491         END DO
492         DO je = 1, iext
493            fbdata%pext(1,jo,1+je) = &
494               & surfdata%rext(jo,pext%ipoint(je))
495         END DO
496      END DO
497
498      ! Write the obfbdata structure
499      CALL write_obfbdata( clfname, fbdata )
500
501      ! Output some basic statistics
502      CALL obs_wri_stats( fbdata )
503
504      CALL dealloc_obfbdata( fbdata )
505
506   END SUBROUTINE obs_wri_surf
507
508   SUBROUTINE obs_wri_stats( fbdata )
509      !!-----------------------------------------------------------------------
510      !!
511      !!                     *** ROUTINE obs_wri_stats  ***
512      !!
513      !! ** Purpose : Output some basic statistics of the data being written out
514      !!
515      !! ** Method  :
516      !!
517      !! ** Action  :
518      !!
519      !!      ! 2014-08  (D. Lea) Initial version
520      !!-----------------------------------------------------------------------
521
522      !! * Arguments
523      TYPE(obfbdata) :: fbdata
524
525      !! * Local declarations
526      INTEGER :: jvar
527      INTEGER :: jo
528      INTEGER :: jk
529      INTEGER :: inumgoodobs
530      INTEGER :: inumgoodobsmpp
531      REAL(wp) :: zsumx
532      REAL(wp) :: zsumx2
533      REAL(wp) :: zomb
534     
535
536      IF (lwp) THEN
537         WRITE(numout,*) ''
538         WRITE(numout,*) 'obs_wri_stats :'
539         WRITE(numout,*) '~~~~~~~~~~~~~~~'
540      ENDIF
541
542      DO jvar = 1, fbdata%nvar
543         zsumx=0.0_wp
544         zsumx2=0.0_wp
545         inumgoodobs=0
546         DO jo = 1, fbdata%nobs
547            DO jk = 1, fbdata%nlev
548               IF ( ( fbdata%pob(jk,jo,jvar) < 9999.0 ) .AND. &
549                  & ( fbdata%pdep(jk,jo) < 9999.0 ) .AND. &
550                  & ( fbdata%padd(jk,jo,1,jvar) < 9999.0 ) ) THEN
551
552                  zomb=fbdata%pob(jk, jo, jvar)-fbdata%padd(jk, jo, 1, jvar)
553                  zsumx=zsumx+zomb
554                  zsumx2=zsumx2+zomb**2
555                  inumgoodobs=inumgoodobs+1
556               ENDIF
557            ENDDO
558         ENDDO
559
560         CALL obs_mpp_sum_integer( inumgoodobs, inumgoodobsmpp )
561         CALL mpp_sum(zsumx)
562         CALL mpp_sum(zsumx2)
563
564         IF (lwp) THEN
565            WRITE(numout,*) 'Type: ',fbdata%cname(jvar),'  Total number of good observations: ',inumgoodobsmpp 
566            WRITE(numout,*) 'Overall mean obs minus model of the good observations: ',zsumx/inumgoodobsmpp
567            WRITE(numout,*) 'Overall RMS obs minus model of the good observations: ',sqrt( zsumx2/inumgoodobsmpp )
568            WRITE(numout,*) ''
569         ENDIF
570
571      ENDDO
572
573   END SUBROUTINE obs_wri_stats
574
575END MODULE obs_write
Note: See TracBrowser for help on using the repository browser.