Subversion Repositories lagranto.ecmwf

Rev

Rev 5 | Go to most recent revision | Details | Last modification | View Log | RSS feed

Rev Author Line No. Line
3 michaesp 1
      PROGRAM trace
2
 
3
C     ********************************************************************
4
C     *                                                                  *
5
C     * Pseudo-lidar plots along trajectories                             *
6
C     *                                                                  *
7
C     * Heini Wernli       first version:       April 1993               *
8
C     * Michael Sprenger   major upgrade:       2008-2009                *
9
C     *                                                                  *
10
C     ********************************************************************
11
 
12
      implicit none
13
 
14
c     --------------------------------------------------------------------
15
c     Declaration of parameters
16
c     --------------------------------------------------------------------
17
 
18
c     Maximum number of levels for input files
19
      integer   nlevmax
20
      parameter (nlevmax=100)
21
 
22
c     Maximum number of input files (dates, length of trajectories)
23
      integer   ndatmax
24
      parameter (ndatmax=500)
25
 
26
c     Numerical epsilon (for float comparison)
27
      real      eps
28
      parameter (eps=0.001)
29
 
30
c     Conversion factors
31
      real      pi180                                   ! deg -> rad
32
      parameter (pi180=3.14159/180.)
33
      real      deg2km                                  ! deg -> km (at equator)
34
      parameter (deg2km=111.)
35
 
36
c     Prefix for primary and secondary fields
37
      character charp
38
      character chars
39
      parameter (charp='P')
40
      parameter (chars='S')
41
 
42
c     --------------------------------------------------------------------
43
c     Declaration of variables
44
c     --------------------------------------------------------------------
45
 
46
c     Input and output format for trajectories (see iotra.f)
47
      integer   inpmode
48
 
49
c     Input parameters
50
      character*80                           inpfile         ! Input trajectory file
51
      character*80                           outfile         ! Output netCDF file
52
      character*80                           outmode         ! Output mode (sum,mean)
53
      integer                                ntra            ! Number of trajectories
54
      integer                                ncol            ! Number of columns (including time, lon, lat, p)
55
      integer                                ntim            ! Number of times per trajectory
56
      integer                                ntrace0         ! Number of trace variables
57
      character*80                           tvar(200)       ! Tracing variable name (only the variable)
58
      character*1                            tfil(200)       ! Filename prefix 
59
      real                                   fac(200)        ! Scaling factor 
60
      integer                                compfl(200)     ! Computation flag (1=compute)
61
      integer                                numdat          ! Number of input files
62
      character*11                           dat(ndatmax)    ! Dates of input files
63
      real                                   timeinc         ! Time increment between input files
64
      real                                   tst             ! Time shift of start relative to first data file
65
      real                                   ten             ! Time shift of end relatiev to first data file  
66
      character*20                           startdate       ! First time/date on trajectory
67
      character*20                           enddate         ! Last time/date on trajectory
68
      character*80                           timecheck       ! Either 'yes' or 'no'
69
      character*80                           intmode         ! Interpolation mode ('normal', 'nearest')
70
	  real                                   pmin,pmax       ! Pressure range for output grid
71
	  integer                                npre            ! Number of pressure levels in output grid
72
	  character*80                           centering       ! Centering around trajectory position ('yes','no')
73
 
74
c     Trajectories
75
      real,allocatable, dimension (:,:,:) :: trainp          ! Input trajectories (ntra,ntim,ncol)
76
      integer                                reftime(6)      ! Reference date
77
      character*80                           varsinp(100)    ! Field names for input trajectory
78
      integer                                fid,fod         ! File identifier for inp and out trajectories
79
      real                                   x0,y0,p0        ! Position of air parcel (physical space)
80
      real                                   reltpos0        ! Relative time of air parcel
81
      real                                   xind,yind,pind  ! Position of air parcel (grid space)
82
      integer                                fbflag          ! Flag for forward (1) or backward (-1) trajectories
83
 
84
c     Meteorological fields from input file
85
      real,allocatable, dimension (:)     :: spt0,spt1       ! Surface pressure
86
      real,allocatable, dimension (:)     :: p3t0,p3t1       ! 3d-pressure 
87
      real,allocatable, dimension (:)     :: f3t0,f3t1       ! 3d field for tracing 
88
      character*80                           svars(100)      ! List of variables on S file
89
      character*80                           pvars(100)      ! List of variables on P file
90
      integer                                n_svars         ! Number of variables on S file
91
      integer                                n_pvars         ! Number of variables on P file
92
 
93
c     Input grid description
94
      real                                   pollon,pollat   ! Longitude/latitude of pole
