Subversion Repositories lagranto.icon

Rev

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