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_vel.F90 in branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/OBS – NEMO

source: branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/OBS/obs_read_vel.F90 @ 3211

Last change on this file since 3211 was 3211, checked in by spickles2, 13 years ago

Stephen Pickles, 11 Dec 2011

Commit to bring the rest of the DCSE NEMO development branch
in line with the latest development version. This includes
array index re-ordering of all OPA_SRC/.

  • Property svn:keywords set to Id
File size: 25.0 KB
Line 
1MODULE obs_read_vel
2   !!======================================================================
3   !!                       ***  MODULE obs_read_vel  ***
4   !! Observation diagnostics: Read the velocity profile observations
5   !!======================================================================
6
7   !!----------------------------------------------------------------------
8   !!   obs_rea_vel_dri : Driver for reading profile obs
9   !!----------------------------------------------------------------------
10
11   !! * Modules used   
12   USE par_kind                 ! Precision variables
13   USE par_oce                  ! Ocean parameters
14   USE in_out_manager           ! I/O manager
15   USE dom_oce                  ! Ocean space and time domain variables
16   USE obs_mpp                  ! MPP support routines for observation diagnostics
17   USE julian                   ! Julian date routines
18   USE obs_utils                ! Observation operator utility functions
19   USE obs_prep                 ! Prepare observation arrays
20   USE obs_grid                 ! Grid search
21   USE obs_sort                 ! Sorting observation arrays
22   USE obs_profiles_def         ! Profile definitions
23   USE obs_conv                 ! Various conversion routines
24   USE obs_types                ! Observation type definitions
25   USE netcdf                   ! NetCDF library
26   USE obs_oper                 ! Observation operators
27   USE obs_vel_io               ! Velocity profile files I/O (non-FB files)
28   USE lib_mpp                  ! For ctl_warn/stop
29
30   IMPLICIT NONE
31
32   !! * Routine accessibility
33   PRIVATE
34
35   PUBLIC obs_rea_vel_dri  ! Read the profile observations
36
37   !! * Control permutation of array indices
38#  include "dom_oce_ftrans.h90"
39
40   !!----------------------------------------------------------------------
41   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
42   !! $Id$
43   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
44   !!----------------------------------------------------------------------
45
46CONTAINS
47 
48   SUBROUTINE obs_rea_vel_dri( kformat, &
49      &                        profdata, knumfiles, cfilenames, &
50      &                        kvars, kextr, kstp, ddobsini, ddobsend, &
51      &                        ldignmis, ldavtimset, ldmod )
52      !!---------------------------------------------------------------------
53      !!
54      !!                   *** ROUTINE obs_rea_pro_dri ***
55      !!
56      !! ** Purpose : Read from file the profile observations
57      !!
58      !! ** Method  : Depending on kformat either ENACT, CORIOLIS or
59      !!              feedback data files are read
60      !!
61      !! ** Action  :
62      !!
63      !! References :
64      !!
65      !! History : 
66      !!      ! :  2009-01 (K. Mogensen) : New merged version of old routines
67      !!----------------------------------------------------------------------
68      !! * Modules used
69   
70      !! * Arguments
71      INTEGER :: kformat   ! Format of input data
72      !                    ! 1: ENACT
73      !                    ! 2: Coriolis
74      TYPE(obs_prof), INTENT(OUT) :: profdata    ! Profile data to be read
75      INTEGER, INTENT(IN) :: knumfiles           ! Number of files to read in
76      CHARACTER(LEN=128), INTENT(IN) ::  cfilenames(knumfiles) ! File names to read in
77      INTEGER, INTENT(IN) :: kvars       ! Number of variables in profdata
78      INTEGER, INTENT(IN) :: kextr       ! Number of extra fields for each var in profdata
79      INTEGER, INTENT(IN) :: kstp        !  Ocean time-step index
80      LOGICAL, INTENT(IN) :: ldignmis    ! Ignore missing files
81      LOGICAL, INTENT(IN) :: ldavtimset  ! Set time to be equal to the end of the day
82      LOGICAL, INTENT(IN) :: ldmod       ! Initialize model from input data
83      REAL(KIND=dp), INTENT(IN) :: ddobsini   ! Obs. ini time in YYYYMMDD.HHMMSS
84      REAL(KIND=dp), INTENT(IN) :: ddobsend   ! Obs. end time in YYYYMMDD.HHMMSS
85
86      !! * Local declarations
87      CHARACTER(LEN=15), PARAMETER :: cpname='obs_rea_vel_dri'
88      INTEGER :: jvar
89      INTEGER :: ji
90      INTEGER :: jj
91      INTEGER :: jk
92      INTEGER :: ij
93      INTEGER :: iflag
94      INTEGER :: inobf
95      INTEGER :: i_file_id
96      INTEGER :: inowin
97      INTEGER :: iyea
98      INTEGER :: imon
99      INTEGER :: iday
100      INTEGER :: ihou
101      INTEGER :: imin
102      INTEGER :: isec
103      INTEGER, DIMENSION(knumfiles) :: &
104         & irefdate
105      INTEGER, DIMENSION(ntyp1770+1) :: &
106         & itypuv,    &
107         & itypuvmpp 
108      INTEGER :: iuv3dtmpp
109      INTEGER, DIMENSION(:), ALLOCATABLE :: &
110         & iobsiu,   &
111         & iobsju,   &
112         & iprocu,   &
113         & iobsiv,   &
114         & iobsjv,   &
115         & iprocv,   &
116         & iindx,    &
117         & ifileidx, &
118         & iprofidx
119      INTEGER :: itype
120      REAL(wp), DIMENSION(:), ALLOCATABLE :: &
121         & zphi, &
122         & zlam
123      REAL(dp), DIMENSION(:), ALLOCATABLE :: &
124         & zdat
125      LOGICAL :: &
126         & llvalprof
127      TYPE(obfbdata), POINTER, DIMENSION(:) :: &
128         & inpfiles
129      REAL(dp), DIMENSION(knumfiles) :: &
130         & djulini, &
131         & djulend
132      INTEGER :: iprof
133      INTEGER :: iproftot
134      INTEGER :: iuv3dt
135      INTEGER, DIMENSION(kvars) :: iv3dt
136      CHARACTER(len=8) :: cl_refdate
137   
138      ! Local initialization
139      iprof = 0
140      iuv3dt = 0
141
142      !-----------------------------------------------------------------------
143      ! Check data the model part is just with feedback data files
144      !-----------------------------------------------------------------------
145      IF ( ldmod .AND. ( kformat /= 0 ) ) THEN
146         CALL ctl_stop( 'Model can only be read from feedback data' )
147         RETURN
148      ENDIF
149
150      !-----------------------------------------------------------------------
151      ! Count the number of files needed and allocate the obfbdata type
152      !-----------------------------------------------------------------------
153     
154      inobf = knumfiles
155
156      ALLOCATE( inpfiles(inobf) )
157
158      prof_files : DO jj = 1, inobf
159         
160         !---------------------------------------------------------------------
161         ! Prints
162         !---------------------------------------------------------------------
163         IF(lwp) THEN
164            WRITE(numout,*)
165            WRITE(numout,*) ' obs_rea_vel_dri : Reading from file = ', &
166               & TRIM( TRIM( cfilenames(jj) ) )
167            WRITE(numout,*) ' ~~~~~~~~~~~~~~~'
168            WRITE(numout,*)
169         ENDIF
170
171         !---------------------------------------------------------------------
172         !  Initialization: Open file and get dimensions only
173         !---------------------------------------------------------------------
174         
175         iflag = nf90_open( TRIM( TRIM( cfilenames(jj) ) ), nf90_nowrite, &
176            &                      i_file_id )
177         
178         IF ( iflag /= nf90_noerr ) THEN
179
180            IF ( ldignmis ) THEN
181               inpfiles(jj)%nobs = 0
182               CALL ctl_warn( 'File ' // TRIM( TRIM( cfilenames(jj) ) ) // &
183                  &           ' not found' )
184            ELSE
185               CALL ctl_stop( 'File ' // TRIM( TRIM( cfilenames(jj) ) ) // &
186                  &           ' not found' )
187            ENDIF
188
189         ELSE 
190           
191            !------------------------------------------------------------------
192            !  Close the file since it is opened in read_proffile
193            !------------------------------------------------------------------
194           
195            iflag = nf90_close( i_file_id )
196
197            !------------------------------------------------------------------
198            !  Read the profile file into inpfiles
199            !------------------------------------------------------------------
200            IF ( kformat == 0 ) THEN
201               CALL init_obfbdata( inpfiles(jj) )
202               IF(lwp) THEN
203                  WRITE(numout,*)
204                  WRITE(numout,*)'Reading from feedback file :', &
205                     &           TRIM( cfilenames(jj) )
206               ENDIF
207               CALL read_obfbdata( TRIM( cfilenames(jj) ), inpfiles(jj), &
208                  &                ldgrid = .TRUE. )
209               IF ( ldmod .AND. ( inpfiles(jj)%nadd == 0 ) ) THEN
210                  CALL ctl_stop( 'Model not in input data' )
211                  RETURN
212               ENDIF
213            ELSEIF ( kformat == 1 ) THEN
214               CALL read_taondbc( TRIM( cfilenames(jj) ), inpfiles(jj), &
215                  &               numout, lwp, .TRUE. )
216            ELSE
217               CALL ctl_stop( 'File format unknown' )
218            ENDIF
219           
220            !------------------------------------------------------------------
221            !  Change longitude (-180,180)
222            !------------------------------------------------------------------
223
224            DO ji = 1, inpfiles(jj)%nobs 
225
226               IF ( inpfiles(jj)%plam(ji) < -180. ) &
227                  &   inpfiles(jj)%plam(ji) = inpfiles(jj)%plam(ji) + 360.
228
229               IF ( inpfiles(jj)%plam(ji) >  180. ) &
230                  &   inpfiles(jj)%plam(ji) = inpfiles(jj)%plam(ji) - 360.
231
232            END DO
233
234            !------------------------------------------------------------------
235            !  Calculate the date  (change eventually)
236            !------------------------------------------------------------------
237            cl_refdate=inpfiles(jj)%cdjuldref(1:8)
238            READ(cl_refdate,'(I8)') irefdate(jj)
239           
240            CALL ddatetoymdhms( ddobsini, iyea, imon, iday, ihou, imin, isec )
241            CALL greg2jul( isec, imin, ihou, iday, imon, iyea, djulini(jj), &
242               &           krefdate = irefdate(jj) )
243            CALL ddatetoymdhms( ddobsend, iyea, imon, iday, ihou, imin, isec )
244            CALL greg2jul( isec, imin, ihou, iday, imon, iyea, djulend(jj), &
245               &           krefdate = irefdate(jj) )
246
247            IF ( ldavtimset ) THEN
248               DO ji = 1, inpfiles(jj)%nobs
249                  !
250                  !  for daily averaged data force the time
251                  !  to be the  end of the day
252                  !
253                  inpfiles(jj)%ptim(ji) = &
254                     & INT(inpfiles(jj)%ptim(ji)) + 1
255               END DO
256            ENDIF
257           
258            IF ( inpfiles(jj)%nobs > 0 ) THEN
259               inpfiles(jj)%iproc = -1
260               inpfiles(jj)%iobsi = -1
261               inpfiles(jj)%iobsj = -1
262            ENDIF
263            inowin = 0
264            DO ji = 1, inpfiles(jj)%nobs
265               IF ( ( inpfiles(jj)%ptim(ji) >  djulini(jj) ) .AND. &
266                  & ( inpfiles(jj)%ptim(ji) <= djulend(jj) )       ) THEN
267                  inowin = inowin + 1
268               ENDIF
269            END DO
270            ALLOCATE( zlam(inowin)  )
271            ALLOCATE( zphi(inowin)  )
272            ALLOCATE( iobsiu(inowin) )
273            ALLOCATE( iobsju(inowin) )
274            ALLOCATE( iprocu(inowin) )
275            ALLOCATE( iobsiv(inowin) )
276            ALLOCATE( iobsjv(inowin) )
277            ALLOCATE( iprocv(inowin) )
278            inowin = 0
279            DO ji = 1, inpfiles(jj)%nobs
280               IF ( ( inpfiles(jj)%ptim(ji) >  djulini(jj) ) .AND. &
281                  & ( inpfiles(jj)%ptim(ji) <= djulend(jj) )       ) THEN
282                  inowin = inowin + 1
283                  zlam(inowin) = inpfiles(jj)%plam(ji)
284                  zphi(inowin) = inpfiles(jj)%pphi(ji)
285               ENDIF
286            END DO
287
288            CALL obs_grid_search( inowin, zlam, zphi, iobsiu, iobsju, iprocu, &
289               & 'U' )
290            CALL obs_grid_search( inowin, zlam, zphi, iobsiv, iobsjv, iprocv, &
291               & 'V' )
292
293            inowin = 0
294            DO ji = 1, inpfiles(jj)%nobs
295               IF ( ( inpfiles(jj)%ptim(ji) >  djulini(jj) ) .AND. &
296                  & ( inpfiles(jj)%ptim(ji) <= djulend(jj) )       ) THEN
297                  inowin = inowin + 1
298                  inpfiles(jj)%iproc(ji,1) = iprocu(inowin)
299                  inpfiles(jj)%iobsi(ji,1) = iobsiu(inowin)
300                  inpfiles(jj)%iobsj(ji,1) = iobsju(inowin)
301                  inpfiles(jj)%iproc(ji,2) = iprocv(inowin)
302                  inpfiles(jj)%iobsi(ji,2) = iobsiv(inowin)
303                  inpfiles(jj)%iobsj(ji,2) = iobsjv(inowin)
304                  IF ( inpfiles(jj)%iproc(ji,1) /= &
305                     & inpfiles(jj)%iproc(ji,2) ) THEN
306                     CALL ctl_stop( 'Error in obs_read_vel:', &
307                        & 'U and V observation on different processors')
308                  ENDIF
309               ENDIF
310            END DO
311            DEALLOCATE( zlam, zphi, iobsiu, iobsju, iprocu,  &
312               & iobsiv, iobsjv, iprocv )
313
314            DO ji = 1, inpfiles(jj)%nobs
315               IF ( ( inpfiles(jj)%ptim(ji) >  djulini(jj) ) .AND. &
316                  & ( inpfiles(jj)%ptim(ji) <= djulend(jj) )       ) THEN
317                  IF ( nproc == 0 ) THEN
318                     IF ( inpfiles(jj)%iproc(ji,1) >  nproc ) CYCLE
319                  ELSE
320                     IF ( inpfiles(jj)%iproc(ji,1) /=  nproc ) CYCLE
321                  ENDIF
322                  llvalprof = .FALSE.
323                  IF ( ( ( nproc == 0 ) .AND. &
324                     & ( inpfiles(jj)%iproc(ji,1) <=  nproc ) ) .OR. &
325                     & ( inpfiles(jj)%iproc(ji,1) ==  nproc ) ) THEN
326                     loop_uv_count : DO ij = 1,inpfiles(jj)%nlev
327                        IF ( inpfiles(jj)%pdep(ij,ji) >= 6000. ) &
328                           & CYCLE
329                        IF ( ( inpfiles(jj)%ivlqc(ij,ji,1) <= 2 ) .AND. &
330                           & ( inpfiles(jj)%ivlqc(ij,ji,2) <= 2 ) .AND. &
331                           & ( inpfiles(jj)%idqc(ij,ji) <= 2 ) ) THEN
332                           iuv3dt = iuv3dt + 1
333                           llvalprof = .TRUE.
334                        ENDIF
335                     END DO loop_uv_count
336                  ENDIF
337                  IF ( llvalprof ) iprof = iprof + 1
338               ENDIF
339            END DO
340
341         ENDIF
342         
343      END DO prof_files
344
345      !-----------------------------------------------------------------------
346      ! Get the time ordered indices of the input data
347      !-----------------------------------------------------------------------
348
349      !---------------------------------------------------------------------
350      !  Loop over input data files to count total number of profiles
351      !---------------------------------------------------------------------
352      iproftot = 0
353      DO jj = 1, inobf
354         DO ji = 1, inpfiles(jj)%nobs
355            IF ( ( inpfiles(jj)%ptim(ji) >  djulini(jj) ) .AND. &
356               & ( inpfiles(jj)%ptim(ji) <= djulend(jj) )       ) THEN
357               iproftot = iproftot + 1
358            ENDIF
359         END DO
360      END DO
361
362      ALLOCATE( iindx(iproftot), ifileidx(iproftot), &
363         &      iprofidx(iproftot), zdat(iproftot) )
364      jk = 0
365      DO jj = 1, inobf
366         DO ji = 1, inpfiles(jj)%nobs
367            IF ( ( inpfiles(jj)%ptim(ji) >  djulini(jj) ) .AND. &
368               & ( inpfiles(jj)%ptim(ji) <= djulend(jj) )       ) THEN
369               jk = jk + 1
370               ifileidx(jk) = jj
371               iprofidx(jk) = ji
372               zdat(jk)     = inpfiles(jj)%ptim(ji)
373            ENDIF
374         END DO
375      END DO
376      CALL sort_dp_indx( iproftot, &
377         &               zdat,     &
378         &               iindx   )
379     
380      iv3dt(:) = -1
381      iv3dt(1:2) = iuv3dt
382      CALL obs_prof_alloc( profdata, kvars, kextr, iprof, iv3dt, &
383         &                 kstp, jpi, jpj, jpk )
384     
385      ! * Read obs/positions, QC, all variable and assign to profdata
386
387      profdata%nprof     = 0
388      profdata%nvprot(:) = 0
389
390      iprof = 0
391
392      iuv3dt = 0
393      itypuv   (:) = 0
394      itypuvmpp(:) = 0
395     
396      DO jk = 1, iproftot
397         
398         jj = ifileidx(iindx(jk))
399         ji = iprofidx(iindx(jk))
400         IF ( ( inpfiles(jj)%ptim(ji) >  djulini(jj) ) .AND.  &
401            & ( inpfiles(jj)%ptim(ji) <= djulend(jj) ) ) THEN
402           
403            IF ( nproc == 0 ) THEN
404               IF ( inpfiles(jj)%iproc(ji,1) >  nproc ) CYCLE
405            ELSE
406               IF ( inpfiles(jj)%iproc(ji,1) /= nproc ) CYCLE
407            ENDIF
408           
409            llvalprof = .FALSE.
410
411            loop_prof : DO ij = 1, inpfiles(jj)%nlev
412               
413               IF ( inpfiles(jj)%pdep(ij,ji) >= 6000. ) &
414                  & CYCLE
415               
416               IF ( ( inpfiles(jj)%ivlqc(ij,ji,1) <= 2 ) .AND. &
417                  & ( inpfiles(jj)%ivlqc(ij,ji,2) <= 2 ) .AND. &
418                  & ( inpfiles(jj)%idqc(ij,ji) <= 2 ) ) THEN
419
420                  llvalprof = .TRUE.
421                  EXIT loop_prof
422
423               ENDIF
424
425            END DO loop_prof
426           
427            ! Set profile information
428           
429            IF ( llvalprof ) THEN
430               
431               iprof = iprof + 1
432
433               CALL jul2greg( isec,                   &
434                  &           imin,                   &
435                  &           ihou,                   &
436                  &           iday,                   &
437                  &           imon,                   &
438                  &           iyea,                   &
439                  &           inpfiles(jj)%ptim(ji), &
440                  &           irefdate(jj) )
441
442
443               ! Profile time coordinates
444               profdata%nyea(iprof) = iyea
445               profdata%nmon(iprof) = imon
446               profdata%nday(iprof) = iday
447               profdata%nhou(iprof) = ihou
448               profdata%nmin(iprof) = imin
449               
450               ! Profile space coordinates
451               profdata%rlam(iprof) = inpfiles(jj)%plam(ji)
452               profdata%rphi(iprof) = inpfiles(jj)%pphi(ji)
453
454               ! Coordinate search parameters
455               profdata%mi  (iprof,1) = inpfiles(jj)%iobsi(ji,1)
456               profdata%mj  (iprof,1) = inpfiles(jj)%iobsj(ji,1)
457               profdata%mi  (iprof,2) = inpfiles(jj)%iobsi(ji,2)
458               profdata%mj  (iprof,2) = inpfiles(jj)%iobsj(ji,2)
459               
460               ! Profile WMO number
461               profdata%cwmo(iprof) = inpfiles(jj)%cdwmo(ji)
462               
463               ! Instrument type
464               READ( inpfiles(jj)%cdtyp(ji), '(I4)' ) itype
465               profdata%ntyp(iprof) = itype
466               
467               ! QC stuff
468
469               profdata%nqc(iprof)     = inpfiles(jj)%ioqc(ji)
470               profdata%nqcf(:,iprof)  = inpfiles(jj)%ioqcf(:,ji)
471               profdata%ipqc(iprof)    = inpfiles(jj)%ipqc(ji)
472               profdata%ipqcf(:,iprof) = inpfiles(jj)%ipqcf(:,ji)
473               profdata%itqc(iprof)    = inpfiles(jj)%itqc(ji)
474               profdata%itqcf(:,iprof) = inpfiles(jj)%itqcf(:,ji)
475               profdata%ivqc(iprof,:)  = inpfiles(jj)%ivqc(ji,:)
476               profdata%ivqcf(:,iprof,:) = inpfiles(jj)%ivqcf(:,ji,:)
477
478               ! Bookkeeping data to match profiles
479               profdata%npidx(iprof) = iprof
480               profdata%npfil(iprof) = iindx(jk)
481
482               ! Observation QC flag (whole profile)
483               profdata%nqc(iprof)  = 0 !TODO
484
485               loop_uv : DO ij = 1, inpfiles(jj)%nlev           
486                   
487                  IF ( inpfiles(jj)%pdep(ij,ji) >= 6000. ) &
488                     & CYCLE
489
490                  IF ( ( inpfiles(jj)%ivlqc(ij,ji,1) <= 2 ) .AND. &
491                     & ( inpfiles(jj)%ivlqc(ij,ji,2) <= 2 ) .AND. &
492                     & ( inpfiles(jj)%idqc(ij,ji) <= 2 ) ) THEN
493                     iuv3dt = iuv3dt + 1
494                  ELSE
495                     CYCLE
496                  ENDIF
497
498                  ! Depth of U observation
499                  profdata%var(1)%vdep(iuv3dt) = &
500                     &                inpfiles(jj)%pdep(ij,ji)
501                 
502                  ! Depth of U observation QC
503                  profdata%var(1)%idqc(iuv3dt) = &
504                     &                inpfiles(jj)%idqc(ij,ji)
505                 
506                  ! Depth of U observation QC flags
507                  profdata%var(1)%idqcf(:,iuv3dt) = &
508                     &                inpfiles(jj)%idqcf(:,ij,ji)
509                 
510                  ! Profile index
511                  profdata%var(1)%nvpidx(iuv3dt) = iprof
512                 
513                  ! Vertical index in original profile
514                  profdata%var(1)%nvlidx(iuv3dt) = ij
515
516                  ! Profile U value
517
518                  profdata%var(1)%vobs(iuv3dt) = &
519                     &                inpfiles(jj)%pob(ij,ji,1)
520                  IF ( ldmod ) THEN
521                     profdata%var(1)%vmod(iuv3dt) = &
522                        &                inpfiles(jj)%padd(ij,ji,1,1)
523                  ENDIF
524                 
525                  ! Profile U qc
526                  profdata%var(1)%nvqc(iuv3dt) = &
527                     & inpfiles(jj)%ivlqc(ij,ji,1)
528                 
529                  ! Profile U qc flags
530                  profdata%var(1)%nvqcf(:,iuv3dt) = &
531                     & inpfiles(jj)%ivlqcf(:,ij,ji,1)
532                 
533                 
534                  ! Depth of V observation
535                  profdata%var(2)%vdep(iuv3dt) = &
536                     &                inpfiles(jj)%pdep(ij,ji)
537                 
538                  ! Depth of V observation QC
539                  profdata%var(2)%idqc(iuv3dt) = &
540                     &                inpfiles(jj)%idqc(ij,ji)
541                 
542                  ! Depth of V observation QC flags
543                  profdata%var(2)%idqcf(:,iuv3dt) = &
544                     &                inpfiles(jj)%idqcf(:,ij,ji)
545                 
546                  ! Profile index
547                  profdata%var(2)%nvpidx(iuv3dt) = iprof
548                 
549                  ! Vertical index in original profile
550                  profdata%var(2)%nvlidx(iuv3dt) = ij
551                 
552                  ! Profile V value
553                  profdata%var(2)%vobs(iuv3dt) = &
554                     &                inpfiles(jj)%pob(ij,ji,2)
555                  IF ( ldmod ) THEN
556                     profdata%var(2)%vmod(iuv3dt) = &
557                        &                inpfiles(jj)%padd(ij,ji,1,2)
558                  ENDIF
559                     
560                  ! Profile V qc
561                  profdata%var(2)%nvqc(iuv3dt) = &
562                     & inpfiles(jj)%ivlqc(ij,ji,2)
563                 
564                  ! Profile V qc flags
565                  profdata%var(2)%nvqcf(:,iuv3dt) = &
566                        & inpfiles(jj)%ivlqcf(:,ij,ji,2)
567                 
568                  ! Observation type
569                  itypuv( profdata%ntyp(iprof) + 1 ) = &
570                     & itypuv( profdata%ntyp(iprof) + 1 ) + 1
571
572
573               END DO loop_uv
574
575            ENDIF
576
577         ENDIF
578
579      END DO
580
581      !-----------------------------------------------------------------------
582      ! Sum up over processors
583      !-----------------------------------------------------------------------
584     
585      CALL obs_mpp_sum_integer ( iuv3dt, iuv3dtmpp )
586     
587      CALL obs_mpp_sum_integers( itypuv, itypuvmpp, ntyp1770 + 1 )
588     
589      !-----------------------------------------------------------------------
590      ! Output number of observations.
591      !-----------------------------------------------------------------------
592      IF(lwp) THEN
593         WRITE(numout,*) 
594         WRITE(numout,'(1X,A)') 'Profile U,V velocity data'
595         WRITE(numout,'(1X,A)') '-------------------------'
596         WRITE(numout,*) 
597         DO ji = 0, ntyp1770
598            IF ( itypuvmpp(ji+1) > 0 ) THEN
599               WRITE(numout,'(1X,A3,A3,I8)') ctypshort(ji), ' = ', &
600                  & itypuvmpp(ji+1)
601            ENDIF
602         END DO
603         WRITE(numout,'(1X,A)') '--------------'
604         WRITE(numout,'(1X,A6,I8)') &
605            & 'Total profile UV data                                 = ',&
606            & iuv3dtmpp
607         WRITE(numout,'(1X,A)') '--------------'
608      ENDIF
609     
610      profdata%nvprot(1)    = iuv3dt
611      profdata%nvprot(2)    = iuv3dt
612      profdata%nvprotmpp(1) = iuv3dtmpp
613      profdata%nvprotmpp(2) = iuv3dtmpp
614      profdata%nprof        = iprof
615
616      !-----------------------------------------------------------------------
617      ! Model level search
618      !-----------------------------------------------------------------------
619      CALL obs_level_search( jpk, gdept_0, &
620         & profdata%nvprot(1), profdata%var(1)%vdep, &
621         & profdata%var(1)%mvk )
622      CALL obs_level_search( jpk, gdept_0, &
623         & profdata%nvprot(2), profdata%var(2)%vdep, &
624         & profdata%var(2)%mvk )
625     
626      !-----------------------------------------------------------------------
627      ! Set model equivalent to 99999
628      !-----------------------------------------------------------------------
629      IF ( .NOT. ldmod ) THEN
630         DO jvar = 1, kvars
631            profdata%var(jvar)%vmod(:) = fbrmdi
632         END DO
633      ENDIF
634      !-----------------------------------------------------------------------
635      ! Deallocate temporary data
636      !-----------------------------------------------------------------------
637      DEALLOCATE( ifileidx, iprofidx, zdat )
638
639      !-----------------------------------------------------------------------
640      ! Deallocate input data
641      !-----------------------------------------------------------------------
642      DO jj = 1, inobf
643         CALL dealloc_obfbdata( inpfiles(jj) )
644      END DO
645      DEALLOCATE( inpfiles )
646
647   END SUBROUTINE obs_rea_vel_dri
648   
649END MODULE obs_read_vel
Note: See TracBrowser for help on using the repository browser.