95
      real                                   ak(100)         ! Vertical layers and levels
96
      real                                   bk(100) 
97
      real                                   xmin,xmax       ! Zonal grid extension
98
      real                                   ymin,ymax       ! Meridional grid extension
99
      integer                                nx,ny,nz        ! Grid dimensions
100
      real                                   dx,dy           ! Horizontal grid resolution
101
      integer                                hem             ! Flag for hemispheric domain
102
      integer                                per             ! Flag for periodic domain
103
      real                                   stagz           ! Vertical staggering
104
      real                                   mdv             ! Missing data value
105
 
106
c	  Output grid  and fields
107
      real                                   levels(1000)    ! Ouput levels
108
      real                                   times (1000)    ! Output times
109
	  real,allocatable, dimension (:,:)   :: out_pos         ! Position of trajectories
110
	  real,allocatable, dimension (:,:)   :: out_val         ! Output lidar field
111
	  real,allocatable, dimension (:,:)   :: out_cnt         ! # output lidar sum ups
112
 
113
 
114
c     Auxiliary variables
115
      integer                                i,j,k,l,n
116
      real                                   rd
117
      character*80                           filename
118
      real                                   time0,time1,reltpos
119
      integer                                itime0,itime1
120
      integer                                stat
121
      real                                   tstart
122
      integer                                iloaded0,iloaded1
123
      real                                   f0
124
      real                                   frac
125
      real                                   tload,tfrac
126
      integer                                isok
127
      character                              ch
128
      integer                                ind
129
      integer                                ind1,ind2,ind3,ind4,ind5
130
      integer                                ind6,ind7,ind8,ind9,ind0
131
      integer                                noutside
132
      real                                   delta
133
      integer                                itrace0
134
      character*80                           string
135
      character*80                           cdfname
136
      character*80                           varname
137
      real                                   time
138
      character*80                           longname
139
      character*80                           unit
140
      integer                                ind_time
141
      integer                                ind_pre
142
 
143
c     Externals 
144
      real                                   int_index4
145
      external                               int_index4
146
 
147
c     --------------------------------------------------------------------
148
c     Start of program, Read parameters, get grid parameters
149
c     --------------------------------------------------------------------
150
 
151
c     Write start message
152
      print*,'========================================================='
153
      print*,'              *** START OF PROGRAM LIDAR ***'
154
      print*
155
 
156
c     Read parameters
157
      open(10,file='trace.param')
158
       read(10,*) inpfile
159
       read(10,*) outfile
160
       read(10,*) outmode
161
       read(10,*) startdate
162
       read(10,*) enddate 
163
       read(10,*) fbflag
164
       read(10,*) numdat
165
       if ( fbflag.eq.1) then
166
          do i=1,numdat
167
             read(10,'(a11)') dat(i)
168
          enddo
169
       else
170
          do i=numdat,1,-1
171
             read(10,'(a11)') dat(i)
172
          enddo
173
       endif
174
       read(10,*) timeinc
175
       read(10,*) tst
176
       read(10,*) ten
177
       read(10,*) ntra
178
       read(10,*) ntim
179
       read(10,*) ncol
180
       read(10,*) ntrace0
181
       do i=1,ntrace0
182
          read(10,*) tvar(i), fac(i), compfl(i), tfil(i)
183
       enddo
184
       read(10,*) n_pvars
185
       do i=1,n_pvars
186
          read(10,*) pvars(i)
187
       enddo
188
       read(10,*) n_svars
189
       do i=1,n_svars
190
          read(10,*) svars(i)
191
       enddo
192
       read(10,*) timecheck
193
       read(10,*) intmode
194
       read(10,*) pmin,pmax,npre
195
       read(10,*) centering
196
      close(10)
197
 
198
c     Remove commented tracing fields
199
      itrace0 = 1
200
      do while ( itrace0.le.ntrace0) 
201
         string = tvar(itrace0)
202
         if ( string(1:1).eq.'#' ) then
203
            do i=itrace0,ntrace0-1
204
               tvar(i)   = tvar(i+1)
205
               fac(i)    = fac(i+1)
206
               compfl(i) = compfl(i+1)
207
               tfil(i)   = tfil(i+1)
208
            enddo
209
            ntrace0 = ntrace0 - 1
210
         else
211
            itrace0 = itrace0 + 1
212
         endif
213
      enddo
214
 
215
c     Set the formats of the input  files
216
      call mode_tra(inpmode,inpfile)
217
      if (inpmode.eq.-1) inpmode=1
218
 
219
C     Convert time shifts <tst,ten> from <hh.mm> into fractional time
220
      call hhmm2frac(tst,frac)
