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_read_sla.F90 in branches/2012/dev_NOC_MERCATOR_2012/NEMOGCM/NEMO/OPA_SRC/OBS – NEMO

source: branches/2012/dev_NOC_MERCATOR_2012/NEMOGCM/NEMO/OPA_SRC/OBS/obs_read_sla.F90 @ 3651

Last change on this file since 3651 was 3651, checked in by cbricaud, 11 years ago

merge dev_MERCATOR_2012_rev3555 into dev_NOC_MERCATOR_2012 ; see ticket 1020

  • Property svn:keywords set to Id
File size: 20.5 KB
Line 
1MODULE obs_read_sla
2   !!======================================================================
3   !!                       ***  MODULE obs_read_sla  ***
4   !! Observation diagnostics: Read the along track SLA data from
5   !!                          the AVISO/CLS database or any SLA data
6   !!                          from feedback files
7   !!======================================================================
8
9   !!----------------------------------------------------------------------
10   !!   obs_rea_sla : Driver for reading SLA data from the AVISO/CLS/feedback
11   !!----------------------------------------------------------------------
12
13   !! * Modules used   
14   USE par_kind                 ! Precision variables
15   USE in_out_manager           ! I/O manager
16   USE dom_oce                  ! Ocean space and time domain variables
17   USE obs_mpp                  ! MPP support routines for observation diagnostics
18   USE julian                   ! Julian date routines
19   USE obs_utils                ! Observation operator utility functions
20   USE obs_grid                 ! Grid search
21   USE obs_sort                 ! Sorting observation arrays
22   USE obs_surf_def             ! Surface observation definitions
23   USE obs_types                ! Observation type definitions
24   USE obs_sla_io               ! I/O for sla files
25   USE netcdf                   ! NetCDF library
26
27   IMPLICIT NONE
28
29   !! * Routine accessibility
30   PRIVATE
31
32   PUBLIC obs_rea_sla  ! Read the SLA observations from the AVISO/SLA database
33
34   !!----------------------------------------------------------------------
35   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
36   !! $Id$
37   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
38   !!----------------------------------------------------------------------
39
40CONTAINS
41
42   SUBROUTINE obs_rea_sla( kformat, &
43      &                    sladata, knumfiles, cfilenames, &
44      &                    kvars, kextr, kstp, ddobsini, ddobsend, &
45      &                    ldignmis, ldmod, ldobstd )
46      !!---------------------------------------------------------------------
47      !!
48      !!                   *** ROUTINE obs_rea_sla ***
49      !!
50      !! ** Purpose : Read from file the SLA data
51      !!
52      !! ** Method  : Depending on kformat either AVISO or
53      !!              feedback data files are read
54      !!
55      !! ** Action  :
56      !!
57      !!
58      !! History : 
59      !!      ! :  2009-01 (K. Mogensen) Initial version based on old versions
60      !!----------------------------------------------------------------------
61      !! * Modules used
62
63      !! * Arguments
64      INTEGER :: kformat    ! Format of input data
65      !                     ! 0: Feedback
66      !                     ! 1: AVISO
67      TYPE(obs_surf), INTENT(INOUT) :: sladata    ! SLA data to be read
68      INTEGER, INTENT(IN) :: knumfiles      ! Number of files to read in
69      CHARACTER(LEN=128), INTENT(IN) ::  &
70         & cfilenames(knumfiles) ! File names to read in
71      INTEGER, INTENT(IN) :: kvars     ! Number of variables in sladata
72      INTEGER, INTENT(IN) :: kextr     ! Number of extra fields for each var in sladata
73      INTEGER, INTENT(IN) :: kstp      ! Ocean time-step index
74      LOGICAL, INTENT(IN) :: ldignmis    ! Ignore missing files
75      LOGICAL, INTENT(IN) :: ldmod       ! Initialize model from input data
76      LOGICAL, INTENT(INOUT), optional :: &
77         & ldobstd        ! Read observation standard deviation from fb. file
78      REAL(KIND=dp), INTENT(IN) :: ddobsini  ! Obs. ini time in YYYYMMDD.HHMMSS
79      REAL(KIND=dp), INTENT(IN) :: ddobsend  ! Obs. end time in YYYYMMDD.HHMMSS
80         
81      !! * Local declarations
82      CHARACTER(LEN=11), PARAMETER :: cpname='obs_rea_sla'
83      INTEGER :: ji
84      INTEGER :: jj
85      INTEGER :: jk
86      INTEGER :: iflag
87      INTEGER :: inobf
88      INTEGER :: i_file_id
89      INTEGER :: inowin
90      INTEGER :: iyea
91      INTEGER :: imon
92      INTEGER :: iday
93      INTEGER :: ihou
94      INTEGER :: imin
95      INTEGER :: isec
96      INTEGER, DIMENSION(knumfiles) :: &
97         & irefdate
98      INTEGER :: iobsmpp
99      INTEGER, DIMENSION(imaxmissions+1) :: &
100         & ityp, &
101         & itypmpp
102      INTEGER, DIMENSION(:), ALLOCATABLE :: &
103         & iobsi,    &
104         & iobsj,    &
105         & iproc,    &
106         & iindx,    &
107         & ifileidx, &
108         & islaidx
109      INTEGER :: itype 
110      REAL(wp), DIMENSION(:), ALLOCATABLE :: &
111         & zphi, &
112         & zlam
113      REAL(dp), DIMENSION(:), ALLOCATABLE :: &
114         & zdat
115      LOGICAL :: llvalprof
116      LOGICAL :: llobstd
117      TYPE(obfbdata), POINTER, DIMENSION(:) :: &
118         & inpfiles
119      REAL(dp), DIMENSION(knumfiles) :: &
120         & djulini, &
121         & djulend
122      INTEGER, DIMENSION(knumfiles) :: &
123         & iobspos, &
124         & iobcpos
125      INTEGER :: iobs
126      INTEGER :: iobstot
127      CHARACTER(len=8) :: cl_refdate
128   
129      ! Local initialization
130      iobs = 0
131      IF ( PRESENT(ldobstd) ) THEN
132         IF (.NOT.ldmod) THEN
133            llobstd = .false.
134         ELSE
135            llobstd = ldobstd
136         ENDIF
137      ELSE
138         llobstd = .FALSE.
139      ENDIF
140 
141      !-----------------------------------------------------------------------
142      ! Check that the model part is just with feedback data files
143      !-----------------------------------------------------------------------
144      IF ( ldmod .AND. ( kformat /= 0 ) ) THEN
145         CALL ctl_stop( 'Model can only be read from feedback data' )
146         RETURN
147      ENDIF
148
149      !-----------------------------------------------------------------------
150      ! Check that the prescribed obs err is just with feedback data files
151      !-----------------------------------------------------------------------
152      IF ( llobstd .AND. ( kformat /= 0 ) ) THEN
153         CALL ctl_stop( 'Observation error can only be read from feedback files' )
154         RETURN
155      ENDIF
156
157      !-----------------------------------------------------------------------
158      ! Count the number of files needed and allocate the obfbdata type
159      !-----------------------------------------------------------------------
160     
161      inobf = knumfiles
162     
163      ALLOCATE( inpfiles(inobf) )
164
165      sla_files : DO jj = 1, inobf
166         
167         !---------------------------------------------------------------------
168         ! Prints
169         !---------------------------------------------------------------------
170         IF(lwp) THEN
171            WRITE(numout,*)
172            WRITE(numout,*) ' obs_rea_sla : Reading from file = ', &
173               & TRIM( TRIM( cfilenames(jj) ) )
174            WRITE(numout,*) ' ~~~~~~~~~~~'
175            WRITE(numout,*)
176         ENDIF
177
178         !---------------------------------------------------------------------
179         !  Initialization: Open file and get dimensions only
180         !---------------------------------------------------------------------
181         
182         iflag = nf90_open( TRIM( TRIM( cfilenames(jj) ) ), nf90_nowrite, &
183            &                      i_file_id )
184         
185         IF ( iflag /= nf90_noerr ) THEN
186
187            IF ( ldignmis ) THEN
188               inpfiles(jj)%nobs = 0
189               CALL ctl_warn( 'File ' // TRIM( TRIM( cfilenames(jj) ) ) // &
190                  &           ' not found' )
191            ELSE
192               CALL ctl_stop( 'File ' // TRIM( TRIM( cfilenames(jj) ) ) // &
193                  &           ' not found' )
194            ENDIF
195
196         ELSE 
197           
198            !------------------------------------------------------------------
199            !  Close the file since it is opened in read_proffile
200            !------------------------------------------------------------------
201           
202            iflag = nf90_close( i_file_id )
203
204            !------------------------------------------------------------------
205            !  Read the profile file into inpfiles
206            !------------------------------------------------------------------
207            IF ( kformat == 0 ) THEN
208               CALL init_obfbdata( inpfiles(jj) )
209               IF(lwp) THEN
210                  WRITE(numout,*)
211                  WRITE(numout,*)'Reading from feedback file :', &
212                     &           TRIM( cfilenames(jj) )
213               ENDIF
214               CALL read_obfbdata( TRIM( cfilenames(jj) ), inpfiles(jj), &
215                  &                ldgrid = .TRUE. )
216               IF ( inpfiles(jj)%nvar < 1 ) THEN
217                  CALL ctl_stop( 'Feedback format error' )
218                  RETURN
219               ENDIF
220               IF ( TRIM(inpfiles(jj)%cname(1)) /= 'SLA' ) THEN
221                  CALL ctl_stop( 'Feedback format error' )
222                  RETURN
223               ENDIF
224               IF ( ldmod .AND. ( ( inpfiles(jj)%nadd == 0 ) .OR.&
225                  &               ( inpfiles(jj)%next == 0 ) ) ) THEN
226                  CALL ctl_stop( 'Model not in input data' )
227                  RETURN
228               ENDIF
229            ELSEIF ( kformat == 1) THEN
230               CALL read_avisofile( TRIM( cfilenames(jj) ), inpfiles(jj), &
231               &                    numout, lwp, .TRUE. )
232            ELSE
233               CALL ctl_stop( 'File format unknown' )
234            ENDIF
235
236            !------------------------------------------------------------------
237            !  Change longitude (-180,180)
238            !------------------------------------------------------------------
239
240            DO ji = 1, inpfiles(jj)%nobs 
241
242               IF ( inpfiles(jj)%plam(ji) < -180. ) &
243                  &   inpfiles(jj)%plam(ji) = inpfiles(jj)%plam(ji) + 360.
244
245               IF ( inpfiles(jj)%plam(ji) >  180. ) &
246                  &   inpfiles(jj)%plam(ji) = inpfiles(jj)%plam(ji) - 360.
247
248            END DO
249
250            !------------------------------------------------------------------
251            !  Calculate the date  (change eventually)
252            !------------------------------------------------------------------
253            cl_refdate=inpfiles(jj)%cdjuldref(1:8)
254            READ(cl_refdate,'(I8)') irefdate(jj)
255           
256            CALL ddatetoymdhms( ddobsini, iyea, imon, iday, ihou, imin, isec )
257            CALL greg2jul( isec, imin, ihou, iday, imon, iyea, djulini(jj), &
258               &           krefdate = irefdate(jj) )
259            CALL ddatetoymdhms( ddobsend, iyea, imon, iday, ihou, imin, isec )
260            CALL greg2jul( isec, imin, ihou, iday, imon, iyea, djulend(jj), &
261               &           krefdate = irefdate(jj) )
262            IF ( inpfiles(jj)%nobs > 0 ) THEN
263               inpfiles(jj)%iproc = -1
264               inpfiles(jj)%iobsi = -1
265               inpfiles(jj)%iobsj = -1
266            ENDIF
267            inowin = 0
268            DO ji = 1, inpfiles(jj)%nobs
269               IF ( inpfiles(jj)%ioqc(ji) > 2 ) CYCLE
270               IF ( inpfiles(jj)%ivqc(ji,1) > 2 ) CYCLE
271               IF ( ( inpfiles(jj)%ptim(ji) >  djulini(jj) ) .AND. &
272                  & ( inpfiles(jj)%ptim(ji) <= djulend(jj) )       ) THEN
273                  inowin = inowin + 1
274               ENDIF
275            END DO
276            ALLOCATE( zlam(inowin)  )
277            ALLOCATE( zphi(inowin)  )
278            ALLOCATE( iobsi(inowin) )
279            ALLOCATE( iobsj(inowin) )
280            ALLOCATE( iproc(inowin) )
281            inowin = 0
282            DO ji = 1, inpfiles(jj)%nobs
283               IF ( inpfiles(jj)%ioqc(ji) > 2 ) CYCLE
284               IF ( inpfiles(jj)%ivqc(ji,1) > 2 ) CYCLE
285               IF ( ( inpfiles(jj)%ptim(ji) >  djulini(jj) ) .AND. &
286                  & ( inpfiles(jj)%ptim(ji) <= djulend(jj) )       ) THEN
287                  inowin = inowin + 1
288                  zlam(inowin) = inpfiles(jj)%plam(ji)
289                  zphi(inowin) = inpfiles(jj)%pphi(ji)
290               ENDIF
291            END DO
292
293            CALL obs_grid_search( inowin, zlam, zphi, iobsi, iobsj, iproc, 'T' )
294
295            inowin = 0
296            DO ji = 1, inpfiles(jj)%nobs
297               IF ( inpfiles(jj)%ioqc(ji) > 2 ) CYCLE
298               IF ( inpfiles(jj)%ivqc(ji,1) > 2 ) CYCLE
299               IF ( ( inpfiles(jj)%ptim(ji) >  djulini(jj) ) .AND. &
300                  & ( inpfiles(jj)%ptim(ji) <= djulend(jj) )       ) THEN
301                  inowin = inowin + 1
302                  inpfiles(jj)%iproc(ji,1) = iproc(inowin)
303                  inpfiles(jj)%iobsi(ji,1) = iobsi(inowin)
304                  inpfiles(jj)%iobsj(ji,1) = iobsj(inowin)
305               ENDIF
306            END DO
307            DEALLOCATE( zlam, zphi, iobsi, iobsj, iproc )
308
309            DO ji = 1, inpfiles(jj)%nobs
310               IF ( inpfiles(jj)%ioqc(ji) > 2 ) CYCLE
311               IF ( inpfiles(jj)%ivqc(ji,1) > 2 ) CYCLE
312               IF ( ( inpfiles(jj)%ptim(ji) >  djulini(jj) ) .AND. &
313                  & ( inpfiles(jj)%ptim(ji) <= djulend(jj) )       ) THEN
314                  IF ( nproc == 0 ) THEN
315                     IF ( inpfiles(jj)%iproc(ji,1) >  nproc ) CYCLE
316                  ELSE
317                     IF ( inpfiles(jj)%iproc(ji,1) /= nproc ) CYCLE
318                  ENDIF
319                  llvalprof = .FALSE.
320                  IF ( ( inpfiles(jj)%ivlqc(1,ji,1) == 1 ) .OR. &
321                     & ( inpfiles(jj)%ivlqc(1,ji,1) == 2 ) ) THEN
322                     iobs = iobs + 1
323                  ENDIF
324               ENDIF
325            END DO
326
327         ENDIF
328     
329      END DO sla_files
330
331      IF (llobstd) THEN
332         
333         DO jj = 1, inobf
334            iobspos(jj) = -1
335            iobcpos(jj) = -1
336            DO ji = 1,inpfiles(jj)%nadd
337               IF ( TRIM(inpfiles(jj)%caddname(ji)) == 'OSTD' ) THEN
338                  iobspos(jj)=ji
339               ENDIF
340               IF ( TRIM(inpfiles(jj)%caddname(ji)) == 'OCNT' ) THEN
341                  iobcpos(jj)=ji
342               ENDIF
343            END DO
344         END DO
345         llobstd = ( ( MINVAL(iobspos) > 0 ) .AND. ( MINVAL(iobcpos) > 0 ) )
346         IF (llobstd) THEN
347            IF (lwp) WRITE(numout,*)'SLA superobs information present.'
348         ELSE
349            IF (lwp) WRITE(numout,*)'SLA superobs information not present.'
350         ENDIF
351
352      ENDIF
353     
354      !-----------------------------------------------------------------------
355      ! Get the time ordered indices of the input data
356      !-----------------------------------------------------------------------
357
358      !---------------------------------------------------------------------
359      !  Loop over input data files to count total number of profiles
360      !---------------------------------------------------------------------
361      iobstot = 0
362      DO jj = 1, inobf
363         DO ji = 1, inpfiles(jj)%nobs
364            IF ( inpfiles(jj)%ioqc(ji) > 2 ) CYCLE
365            IF ( inpfiles(jj)%ivqc(ji,1) > 2 ) CYCLE
366            IF ( ( inpfiles(jj)%ptim(ji) >  djulini(jj) ) .AND. &
367               & ( inpfiles(jj)%ptim(ji) <= djulend(jj) )       ) THEN
368               iobstot = iobstot + 1
369            ENDIF
370         END DO
371      END DO
372
373      ALLOCATE( iindx(iobstot), ifileidx(iobstot), &
374         &      islaidx(iobstot), zdat(iobstot) )
375      jk = 0
376      DO jj = 1, inobf
377         DO ji = 1, inpfiles(jj)%nobs
378            IF ( inpfiles(jj)%ioqc(ji) > 2 ) CYCLE
379            IF ( inpfiles(jj)%ivqc(ji,1) > 2 ) CYCLE
380            IF ( ( inpfiles(jj)%ptim(ji) >  djulini(jj) ) .AND. &
381               & ( inpfiles(jj)%ptim(ji) <= djulend(jj) )       ) THEN
382               jk = jk + 1
383               ifileidx(jk) = jj
384               islaidx(jk) = ji
385               zdat(jk)     = inpfiles(jj)%ptim(ji)
386            ENDIF
387         END DO
388      END DO
389      CALL sort_dp_indx( iobstot, &
390         &               zdat,     &
391         &               iindx   )
392     
393      CALL obs_surf_alloc( sladata, iobs, kvars, kextr, &
394         &                 jpi, jpj, kstp )
395     
396      ! * Read obs/positions, QC, all variable and assign to sladata
397 
398      iobs = 0
399
400      ityp   (:) = 0
401      itypmpp(:) = 0
402     
403      DO jk = 1, iobstot
404         
405         jj = ifileidx(iindx(jk))
406         ji = islaidx(iindx(jk))
407
408         IF ( inpfiles(jj)%ioqc(ji) > 2 ) CYCLE
409         IF ( inpfiles(jj)%ivqc(ji,1) > 2 ) CYCLE
410
411         IF ( ( inpfiles(jj)%ptim(ji) >  djulini(jj) ) .AND.  &
412            & ( inpfiles(jj)%ptim(ji) <= djulend(jj) ) ) THEN
413           
414            IF ( nproc == 0 ) THEN
415               IF ( inpfiles(jj)%iproc(ji,1) >  nproc ) CYCLE
416            ELSE
417               IF ( inpfiles(jj)%iproc(ji,1) /= nproc ) CYCLE
418            ENDIF
419           
420            ! Set observation information
421           
422            IF ( ( inpfiles(jj)%ivlqc(1,ji,1) == 1 ) .OR. &
423               & ( inpfiles(jj)%ivlqc(1,ji,1) == 2 ) ) THEN
424
425               iobs = iobs + 1
426
427               CALL jul2greg( isec,                   &
428                  &           imin,                   &
429                  &           ihou,                   &
430                  &           iday,                   &
431                  &           imon,                   &
432                  &           iyea,                   &
433                  &           inpfiles(jj)%ptim(ji), &
434                  &           irefdate(jj) )
435
436
437               ! SLA time coordinates
438               sladata%nyea(iobs) = iyea
439               sladata%nmon(iobs) = imon
440               sladata%nday(iobs) = iday
441               sladata%nhou(iobs) = ihou
442               sladata%nmin(iobs) = imin
443               
444               ! SLA space coordinates
445               sladata%rlam(iobs) = inpfiles(jj)%plam(ji)
446               sladata%rphi(iobs) = inpfiles(jj)%pphi(ji)
447
448               ! Coordinate search parameters
449               sladata%mi  (iobs) = inpfiles(jj)%iobsi(ji,1)
450               sladata%mj  (iobs) = inpfiles(jj)%iobsj(ji,1)
451               
452               ! Instrument type
453               READ( inpfiles(jj)%cdtyp(ji), '(I4)' ) itype
454               sladata%ntyp(iobs) = itype
455               ityp(itype+1) = ityp(itype+1) + 1
456
457               ! Identifier
458               sladata%cwmo(iobs) = inpfiles(jj)%cdwmo(ji)
459
460               ! Bookkeeping data to match observations
461               sladata%nsidx(iobs) = iobs
462               sladata%nsfil(iobs) = iindx(jk)
463
464               ! QC flags
465               sladata%nqc(iobs) = inpfiles(jj)%ivqc(ji,1)
466
467               ! Observed value
468               sladata%robs(iobs,1) = inpfiles(jj)%pob(1,ji,1)
469
470
471               ! Model and MDT is set to fbrmdi unless read from file
472               IF ( ldmod ) THEN
473                  sladata%rmod(iobs,1) = inpfiles(jj)%padd(1,ji,1,1)
474                  sladata%rext(iobs,1) = inpfiles(jj)%padd(1,ji,2,1)
475                  sladata%rext(iobs,2) = inpfiles(jj)%pext(1,ji,1)
476                  IF (llobstd) THEN
477                     sladata%rext(iobs,3) = &
478                        & inpfiles(jj)%padd(1,ji,iobspos(jj),1)
479                     sladata%rext(iobs,4) = &
480                        & inpfiles(jj)%padd(1,ji,iobcpos(jj),1)
481                  ENDIF
482               ELSE
483                  sladata%rmod(iobs,1) = fbrmdi
484                  sladata%rext(iobs,:) = fbrmdi
485               ENDIF
486
487            ENDIF
488         ENDIF
489
490      END DO
491
492      !-----------------------------------------------------------------------
493      ! Sum up over processors
494      !-----------------------------------------------------------------------
495     
496      CALL obs_mpp_sum_integer( iobs, iobsmpp )
497      CALL obs_mpp_sum_integers( ityp, itypmpp, imaxmissions + 1 )
498     
499      !-----------------------------------------------------------------------
500      ! Output number of observations.
501      !-----------------------------------------------------------------------
502      IF (lwp) THEN
503
504         WRITE(numout,*)
505         WRITE(numout,'(1X,A)')'Altimeter satellites'
506         WRITE(numout,'(1X,A)')'--------------------'
507         DO jj = 1,8
508            IF ( itypmpp(jj) > 0 ) THEN
509               WRITE(numout,'(1X,A38,A2,I10)')calttyp(jj-1),'= ',itypmpp(jj)
510            ENDIF
511         END DO
512         WRITE(numout,'(1X,A50)')'--------------------------------------------------'
513         WRITE(numout,'(1X,A40,I10)')'Total                                 = ',iobsmpp
514         WRITE(numout,*)
515
516      ENDIF
517
518      !-----------------------------------------------------------------------
519      ! Deallocate temporary data
520      !-----------------------------------------------------------------------
521      DEALLOCATE( ifileidx, islaidx, zdat )
522
523      !-----------------------------------------------------------------------
524      ! Deallocate input data
525      !-----------------------------------------------------------------------
526      DO jj = 1, inobf
527         CALL dealloc_obfbdata( inpfiles(jj) )
528      END DO
529      DEALLOCATE( inpfiles )
530
531      !-----------------------------------------------------------------------
532      ! Reset ldobstd if the data is present
533      !-----------------------------------------------------------------------
534      IF ( PRESENT(ldobstd) ) THEN
535         ldobstd = llobstd
536      ENDIF
537
538   END SUBROUTINE obs_rea_sla
539
540END MODULE obs_read_sla
Note: See TracBrowser for help on using the repository browser.