221
      tst = frac
222
      call hhmm2frac(ten,frac)
223
      ten = frac
224
 
225
c     Set the time for the first data file (depending on forward/backward mode)
226
      if (fbflag.eq.1) then
227
        tstart = -tst
228
      else
229
        tstart = tst
230
      endif
231
 
232
c     Read the constant grid parameters (nx,ny,nz,xmin,xmax,ymin,ymax,pollon,pollat)
233
c     The negative <-fid> of the file identifier is used as a flag for parameter retrieval  
234
      filename = charp//dat(1)
235
      varname  = 'U'
236
      call input_open (fid,filename)
237
      call input_grid (-fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny,
238
     >                 tstart,pollon,pollat,rd,rd,nz,rd,rd,rd,timecheck)
239
      call input_close(fid)
240
 
241
C     Allocate memory for some meteorological arrays
242
      allocate(spt0(nx*ny),stat=stat)
243
      if (stat.ne.0) print*,'*** error allocating array spt0 ***'   ! Surface pressure
244
      allocate(spt1(nx*ny),stat=stat)
245
      if (stat.ne.0) print*,'*** error allocating array spt1 ***'
246
      allocate(p3t0(nx*ny*nz),stat=stat)
247
      if (stat.ne.0) print*,'*** error allocating array p3t0 ***'   ! Pressure
248
      allocate(p3t1(nx*ny*nz),stat=stat)
249
      if (stat.ne.0) print*,'*** error allocating array p3t1 ***'
250
      allocate(f3t0(nx*ny*nz),stat=stat)
251
      if (stat.ne.0) print*,'*** error allocating array p3t0 ***'   ! Lidar field
252
      allocate(f3t1(nx*ny*nz),stat=stat)
253
      if (stat.ne.0) print*,'*** error allocating array p3t1 ***'
254
 
255
c	  Allocate memory for output field
256
	  allocate(out_pos(ntim,npre),stat=stat)
257
      if (stat.ne.0) print*,'*** error allocating array out_pos ***'
258
	  allocate(out_val(ntim,npre),stat=stat)
259
      if (stat.ne.0) print*,'*** error allocating array out_val ***'
260
	  allocate(out_cnt(ntim,npre),stat=stat)
261
      if (stat.ne.0) print*,'*** error allocating array out_cnt ***'
262
 
263
C     Get memory for trajectory arrays
264
      allocate(trainp(ntra,ntim,ncol),stat=stat)
265
      if (stat.ne.0) print*,'*** error allocating array tra      ***' 
266
 
267
c     Set the flags for periodic domains
268
      if ( abs(xmax-xmin-360.).lt.eps ) then
269
         per = 1
270
      elseif ( abs(xmax-xmin-360.+dx).lt.eps ) then
271
         per = 2
272
      else
273
         per = 0
274
      endif
275
 
276
C     Set logical flag for periodic data set (hemispheric or not)
277
      hem = 0
278
      if (per.eq.0.) then
279
         delta=xmax-xmin-360.
280
         if (abs(delta+dx).lt.eps) then               ! Program aborts: arrays must be closed
281
            print*,' ERROR: arrays must be closed... Stop'
282
         else if (abs(delta).lt.eps) then ! Periodic and hemispheric
283
           hem=1
284
           per=360.
285
        endif
286
      else                                            ! Periodic and hemispheric
287
         hem=1
288
      endif
289
 
290
c     Write some status information
291
      print*,'---- INPUT PARAMETERS -----------------------------------'
292
      print*
293
      print*,'  Input trajectory file  : ',trim(inpfile)
294
      print*,'  Format of input file   : ',inpmode
295
      print*,'  Output netCDF    file  : ',trim(outfile)
296
      print*,'  Format of output file  : ',trim(outmode)
297
      print*,'  Forward/backward       : ',fbflag
298
      print*,'  #tra                   : ',ntra
299
      print*,'  #col                   : ',ncol
300
      print*,'  #tim                   : ',ntim
301
      print*,'  No time check          : ',trim(timecheck)
302
      print*,'  Interpolation mode     : ',trim(intmode)
303
      do i=1,ntrace0
304
         if (compfl(i).eq.0) then
305
            print*,'  Tracing field          : ',
306
     >                 trim(tvar(i)), fac(i), ' 0 ', tfil(i)
307
         else
308
            print*,'  Tracing field          : ',
309
     >                trim(tvar(i)),' : online calc not supported'
310
         endif
311
      enddo
312
      print*,'  Output (pmin,pmax,n)   : ',pmin,pmax,npre
313
      print*,'  Centering              : ',trim(centering)
314
      print*
315
      print*,'---- INPUT DATA FILES -----------------------------------'
316
      print*
317
      call frac2hhmm(tstart,tload)
318
      print*,'  Time of 1st data file  : ',tload
319
      print*,'  #input files           : ',numdat
320
      print*,'  time increment         : ',timeinc
321
      call frac2hhmm(tst,tload)
322
      print*,'  Shift of start         : ',tload
323
      call frac2hhmm(ten,tload)
324
      print*,'  Shift of end           : ',tload
325
      print*,'  First/last input file  : ',trim(dat(1)),
326
     >                                     ' ... ',
327
     >                                     trim(dat(numdat))
328
      print*,'  Primary variables      : ',trim(pvars(1))
329
      do i=2,n_pvars
330
         print*,'                         : ',trim(pvars(i))
331
      enddo
332
      if ( n_svars.ge.1 ) then
333
         print*,'  Secondary variables    : ',trim(svars(1))
334
         do i=2,n_svars
335
            print*,'                         : ',trim(svars(i))
336
         enddo
337
      endif
338
      print*
339
      print*,'---- CONSTANT GRID PARAMETERS ---------------------------'
340
      print*
341
      print*,'  xmin,xmax     : ',xmin,xmax
342
      print*,'  ymin,ymax     : ',ymin,ymax
343
      print*,'  dx,dy         : ',dx,dy
344
      print*,'  pollon,pollat : ',pollon,pollat
345
      print*,'  nx,ny,nz      : ',nx,ny,nz
346
      print*,'  per, hem      : ',per,hem
347
      print*
348
 
349
c     --------------------------------------------------------------------
350
c     Load the input trajectories
351
c     --------------------------------------------------------------------
352
 
353
c     Read the input trajectory file
354
      call ropen_tra(fid,inpfile,ntra,ntim,ncol,reftime,varsinp,inpmode)
355
      call read_tra (fid,trainp,ntra,ntim,ncol,inpmode)
356
      call close_tra(fid,inpmode)
357
 
358
c     Check that first four columns correspond to time,lon,lat,p
359
      if ( (varsinp(1).ne.'time' ).or.
360
     >     (varsinp(2).ne.'xpos' ).and.(varsinp(2).ne.'lon' ).or.
361
     >     (varsinp(3).ne.'ypos' ).and.(varsinp(3).ne.'lat' ).or.
362
     >     (varsinp(4).ne.'ppos' ).and.(varsinp(4).ne.'p'   ) )
363
     >then
364
         print*,' ERROR: problem with input trajectories ...'
365
         stop
366
      endif
367
      varsinp(1) = 'time'
368
      varsinp(2) = 'lon'
369
      varsinp(3) = 'lat'
370
      varsinp(4) = 'p'
371
 
372
c     Write some status information of the input trajectories
373
      print*,'---- INPUT TRAJECTORIES ---------------------------------'
374
      print*
375
      print*,' Start date             : ',trim(startdate)
376
      print*,' End date               : ',trim(enddate)
377
      print*,' Reference time (year)  : ',reftime(1)
378
      print*,'                (month) : ',reftime(2)
379
      print*,'                (day)   : ',reftime(3)
380
      print*,'                (hour)  : ',reftime(4)
381
      print*,'                (min)   : ',reftime(5)
382
      print*,' Time range (min)       : ',reftime(6)
383
      do i=1,ncol
384
         print*,' Var                    :',i,trim(varsinp(i))
385
      enddo
386
      print*
387
 
388
c     Check that first time is 0 - otherwise the tracing will produce
389
c     wrong results because later in the code only absolute times are
390
c     considered: <itime0   = int(abs(tfrac-tstart)/timeinc) + 1>. This 
391
c     will be changed in a future version.
392
      if ( abs( trainp(1,1,1) ).gt.eps ) then
393
         print*,' ERROR: First time of trajectory must be 0, i.e. '
394
         print*,'     correspond to the reference date. Otherwise'
395
         print*,'     the tracing will give wrong results... STOP'
396
         stop
397
      endif
398
 
399
c     --------------------------------------------------------------------
400
c     Trace the fields (fields available on input files)
401
c     --------------------------------------------------------------------
402
 
403
      print*
404
      print*,'---- LIDAR FROM PRIMARY AND SECONDARY DATA FILES ------'
405
 
406
c     Loop over all tracing fields
407
      do i=1,ntrace0
408
 
409
c	      Skip all fields marked for online calculation
410
          if ( compfl(i).eq.1 ) goto 110
411
 
412
c	      Init the output fields: position and lidar field
413
	      do k=1,ntim
414
	      	do l=1,npre
415
	      		out_pos(k,l) = 0.
416
	      	    out_val(k,l) = 0.
417
	      	    out_cnt(k,l) = 0.
418
	      	 enddo
419
	      enddo
420
 
421
c         Write some status information
422
          print*
423
          print*,' Now lidaring           : ',
424
     >         trim(tvar(i)),compfl(i),' ',trim(tfil(i))
425
 
426
c         Reset flags for load manager
427
          iloaded0 = -1
428
          iloaded1 = -1
429
 
430
c         Reset the counter for fields outside domain
431
          noutside = 0
432
 
433
c         Loop over all times
434
          do j=1,ntim
435
 
436
c            Convert trajectory time from hh.mm to fractional time
437
             call hhmm2frac(trainp(1,j,1),tfrac)
438
 
439
c            Get the times which are needed
440
             itime0   = int(abs(tfrac-tstart)/timeinc) + 1
441
             time0    = tstart + fbflag * real(itime0-1) * timeinc
442
             itime1   = itime0 + 1
443
             time1    = time0 + fbflag * timeinc
444
             if ( itime1.gt.numdat ) then
445
                itime1 = itime0
446
                time1  = time0
447
             endif
448
 
449
c            Load manager: Check whether itime0 can be copied from itime1
450
             if ( itime0.eq.iloaded1 ) then
451
                f3t0     = f3t1
452
                p3t0     = p3t1
453
                spt0     = spt1
454
                iloaded0 = itime0
455
             endif
456
 
457
c            Load manager: Check whether itime1 can be copied from itime0
458
             if ( itime1.eq.iloaded0 ) then
459
                f3t1     = f3t0
460
                p3t1     = p3t0
461
                spt1     = spt0
462
                iloaded1 = itime1
463
             endif
464
 
465
c            Load manager:  Load first time (tracing variable and grid)
466
             if ( itime0.ne.iloaded0 ) then
467
 
468
                filename = tfil(i)//dat(itime0)
469
                call frac2hhmm(time0,tload)
470
                varname  = tvar(i) 
471
                write(*,'(a23,a20,a3,a5,f7.2)') 
472
     >               '    ->  loading          : ',
473
     >               trim(filename),'  ',trim(varname),tload
474
                call input_open (fid,filename)                
475
                call input_wind 
476
     >               (fid,varname,f3t0,tload,stagz,mdv,
477
     >               xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck) 
478
 
479
                call input_grid      
480
     >               (fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny,
481
     >                tload,pollon,pollat,p3t0,spt0,nz,ak,bk,stagz,
482
     >                timecheck)
483
                call input_close(fid)
484
 
485
                iloaded0 = itime0
486
 
487
             endif
488
 
489
c            Load manager: Load second time (tracing variable and grid)
490
             if ( itime1.ne.iloaded1 ) then
491
 
492
                filename = tfil(i)//dat(itime1)
493
                call frac2hhmm(time1,tload)
494
                varname  = tvar(i) 
495
                write(*,'(a23,a20,a3,a5,f7.2)') 
496
     >               '    ->  loading          : ',
497
     >               trim(filename),'  ',trim(varname),tload
498
                call input_open (fid,filename)
499
                call input_wind 
500
     >               (fid,varname,f3t1,tload,stagz,mdv,
501
     >               xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)      
502
                call input_grid      
503
     >               (fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny,
504
     >               tload,pollon,pollat,p3t1,spt1,nz,ak,bk,stagz,
505
     >               timecheck)
506
                call input_close(fid)
507
 
508
                iloaded1 = itime1
509
 
510
             endif
511
 
512
c            Loop over all trajectories
513
             do k=1,ntra
514
 
515
c               Set the horizontal position where to interpolate to
516
                x0       = trainp(k,j,2)                          ! Longitude
517
                y0       = trainp(k,j,3)                          ! Latitude
518
 
519
c               Set the relative time
520
                call hhmm2frac(trainp(k,j,1),tfrac)
521
                reltpos0 = fbflag * (tfrac-time0)/timeinc    
522
 
523
c               Handle periodic boundaries in zonal direction
524
                if ( (x0.gt.xmax).and.(per.ne.0) ) x0 = x0 - 360.
525
                if ( (x0.lt.xmin).and.(per.ne.0) ) x0 = x0 + 360.
526
 
527
c               Handle pole problems for hemispheric data (taken from caltra.f)
528
                if ((hem.eq.1).and.(y0.gt.90.)) then
529
                   y0=180.-y0
530
                   x0=x0+per/2.
531
                endif
532
                if ((hem.eq.1).and.(y0.lt.-90.)) then
533
                   y0=-180.-y0
534
                   x0=x0+per/2.
535
                endif
536
                if (y0.gt.89.99) then
537
                   y0=89.99
538
                endif
539
 
540
c               Loop over pressure profile
541
	            do l=1,npre
542
 
543
c                 Set the pressure where to interpolate to
544
	              p0 = pmin + real(l-1)/real(npre-1) * (pmax-pmin)
545
	              if ( centering.eq.'yes' )then
546
	              	  p0 = p0 + trainp(k,j,4)
547
	              endif
548
 
549
C                 Get the index where to interpolate (x0,y0,p0)
550
                  if ( (abs(x0-mdv).gt.eps).and.
551
     >                 (abs(y0-mdv).gt.eps) )
552
     >            then
553
                     call get_index4 (xind,yind,pind,x0,y0,p0,reltpos0,
554
     >                                p3t0,p3t1,spt0,spt1,3,
555
     >                                nx,ny,nz,xmin,ymin,dx,dy,mdv)
556
                  else
557
                     xind = mdv
558
                     yind = mdv
559
                     pind = mdv
560
                  endif
561
 
562
c                 If requested, apply nearest-neighbor interpolation
563
                  if ( intmode.eq.'nearest') then
564
 
565
                     xind = real( nint(xind) )
566
                     yind = real( nint(yind) )
567
                     pind = real( nint(pind) )
568
 
569
                     if ( xind.lt.1.  ) xind = 1.
570
                     if ( xind.gt.nx  ) xind = real(nx)
571
                     if ( yind.lt.1.  ) yind = 1.
572
                     if ( yind.gt.ny  ) yind = real(ny)
573
 
574
                     if ( pind.lt.1.  ) pind = 1.
575
                     if ( pind.gt.nz  ) pind = real(nz)
576
 
577
                  endif
578
 
579
c                 Do the interpolation: everthing is ok
580
                  if ( (xind.ge.1.).and.(xind.le.real(nx)).and.
581
     >                 (yind.ge.1.).and.(yind.le.real(ny)).and.
582
     >                 (pind.ge.1.).and.(pind.le.real(nz)) )
583
     >            then
584
                     f0 = int_index4(f3t0,f3t1,nx,ny,nz,
585
     >                               xind,yind,pind,reltpos0,mdv)
586
 
587
c                 Set to missing data
588
                  else
589
                     f0       = mdv
590
                  endif
591
 
592
c	              Save result to output array
593
                  if (abs(f0-mdv).gt.eps) then
594
                     out_val(j,l) = out_val(j,l) + f0 * fac(i)
595
                     out_cnt(j,l) = out_cnt(j,l) + 1.
596
 
597
	              endif
598
 
599
c              End loop over all pressure levels
600
	           enddo
601
 
602
c	           Save the output trajectory position
603
	           ind_time = j
604
	           if ( centering.eq.'no' ) then
605
	              ind_pre  = nint( real(npre) *
606
     >          	       ( (trainp(k,j,4) - pmin)/(pmax-pmin) ) + 1.)
607
	           else
608
	           	  ind_pre  = nint( real(npre) *
609
     >          	       ( (0.            - pmin)/(pmax-pmin) ) + 1.)
610
	           endif
611
 
612
	           if ( (ind_time.ge.1).and.(ind_time.le.ntim).and.
613
     >              (ind_pre .ge.1).and.(ind_pre .le.npre) )
614
     >         then
615
                    out_pos(ind_time,ind_pre) =
616
     >                          	out_pos(ind_time,ind_pre) + 1.
617
	           endif
618
 
619
c	         End loop over all trajectories
620
             enddo
621
 
622
c	      End loop over all times
623
          enddo
624
 
625
c	      Write the trajectory position to netCDF file - only once
626
	      if ( i.eq.1 ) then
627
	      	  cdfname  = outfile
628
	      	  varname  = 'POSITION'
629
	      	  longname = 'position of trajectory points'
630
	      	  unit     = 'none'
631
	      	  time     = 0.
632
              do k=1,npre
633
              	levels(k) = pmin + real(k-1)/real(npre-1) * (pmax-pmin)
634
              enddo
635
              do k=1,ntim
636
                 times(k) = trainp(1,k,1)
637
              enddo
638
              call writecdf2D_cf
639
     >            (cdfname,varname,longname,unit,out_pos,time,levels,
640
     >             times,npre,ntim,1,1)
641
	      endif
642
 
643
c	      If no valid lidar count: set the field to missing data
644
          do k=1,ntim
645
          	do l=1,npre
646
          		if (abs(out_cnt(k,l)).lt.eps) then
647
          			out_val(k,l) = mdv
648
          	    endif
649
          	 enddo
650
          enddo
651
 
652
c	      If requested, calculate the mean of the lidar field
653
	      if ( outmode.eq.'mean' ) then
654
	      	do k=1,ntim
655
	      		do l=1,npre
656
	      			if ( (abs(out_val(k,l)-mdv).gt.eps).and.
657
     >                   (abs(out_cnt(k,l)    ).gt.0. ) )
658
     >              then
659
	      				out_val(k,l) = out_val(k,l) / out_cnt(k,l)
660
	      		    endif
661
	      		 enddo
662
	          enddo
663
	      endif
664
 
665
c	      Write the lidar field and count
666
	      cdfname  = outfile
667
	      if (outmode.eq.'sum' ) then
668
	         varname  = trim(tvar(i))//'_SUM'
669
	      elseif (outmode.eq.'mean' ) then
670
	         varname  = trim(tvar(i))//'_MEAN'
671
	      endif
672
	      longname = 'sum over all '//trim(tvar(i))//' profiles'
673
	      unit     = 'not given'
674
	      time     = 0.
675
          call writecdf2D_cf
676
     >            (cdfname,varname,longname,unit,out_val,time,levels,
677
     >             times,npre,ntim,0,1)
678
 
679
	  	  cdfname  = outfile
680
	      varname  = trim(tvar(i))//'_CNT'
681
	      longname = 'counts of all '//trim(tvar(i))//' profiles'
682
	      unit     = 'not given'
683
	      time     = 0.
684
          call writecdf2D_cf
685
     >            (cdfname,varname,longname,unit,out_cnt,time,levels,
686
     >             times,npre,ntim,0,1)
687
 
688
c         Exit point for loop over all tracing variables
689
 110      continue
690
 
691
c	   End loop over all lidar variables
692
       enddo
693
 
694
 
695
c     --------------------------------------------------------------------
696
c     Write output to netCDF file
697
c     --------------------------------------------------------------------
698
 
699
c     Write status information
700
      print*
701
      print*,'---- WRITE OUTPUT LIDAR FIELDS --------------------------'
702
      print*
703
 
704
 
705
c     Write some status information, and end of program message
706
      print*  
707
      print*,'---- STATUS INFORMATION --------------------------------'
708
      print*
709
      print*,' ok'
710
      print*
711
      print*,'              *** END OF PROGRAM LIDAR ***'
712
      print*,'========================================================='
713
 
714
 
715
      end 
716
 
717
 
718
c     ********************************************************************
719
c     * INPUT / OUTPUT SUBROUTINES                                       *
720
c     ********************************************************************
721
 
722
c     --------------------------------------------------------------------
723
c     Subroutines to write 2D CF netcdf output file
724
c     --------------------------------------------------------------------
725
 
726
      subroutine writecdf2D_cf
727
     >          (cdfname,varname,longname,unit,arr,time,levels,times,
728
     >           npre,ntim,crefile,crevar)
729
 
730
c     Create and write to the CF netcdf file <cdfname>. The variable
731
c     with name <varname> and with time <time> is written. The data
732
c     are in the two-dimensional array <arr>.  The flags <crefile> and
733
c     <crevar> determine whether the file and/or the variable should
734
c     be created.
735
 
736
      USE netcdf
737
 
738
      IMPLICIT NONE
739
 
740
c     Declaration of input parameters
741
      character*80 cdfname
742
      character*80 varname,longname,unit
743
      integer      npre,ntim
744
      real         arr(ntim,npre)
745
      real         levels(npre)
746
      real         times (ntim)
747
      real         time
748
      integer      crefile,crevar
749
 
750
c     Numerical epsilon
751
      real         eps
752
      parameter    (eps=1.e-5)
753
 
754
c     Local variables
755
      integer      ierr
756
      integer      ncID
757
      integer      LevDimId,    varLevID
758
      integer      TimeDimID,   varTimeID
759
      real         timeindex
760
      integer      i,j
761
      integer      nvars,varids(100)
762
      integer      ndims,dimids(100)
763
      real         timelist(1000)
764
      integer      ntimes
765
      integer      ind
766
      integer      varID
767
 
768
c     Quick an dirty solution for fieldname conflict
769
      if ( varname.eq.'time' ) varname = 'TIME'
770
 
771
c     Initially set error to indicate no errors.
772
      ierr = 0
773
 
774
c     ---- Create the netCDF - skip if <crefile=0> ----------------------
775
      if ( crefile.ne.1 ) goto 100
776
 
777
c     Create the file
778
      ierr = nf90_create(trim(cdfname), NF90_CLOBBER, ncID)
779
 
780
c     Define dimensions
781
      ierr=nf90_def_dim(ncID,'level',npre, LevDimID )
782
      ierr=nf90_def_dim(ncID,'time' ,ntim, TimeDimID)
783
 
784
c     Define coordinate Variables
785
      ierr = nf90_def_var(ncID,'level',NF90_FLOAT,
786
     >     (/ LevDimID /),varLevID)
787
      ierr = nf90_put_att(ncID, varLevID, "standard_name","level")
788
      ierr = nf90_put_att(ncID, varLevID, "units"        ,"hPa")
789
 
790
      ierr = nf90_def_var(ncID,'time',NF90_FLOAT,
791
     >     (/ TimeDimID /), varTimeID)
792
      ierr = nf90_put_att(ncID, varTimeID, "long_name",    "time")
793
      ierr = nf90_put_att(ncID, varTimeID, "units",       "hours")
794
 
795
c     Write global attributes
796
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'Conventions', 'CF-1.0')
797
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'title',
798
     >     'pseudo-lidar from trajectory file')
799
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'source',
800
     >     'Lagranto Trajectories')
801
      ierr = nf90_put_att(ncID, NF90_GLOBAL, 'institution',
802
     >     'ETH Zurich, IACETH')
803
 
804
c     Check whether the definition was successful
805
      ierr = nf90_enddef(ncID)
806
      if (ierr.gt.0) then
807
         print*, 'An error occurred while attempting to ',
808
     >        'finish definition mode.'
809
         stop
810
      endif
811
 
812
c     Write coordinate data
813
      ierr = nf90_put_var(ncID,varLevID  ,levels)
814
      ierr = nf90_put_var(ncID,varTimeID ,times )
815
 
816
c     Close netCDF file
817
      ierr = nf90_close(ncID)
818
 
819
 100  continue
820
 
821
c     ---- Define a new variable - skip if <crevar=0> -----------------------
822
 
823
      if ( crevar.ne.1 ) goto 110
824
 
825
c     Open the file for read/write access
826
      ierr = nf90_open  (trim(cdfname), NF90_WRITE  , ncID)
827
 
828
c     Get the IDs for dimensions
829
      ierr = nf90_inq_dimid(ncID,'level', LevDimID )
830
      ierr = nf90_inq_dimid(ncID,'time' , TimeDimID)
831
 
832
c     Enter define mode
833
      ierr = nf90_redef(ncID)
834
 
835
c     Write definition and add attributes
836
      ierr = nf90_def_var(ncID,varname,NF90_FLOAT,
837
     >                   (/ varTimeID, varLevID /),varID)
838
      ierr = nf90_put_att(ncID, varID, "long_name" , longname )
839
      ierr = nf90_put_att(ncID, varID, "units"     , unit     )
840
      ierr = nf90_put_att(ncID, varID, '_FillValue', -999.99  )
841
 
842
c     Check whether definition was successful
843
      ierr = nf90_enddef(ncID)
844
      if (ierr.gt.0) then
845
         print*, 'An error occurred while attempting to ',
846
     >           'finish definition mode.'
847
         stop
848
      endif
849
 
850
c     Close netCDF file
851
      ierr = nf90_close(ncID)
852
 
853
 110  continue
854
 
855
c     ---- Write data --------------------------------------------------
856
 
857
c     Open the file for read/write access
858
      ierr = nf90_open  (trim(cdfname), NF90_WRITE , ncID)
859
 
860
c     Get the varID
861
      ierr = nf90_inq_varid(ncID,varname, varID )
862
      if (ierr.ne.0) then
863
         print*,'Variable ',trim(varname),' is not defined on ',
864
     >          trim(cdfname)
865
         stop
866
      endif
867
 
868
c     Write data block
869
      ierr = nf90_put_var(ncID,varID,arr,
870
     >                    start = (/ 1, 1 /),
871
     >                    count = (/ ntim, npre/) )
872
 
873
c     Check whether writing was successful
874
      ierr = nf90_close(ncID)
875
      if (ierr.ne.0) then
876
         write(*,*) trim(nf90_strerror(ierr))
877
         write(*,*) 'An error occurred while attempting to ',
878
     >              'close the netcdf file.'
879
         write(*,*) 'in clscdf_CF'
880
      endif
881
 
882
      end
883
 
884
 
885
 
886
 
887
 
888
 
889
 
890
 
891