Subversion Repositories lagranto.icon

Rev

Details | Last modification | View Log | RSS feed

Rev Author Line No. Line
3 michaesp 1
      PROGRAM caltra
2
 
3
C     ********************************************************************
4
C     *                                                                  *
5
C     * Calculates 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     Distance in m between 2 lat circles 
31
      real	deltay
32
      parameter	(deltay=1.112E5)
33
 
34
c     Numerical method for the integration (0=iterative Euler, 1=Runge-Kutta)
35
      integer   imethod
36
      parameter (imethod=1)
37
 
38
c     Number of iterations for iterative Euler scheme
39
      integer   numit
40
      parameter (numit=3)
41
 
42
c     Input and output format for trajectories (see iotra.f)
43
      integer   inpmode
44
      integer   outmode
45
 
46
c     Filename prefix (typically 'P')
47
      character*1 prefix
48
      parameter   (prefix='P')
49
 
50
c     --------------------------------------------------------------------
51
c     Declaration of variables
52
c     --------------------------------------------------------------------
53
 
54
c     Input parameters
55
      integer                                fbflag          ! Flag for forward/backward mode
56
      integer                                numdat          ! Number of input files
57
      character*13                           dat(ndatmax)    ! Dates of input files
58
      real                                   timeinc         ! Time increment between input files
59
      real                                   per             ! Periodicity (=0 if none)
60
      integer                                ntra            ! Number of trajectories
61
      character*80                           cdfname         ! Name of output files
62
      real                                   ts              ! Time step
63
      real                                   tst,ten         ! Shift of start and end time relative to first data file
64
      integer                                deltout         ! Output time interval (in minutes)
65
      integer                                jflag           ! Jump flag (if =1 ground-touching trajectories reenter atmosphere)
66
      real                                   wfactor         ! Factor for vertical velocity field
67
      character*80                           strname         ! File with start positions
68
      character*80                           timecheck       ! Either 'yes' or 'no'
69
      real                                   balloon         ! vertical velocity for balloon sounding
70
 
71
c     Trajectories
72
      integer                                ncol            ! Number of columns for insput trajectories
73
      real,allocatable, dimension (:,:,:) :: trainp          ! Input start coordinates (ntra,1,ncol)
74
      real,allocatable, dimension (:,:,:) :: traout          ! Output trajectories (ntra,ntim,4)
75
      integer                                reftime(6)      ! Reference date
76
      character*80                           vars(200)       ! Field names
77
      real,allocatable, dimension (:)     :: xx0,yy0,pp0     ! Position of air parcels
78
      integer,allocatable, dimension (:)  :: leftflag        ! Flag for domain-leaving
79
      real                                   xx1,yy1,pp1     ! Updated position of air parcel
80
      integer                                leftcount       ! Number of domain leaving trajectories
81
      integer                                ntim            ! Number of output time steps
82
 
83
c     Meteorological fields
84
      real,allocatable, dimension (:)     :: uut0,uut1       ! Zonal wind
85
      real,allocatable, dimension (:)     :: vvt0,vvt1       ! Meridional wind
86
      real,allocatable, dimension (:)     :: wwt0,wwt1       ! Vertical wind
87
 
88
c     Grid description
89
      real,allocatable, dimension (:)     :: zuv,zw         ! 3d-model height for U, V, W grid
90
      real,allocatable, dimension (:)     :: zb             ! Surface height
91
 
92
      real                                   pollon,pollat   ! Longitude/latitude of pole
93
      real                                   polgam          ! Rotation of grid
94
      real                                   xmin,xmax       ! Zonal grid extension
95
      real                                   ymin,ymax       ! Meridional grid extension
96
      integer                                nx,ny,nz        ! Grid dimensions
97
      real                                   dx,dy           ! Horizontal grid resolution
98
      integer                                hem             ! Flag for hemispheric domain
99
      real                                   mdv             ! Missing data value
100
      integer                                fflag           ! Fplane flag (0/1)
101
      real                                   lat_fixed       ! fixed latitude for fplane
102
 
103
c     Auxiliary variables                 
104
      real                                   delta,rd
105
      integer	                             itm,iloop,i,j,k,filo,lalo
106
      integer                                ierr,stat
107
      integer                                cdfid,fid
108
      real	                                 tstart,time0,time1,time
109
      real                                   reltpos0,reltpos1
110
      real                                   xind,yind,pind,pp,sp,stagz
111
      character*80                           filename,varname
112
      integer                                reftmp(6)
113
      character                              ch
114
      real                                   frac,tload
115
      integer                                itim
116
      real                                   lon,lat,rlon,rlat
117
 
118
c     Externals 
119
      real                                   lmtolms        ! Grid rotation
120
      real                                   phtophs    
121
      real                                   lmstolm
122
      real                                   phstoph        
123
      external                               lmtolms,phtophs
124
      external                               lmstolm,phstoph
125
      real                                   phi2phirot
126
      real                                   phirot2phi
127
      real                                   rla2rlarot
128
      real                                   rlarot2rla
129
      external                               phi2phirot,phirot2phi
130
      external                               rla2rlarot,rlarot2rla
131
 
132
c     --------------------------------------------------------------------
133
c     Start of program, Read parameters
134
c     --------------------------------------------------------------------
135
 
136
c     Write start message
137
      print*,'========================================================='
138
      print*,'              *** START OF PROGRAM CALTRA ***'
139
      print*
140
 
141
c     Open the parameter file
142
      open(9,file='caltra.param')
143
 
144
c     Read flag for forward/backward mode (fbflag)
145
      read(9,*) fbflag
146
 
147
c     Read number of input files (numdat)
148
      read(9,*) numdat
149
      if (numdat.gt.ndatmax) then
150
        print*,' ERROR: too many input files ',numdat,ndatmax
151
        goto 993
152
      endif
153
 
154
c     Read list of input dates (dat, sort depending on forward/backward mode)
155
      if (fbflag.eq.1) then
156
        do itm=1,numdat
157
          read(9,'(a)') dat(itm)
158
        enddo
159
      else
160
        do itm=numdat,1,-1
161
          read(9,'(a)') dat(itm)
162
        enddo
163
      endif
164
 
165
c     Read time increment between input files (timeinc)
166
      read(9,*) timeinc
167
 
168
C     Read if data domain is periodic and its periodicity
169
      read(9,*) per
170
 
171
c     Read the number of trajectories and name of position file
172
      read(9,*) strname
173
      read(9,*) ntra
174
      read(9,*) ncol 
175
      if (ntra.eq.0) goto 991
176
 
177
C     Read the name of the output trajectory file and set the constants file
178
      read(9,*) cdfname
179
 
180
C     Read the timestep for trajectory calculation (convert from minutes to hours)
181
      read(9,*) ts
182
      ts=ts/60.        
183
 
184
C     Read shift of start and end time relative to first data file
185
      read(9,*) tst
186
      read(9,*) ten
187
 
188
C     Read output time interval (in minutes)
189
      read(9,*) deltout
190
 
191
C     Read jumpflag (if =1 ground-touching trajectories reenter the atmosphere)
192
      read(9,*) jflag
193
 
194
C     Read factor for vertical velocity field
195
      read(9,*) wfactor
196
 
197
c     Read the reference time and the time range
198
      read(9,*) reftime(1)                  ! year
199
      read(9,*) reftime(2)                  ! month
200
      read(9,*) reftime(3)                  ! day
201
      read(9,*) reftime(4)                  ! hour
202
      read(9,*) reftime(5)                  ! min
203
      read(9,*) reftime(6)                  ! time range (in min)
204
 
205
c     Read flag for 'no time check'
206
      read(9,*) timecheck
207
 
208
c     Read in fplane option and latitude
209
      read(9,*) fflag                       ! fplane flag (1/0)                         ! if fflag=1 -> fplane
210
      read(9,*) lat_fixed                   ! fixed latitude for f-plane
211
 
212
c     Read in vertical velocity for balloon sounding
213
      read(9,*) balloon                     ! will be added to vertical wind, in [m/s]
214
 
215
c     Close the input file
216
      close(9)
217
 
218
c     Calculate the number of output time steps
219
      ntim = abs(reftime(6)/deltout) + 1
220
 
221
c     Set the formats of the input and output files
222
      call mode_tra(inpmode,strname)
223
      call mode_tra(outmode,cdfname)
224
      if (outmode.eq.-1) outmode=1
225
 
226
c     Write some status information
227
      print*,'---- INPUT PARAMETERS -----------------------------------'
228
      print* 
229
      print*,'  Forward/Backward       : ',fbflag
230
      print*,'  #input files           : ',numdat
231
      print*,'  First/last input file  : ',trim(dat(1)),' ... ',
232
     >                                     trim(dat(numdat))
233
      print*,'  time increment         : ',timeinc
234
      print*,'  Output file            : ',trim(cdfname)
235
      print*,'  Time step (min)        : ',60.*ts
236
      write(*,'(a27,f7.2,f7.2)') '   Time shift (start,end) : ',tst,ten
237
      print*,'  Output time interval   : ',deltout
238
      print*,'  Jump flag              : ',jflag
239
      print*,'  Vertical wind (scale)  : ',wfactor
240
      print*,'  Trajectory pos file    : ',trim(strname)
241
      print*,'  # of trajectories      : ',ntra
242
      print*,'  # of output timesteps  : ',ntim
243
      if ( inpmode.eq.-1) then
244
         print*,'  Input format           : (lon,lat,p)-list'
245
      else
246
         print*,'  Input format           : ',inpmode
247
      endif
248
      print*,'  Output format          : ',outmode
249
      print*,'  Periodicity            : ',per
250
      print*,'  Time check             : ',trim(timecheck)
251
      if ( abs(balloon).gt.eps ) then
252
         print*
253
         print*,'  WARNING -> Lagranto in Sounding Mode'
254
         print*
255
         print*,'  Balloon velocity [m/s] : ',balloon
256
      endif
257
      print*
258
 
259
      print*,'---- FIXED NUMERICAL PARAMETERS -------------------------'
260
      print*
261
      print*,'  Numerical scheme       : ',imethod
262
      print*,'  Number of iterations   : ',numit
263
      print*,'  Filename prefix        : ',prefix
264
      print*,'  Missing data value     : ',mdv
265
      print*
266
 
267
c     Run a simple consistency check
268
      call do_timecheck(dat,numdat,timeinc,fbflag,reftime)
269
 
270
c     --------------------------------------------------------------------
271
c     Read grid parameters, checks and allocate memory
272
c     --------------------------------------------------------------------
273
 
274
c     Read the constant grid parameters (nx,ny,nz,xmin,xmax,ymin,ymax,pollon,pollat)
275
c     The negative <-fid> of the file identifier is used as a flag for parameter retrieval  
276
      filename = prefix//dat(1)
277
      varname  = 'U'
278
      nx       = 1
279
      ny       = 1
280
      nz       = 1
281
      tload   = -tst
282
      call input_open (fid,filename)
283
      call input_grid (-fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny,
284
     >                 tload,pollon,pollat,polgam,rd,rd,nz,rd,timecheck)
285
      call input_close(fid)
286
 
287
C     Check if the number of levels is too large
288
      if (nz.gt.nlevmax) goto 993
289
 
290
C     Set logical flag for periodic data set (hemispheric or not)
291
      hem = 0
292
      if (per.eq.0.) then
293
         delta=xmax-xmin-360.
294
         if (abs(delta+dx).lt.eps) then               ! Program aborts: arrays must be closed
295
            goto 992
296
        else if (abs(delta).lt.eps) then              ! Periodic and hemispheric
297
           hem=1
298
           per=360.
299
        endif
300
      else                                            ! Periodic and hemispheric
301
         hem=1
302
         per=xmax-xmin+dx
303
      endif
304
 
305
C     Allocate memory for some meteorological arrays
306
      allocate(uut0(nx*ny*nz),stat=stat)
307
      if (stat.ne.0) print*,'*** error allocating array uut0 ***'   ! Zonal wind
308
      allocate(uut1(nx*ny*nz),stat=stat)
309
      if (stat.ne.0) print*,'*** error allocating array uut1 ***'
310
      allocate(vvt0(nx*ny*nz),stat=stat)
311
      if (stat.ne.0) print*,'*** error allocating array vvt0 ***'   ! Meridional wind
312
      allocate(vvt1(nx*ny*nz),stat=stat)
313
      if (stat.ne.0) print*,'*** error allocating array vvt1 ***'
314
      allocate(wwt0(nx*ny*nz),stat=stat)
315
      if (stat.ne.0) print*,'*** error allocating array wwt0 ***'   ! Vertical wind
316
      allocate(wwt1(nx*ny*nz),stat=stat)
317
      if (stat.ne.0) print*,'*** error allocating array wwt1 ***'
318
 
319
c	  Allocate memory for model grid
320
      allocate(zuv(nx*ny*nz),stat=stat)
321
      if (stat.ne.0) print*,'*** error allocating array zu ***'     ! Z for U and V grid
322
      allocate(zw(nx*ny*nz),stat=stat)
323
      if (stat.ne.0) print*,'*** error allocating array zw ***'     ! Z for W grid
324
      allocate(zb(nx*ny),stat=stat)
325
      if (stat.ne.0) print*,'*** error allocating array zb ***'     ! surface height
326
      allocate(uut0(nx*ny*nz),stat=stat)
327
 
328
C     Get memory for trajectory arrays
329
      allocate(trainp(ntra,1,ncol),stat=stat)
330
      if (stat.ne.0) print*,'*** error allocating array trainp   ***' ! Input start coordinates
331
      allocate(traout(ntra,ntim,4),stat=stat)
332
      if (stat.ne.0) print*,'*** error allocating array traout   ***' ! Output trajectories
333
      allocate(xx0(ntra),stat=stat)
334
      if (stat.ne.0) print*,'*** error allocating array xx0      ***' ! X position (longitude)
335
      allocate(yy0(ntra),stat=stat)
336
      if (stat.ne.0) print*,'*** error allocating array yy0      ***' ! Y position (latitude)
337
      allocate(pp0(ntra),stat=stat)
338
      if (stat.ne.0) print*,'*** error allocating array pp0      ***' ! Pressure
339
      allocate(leftflag(ntra),stat=stat)
340
      if (stat.ne.0) print*,'*** error allocating array leftflag ***' ! Leaving-domain flag
341
 
342
c     Write some status information
343
      print*,'---- CONSTANT GRID PARAMETERS ---------------------------'
344
      print*
345
      print*,'  xmin,xmax     : ',xmin,xmax
346
      print*,'  ymin,ymax     : ',ymin,ymax
347
      print*,'  dx,dy         : ',dx,dy
348
      print*,'  pollon,pollat : ',pollon,pollat
349
      print*,'  polgam        : ',polgam
350
      print*,'  nx,ny,nz      : ',nx,ny,nz
351
      print*,'  per, hem      : ',per,hem
352
      print*
353
 
354
c     --------------------------------------------------------------------
355
c     Initialize the trajectory calculation
356
c     --------------------------------------------------------------------
357
 
358
c     Read start coordinates from file - Format (lon,lat,lev)
359
      if (inpmode.eq.-1) then
360
         open(fid,file=strname)
361
          do i=1,ntra
362
             read(fid,*) xx0(i),yy0(i),pp0(i)
363
          enddo
364
         close(fid)
365
 
366
c     Read start coordinates from trajectory file - check consistency of ref time
367
      else
368
         call ropen_tra(cdfid,strname,ntra,1,ncol,reftmp,vars,inpmode)
369
         call read_tra (cdfid,trainp,ntra,1,ncol,inpmode)
370
         do i=1,ntra
371
            time   = trainp(i,1,1)
372
            xx0(i) = trainp(i,1,2) 
373
            yy0(i) = trainp(i,1,3) 
374
            pp0(i) = trainp(i,1,4) 
375
         enddo
376
         call close_tra(cdfid,inpmode)
377
 
378
         if ( ( reftime(1).ne.reftmp(1) ).or.
379
     >        ( reftime(2).ne.reftmp(2) ).or.
380
     >        ( reftime(3).ne.reftmp(3) ).or.
381
     >        ( reftime(4).ne.reftmp(4) ).or.
382
     >        ( reftime(5).ne.reftmp(5) ) )
383
     >   then
384
            print*,' WARNING: Inconsistent reference times'
385
            write(*,'(5f8.2)') (reftime(i),i=1,5)
386
            write(*,'(5f8.2)') (reftmp (i),i=1,5)
387
            print*,'Enter a key to proceed...'
388
            stop
389
         endif
390
      endif
391
 
392
c     Rotate coordinates
393
      do i=1,ntra
394
 
395
         lon = xx0(i)
396
         lat = yy0(i)
397
 
398
         if ( abs(polgam).gt.eps ) then
399
 
400
           rlat = phi2phirot ( lat, lon, pollat, pollon )
401
           rlon = rla2rlarot ( lat, lon, pollat, pollon, polgam )
402
 
403
         else
404
 
405
            if ( abs(pollat-90.).gt.eps) then
406
               rlon = lmtolms(lat,lon,pollat,pollon)
407
               rlat = phtophs(lat,lon,pollat,pollon)
408
            else
409
               rlon = lon
410
               rlat = lat
411
            endif
412
 
413
         endif
414
 
415
         xx0(i) = rlon
416
         yy0(i) = rlat
417
 
418
      enddo
419
 
420
c     Set sign of time range
421
      reftime(6) = fbflag * reftime(6)
422
 
423
c     Write some status information
424
      print*,'---- REFERENCE DATE---------- ---------------------------'
425
      print*
426
      print*,' Reference time (year)  :',reftime(1)
427
      print*,'                (month) :',reftime(2)
428
      print*,'                (day)   :',reftime(3)
429
      print*,'                (hour)  :',reftime(4)
430
      print*,'                (min)   :',reftime(5)
431
      print*,' Time range             :',reftime(6),' min'
432
      print*
433
 
434
C     Save starting positions 
435
      itim = 1
436
      do i=1,ntra
437
         traout(i,itim,1) = 0.
438
         traout(i,itim,2) = xx0(i)
439
         traout(i,itim,3) = yy0(i)
440
         traout(i,itim,4) = pp0(i)
441
      enddo
442
 
443
c     Init the flag and the counter for trajectories leaving the domain
444
      leftcount=0
445
      do i=1,ntra
446
         leftflag(i)=0
447
      enddo
448
 
449
C     Convert time shifts <tst,ten> from <hh.mm> into fractional time
450
      call hhmm2frac(tst,frac)
451
      tst = frac
452
      call hhmm2frac(ten,frac)
453
      ten = frac
454
 
455
c     Read models heights for U, V, W grid
456
	  filename = prefix//dat(1)
457
	  call input_open (fid,filename)
458
	  varname = 'U'
459
	  stagz   =  0.
460
	  call input_grid
461
     >    (fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny,
462
     >     tload,pollon,pollat,polgam,zuv,zb,nz,stagz,timecheck)
463
	  varname = 'W'
464
	  stagz   = -0.5
465
	  call input_grid
466
     >    (fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny,
467
     >     tload,pollon,pollat,polgam,zw,zb,nz,stagz,timecheck)
468
 
469
      call input_close(fid)
470
 
471
c     -----------------------------------------------------------------------
472
c     Loop to calculate trajectories
473
c     -----------------------------------------------------------------------
474
 
475
c     Write some status information
476
      print*,'---- TRAJECTORIES ----------- ---------------------------'
477
      print*    
478
 
479
C     Set the time for the first data file (depending on forward/backward mode)
480
      if (fbflag.eq.1) then
481
        tstart = -tst
482
      else
483
        tstart = tst
484
      endif
485
 
486
c     Read wind fields and vertical grid from first file
487
      filename = prefix//dat(1)
488
 
489
      call frac2hhmm(tstart,tload)
490
 
491
      write(*,'(a16,a20,f7.2)') '  (file,time) : ',
492
     >                       trim(filename),tload
493
 
494
      call input_open (fid,filename)
495
      varname='U'                                      ! U
496
      call input_wind 
497
     >    (fid,varname,uut1,tload,stagz,mdv,
498
     >     xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
499
      varname='V'                                      ! V
500
      call input_wind 
501
     >    (fid,varname,vvt1,tload,stagz,mdv,
502
     >     xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
503
      varname='W'                                      ! W
504
      call input_wind 
505
     >    (fid,varname,wwt1,tload,stagz,mdv,
506
     >     xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)      
507
      call input_close(fid)
508
 
509
c     Loop over all input files (time step is <timeinc>)
510
      do itm=1,numdat-1
511
 
512
c       Calculate actual and next time
513
        time0 = tstart+real(itm-1)*timeinc*fbflag
514
        time1 = time0+timeinc*fbflag
515
 
516
c       Copy old velocities to new ones
517
        do i=1,nx*ny*nz
518
           uut0(i)=uut1(i)
519
           vvt0(i)=vvt1(i)
520
           wwt0(i)=wwt1(i)
521
        enddo
522
 
523
c       Read wind fields and surface pressure at next time
524
        filename = prefix//dat(itm+1)
525
 
526
        call frac2hhmm(time1,tload)
527
        write(*,'(a16,a20,f7.2)') '  (file,time) : ',
528
     >                          trim(filename),tload
529
 
530
        call input_open (fid,filename)
531
        varname='U'                                     ! U
532
        call input_wind 
533
     >       (fid,varname,uut1,tload,stagz,mdv,
534
     >        xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
535
        varname='V'                                     ! V
536
        call input_wind 
537
     >       (fid,varname,vvt1,tload,stagz,mdv,
538
     >        xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
539
        varname='W'                                     ! W
540
        call input_wind 
541
     >       (fid,varname,wwt1,tload,stagz,mdv,
542
     >        xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)      
543
        call input_close(fid)
544
 
545
C       Determine the first and last loop indices
546
        if (numdat.eq.2) then
547
          filo = nint(tst/ts)+1
548
          lalo = nint((timeinc-ten)/ts) 
549
        elseif ( itm.eq.1 ) then
550
          filo = nint(tst/ts)+1
551
          lalo = nint(timeinc/ts)
552
        else if (itm.eq.numdat-1) then
553
          filo = 1
554
          lalo = nint((timeinc-ten)/ts)
555
        else
556
          filo = 1
557
          lalo = nint(timeinc/ts)
558
        endif
559
 
560
c       Split the interval <timeinc> into computational time steps <ts>
561
        do iloop=filo,lalo
562
 
563
C         Calculate relative time position in the interval timeinc (0=beginning, 1=end)
564
          reltpos0 = ((real(iloop)-1.)*ts)/timeinc
565
          reltpos1 = real(iloop)*ts/timeinc
566
 
567
C         Initialize counter for domain leaving trajectories
568
          leftcount=0
569
 
570
c         Timestep for all trajectories
571
          do i=1,ntra
572
 
573
C           Check if trajectory has already left the data domain
574
            if (leftflag(i).ne.1) then	
575
 
576
c             Iterative Euler timestep (x0,y0,p0 -> x1,y1,p1)
577
              if (imethod.eq.1) then
578
                 call euler(
579
     >               xx1,yy1,pp1,leftflag(i),
580
     >               xx0(i),yy0(i),pp0(i),reltpos0,reltpos1,
581
     >               ts*3600,numit,jflag,mdv,wfactor,fbflag,
582
     >               zuv,zw,zb,uut0,uut1,vvt0,vvt1,wwt0,wwt1,
583
     >               xmin,ymin,dx,dy,per,hem,nx,ny,nz,fflag,
584
     >               lat_fixed,balloon)
585
 
586
c             Runge-Kutta timestep (x0,y0,p0 -> x1,y1,p1)
587
              else if (imethod.eq.2) then
588
                 call runge(
589
     >               xx1,yy1,pp1,leftflag(i),
590
     >               xx0(i),yy0(i),pp0(i),reltpos0,reltpos1,
591
     >               ts*3600,numit,jflag,mdv,wfactor,fbflag,
592
     >               zuv,zw,zb,uut0,uut1,vvt0,vvt1,wwt0,wwt1,
593
     >               xmin,ymin,dx,dy,per,hem,nx,ny,nz)
594
 
595
              endif
596
 
597
c             Update trajectory position, or increase number of trajectories leaving domain
598
              if (leftflag(i).eq.1) then
599
                leftcount=leftcount+1
600
                print*,'     -> Trajectory ',i,' leaves domain'
601
                xx0(i)=-999.
602
                yy0(i)=-999.
603
                pp0(i)=-999.
604
              else
605
                xx0(i)=xx1
606
                yy0(i)=yy1
607
                pp0(i)=pp1
608
              endif
609
 
610
c          Trajectory has already left data domain (mark as <mdv>)
611
           else
612
              xx0(i)=-999.
613
              yy0(i)=-999.
614
              pp0(i)=-999.
615
 
616
           endif
617
 
618
          enddo
619
 
620
C         Save positions only every deltout minutes
621
          delta = aint(iloop*60*ts/deltout)-iloop*60*ts/deltout
622
          if (abs(delta).lt.eps) then
623
            time = time0+reltpos1*timeinc*fbflag
624
            itim = itim + 1
625
            do i=1,ntra
626
               call frac2hhmm(time,tload)
627
               traout(i,itim,1) = tload
628
               if ( abs( pp0(i) + 999. ).lt.eps ) then
629
                   traout(i,itim,2) = traout(i,itim-1,2)
630
                   traout(i,itim,3) = traout(i,itim-1,3)
631
                   traout(i,itim,4) = -999.
632
               else
633
                   traout(i,itim,2) = xx0(i)
634
                   traout(i,itim,3) = yy0(i)
635
                   traout(i,itim,4) = pp0(i)
636
               endif
637
            enddo
638
          endif
639
 
640
        enddo
641
 
642
      enddo
643
 
644
c     Coordinate rotation
645
      do i=1,ntra
646
         do j=1,ntim
647
 
648
            rlon = traout(i,j,2)
649
            rlat = traout(i,j,3)
650
 
651
            if ( abs(polgam).gt.eps ) then
652
 
653
              lon = rlarot2rla (rlat, rlon, pollat, pollon, polgam)
654
              lat = phirot2phi (rlat, rlon, pollat, pollon, polgam)
655
 
656
            else
657
 
658
              if ( abs(pollat-90.).gt.eps) then
659
                 lon = lmstolm(rlat,rlon,pollat,pollon)
660
                 lat = phstoph(rlat,rlon,pollat,pollon)
661
              else
662
                 lon = rlon
663
                 lat = rlat
664
              endif
665
 
666
            endif
667
 
668
            traout(i,j,2) = lon
669
            traout(i,j,3) = lat
670
 
671
         enddo
672
      enddo
673
 
674
c     Write trajectory file
675
      vars(1)  ='time'
676
      vars(2)  ='lon'
677
      vars(3)  ='lat'
678
      vars(4)  ='z'
679
      call wopen_tra(cdfid,cdfname,ntra,ntim,4,reftime,vars,outmode)
680
      call write_tra(cdfid,traout,ntra,ntim,4,outmode)
681
      call close_tra(cdfid,outmode)   
682
 
683
c     Write some status information, and end of program message
684
      print*  
685
      print*,'---- STATUS INFORMATION --------------------------------'
686
      print*
687
      print*,'  #leaving domain    ', leftcount
688
      print*,'  #staying in domain ', ntra-leftcount
689
      print*
690
      print*,'              *** END OF PROGRAM CALTRA ***'
691
      print*,'========================================================='
692
 
693
      stop
694
 
695
c     ------------------------------------------------------------------
696
c     Exception handling
697
c     ------------------------------------------------------------------
698
 
699
 991  write(*,*) '*** ERROR: all start points outside the data domain'
700
      call exit(1)
701
 
702
 992  write(*,*) '*** ERROR: close arrays on files (prog. closear)'
703
      call exit(1)
704
 
705
 993  write(*,*) '*** ERROR: problems with array size'
706
      call exit(1)
707
 
708
      end 
709
 
710
 
711
c     *******************************************************************
712
c     * Time step : either Euler or Runge-Kutta                         *
713
c     *******************************************************************
714
 
715
C     Time-step from (x0,y0,p0) to (x1,y1,p1)
716
C
717
C     (x0,y0,p0) input	coordinates (long,lat,p) for starting point
718
C     (x1,y1,p1) output	coordinates (long,lat,p) for end point
719
C     deltat	 input	timestep in seconds
720
C     numit	 input	number of iterations
721
C     jump	 input  flag (=1 trajectories don't enter the ground)
722
C     left	 output	flag (=1 if trajectory leaves data domain)
723
 
724
c     -------------------------------------------------------------------
725
c     Iterative Euler time step
726
c     -------------------------------------------------------------------
727
 
728
      subroutine euler(x1,y1,p1,left,x0,y0,p0,reltpos0,reltpos1,
729
     >                 deltat,numit,jump,mdv,wfactor,fbflag,
730
     >		           zuv,zw,zb,uut0,uut1,vvt0,vvt1,wwt0,wwt1,
731
     >                 xmin,ymin,dx,dy,per,hem,nx,ny,nz,
732
     >                 fflag,lat_fixed,balloon)
733
 
734
      implicit none
735
 
736
c     Declaration of subroutine parameters
737
      integer      nx,ny,nz
738
      real         x1,y1,p1
739
      integer      left
740
      real	       x0,y0,p0
741
      real         reltpos0,reltpos1
742
      real   	   deltat
743
      integer      numit
744
      integer      jump
745
      real         wfactor
746
      integer      fbflag
747
      real     	   zb(nx*ny)
748
      real         uut0(nx*ny*nz),uut1(nx*ny*nz)
749
      real 	       vvt0(nx*ny*nz),vvt1(nx*ny*nz)
750
      real         wwt0(nx*ny*nz),wwt1(nx*ny*nz)
751
      real         zuv(nx*ny*nz),zw(nx*ny*nz)
752
      real         xmin,ymin,dx,dy
753
      real         per
754
      integer      hem
755
      real         mdv
756
      integer      fflag
757
      real         lat_fixed
758
      real         balloon
759
 
760
c     Numerical and physical constants
761
      real         deltay
762
      parameter    (deltay=1.112E5)  ! Distance in m between 2 lat circles
763
      real         pi                       
764
      parameter    (pi=3.1415927)    ! Pi
765
 
766
c     Test flag
767
      integer      test
768
      parameter    (test=0)
769
      integer      maxindstep
770
      parameter    (maxindstep=5.)
771
 
772
c     Auxiliary variables
773
      real         xmax,ymax
774
      real	       xind,yind,pind
775
      real	       u0,v0,w0,u1,v1,w1,u,v,w,sp
776
      integer	   icount
777
      character    ch
778
      real         xind0,yind0,zind0,xind1,yind1,zind1
779
 
780
c     Externals    
781
      real         int_index4
782
      external     int_index4
783
 
784
c     Reset the flag for domain-leaving
785
      left=0
786
 
787
c     Set the esat-north bounray of the domain
788
      xmax = xmin+real(nx-1)*dx
789
      ymax = ymin+real(ny-1)*dy
790
 
791
C     Interpolate wind fields to starting position (x0,y0,p0)
792
      call get_index4 (xind,yind,pind,x0,y0,p0,reltpos0,
793
     >                 zuv,zuv,zb,zb,3,nx,ny,nz,xmin,ymin,dx,dy,mdv)
794
      u0 = int_index4(uut0,uut1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
795
      v0 = int_index4(vvt0,vvt1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
796
 
797
      call get_index4 (xind,yind,pind,x0,y0,p0,reltpos0,
798
     >                 zw,zw,zb,zb,3,nx,ny,nz,xmin,ymin,dx,dy,mdv)
799
      w0 = int_index4(wwt0,wwt1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
800
 
801
c     In testmode: remember the indices at starting position
802
      if ( test.eq.1 ) then
803
           zind0 = pind
804
           xind0 = xind
805
           yind0 = yind
806
      endif
807
 
808
c     Force the near-surface wind to zero
809
      if (pind.lt.1.) w0=w0*pind
810
 
811
C     For first iteration take ending position equal to starting position
812
      x1=x0
813
      y1=y0
814
      p1=p0
815
 
816
C     Iterative calculation of new position
817
      do icount=1,numit
818
 
819
C        Calculate new winds for advection
820
         call get_index4 (xind,yind,pind,x1,y1,p1,reltpos1,
821
     >                    zuv,zuv,zb,zb,3,nx,ny,nz,xmin,ymin,dx,dy,mdv)
822
         u1 = int_index4(uut0,uut1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
823
         v1 = int_index4(vvt0,vvt1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
824
 
825
         call get_index4 (xind,yind,pind,x1,y1,p1,reltpos1,
826
     >                    zw,zw,zb,zb,3,nx,ny,nz,xmin,ymin,dx,dy,mdv)
827
         w1 = int_index4(wwt0,wwt1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
828
 
829
c        Force the near-surface wind to zero
830
         if (pind.lt.1.) w1=w1*pind
831
 
832
c        Get the new velocity in between
833
         u=(u0+u1)/2.
834
         v=(v0+v1)/2.
835
         w=(w0+w1)/2. + balloon
836
 
837
C        Calculate new positions (fplane or normal)
838
         if ( fflag.eq.1 ) then
839
            x1 = x0 + fbflag*u*deltat/(deltay*cos(lat_fixed*pi/180.))
840
            y1 = y0 + fbflag*v*deltat/deltay
841
            p1 = p0 + fbflag*wfactor*w*deltat
842
         else
843
            x1 = x0 + fbflag*u*deltat/(deltay*cos(y0*pi/180.))
844
            y1 = y0 + fbflag*v*deltat/deltay
845
            p1 = p0 + fbflag*wfactor*w*deltat
846
         endif
847
 
848
c       Handle pole problems (crossing and near pole trajectory)
849
        if ((hem.eq.1).and.(y1.gt.90.)) then
850
          y1=180.-y1
851
          x1=x1+per/2.
852
        endif
853
        if ((hem.eq.1).and.(y1.lt.-90.)) then
854
          y1=-180.-y1
855
          x1=x1+per/2.
856
        endif
857
        if (y1.gt.89.99) then
858
           y1=89.99
859
        endif
860
 
861
c       Handle crossings of the dateline
862
        if ((hem.eq.1).and.(x1.gt.xmin+per-dx)) then
863
           x1=xmin+amod(x1-xmin+dx,per)
864
        endif
865
        if ((hem.eq.1).and.(x1.lt.xmin)) then
866
           x1=xmin+per+amod(x1-xmin-dx,per)
867
        endif
868
 
869
C       Interpolate surface height to actual position
870
        call get_index4 (xind,yind,pind,x1,y1,0.,reltpos1,
871
     >                   zuv,zuv,zb,zb,3,
872
     >                   nx,ny,nz,xmin,ymin,dx,dy,mdv)
873
        sp = int_index4 (zb,zb,nx,ny,1,xind,yind,1.,reltpos1,mdv)
874
 
875
c       Handle trajectories which cross the lower boundary (jump flag)
876
        if ((jump.eq.1).and.(p1.lt.sp)) p1=sp+10.
877
 
878
C       Check if trajectory leaves data domain
879
        if ( ( (hem.eq.0).and.(x1.lt.xmin)    ).or.
880
     >       ( (hem.eq.0).and.(x1.gt.xmax-dx) ).or.
881
     >         (y1.lt.ymin).or.(y1.gt.ymax).or.(p1.lt.sp) )
882
     >  then
883
          left=1
884
          goto 100
885
        endif
886
 
887
      enddo
888
 
889
c     In test mode, get the indices at the final position
890
      if ( test.eq.1 ) then
891
         call get_index4 (xind,yind,pind,x1,y1,p1,reltpos1,
892
     >                    zuv,zuv,zb,zb,3,nx,ny,nz,xmin,ymin,dx,dy,mdv)
893
         xind1 = xind
894
         yind1 = yind
895
         call get_index4 (xind,yind,pind,x1,y1,p1,reltpos1,
896
     >                    zw,zw,zb,zb,3,nx,ny,nz,xmin,ymin,dx,dy,mdv)
897
         zind1 = pind
898
 
899
         if ( ( abs(xind1-xind0).gt.maxindstep ).or.
900
     >        ( abs(yind1-yind0).gt.maxindstep ).or.
901
     >        ( abs(zind1-zind0).gt.maxindstep ) )
902
     >   then
903
             print*
904
             print*,' ERROR: Index step too large : convergence problem'
905
             print*
906
             print*,'xind : ',xind0,' -> ',xind1,' : ',abs(xind1-xind0)
907
             print*,'yind : ',yind0,' -> ',yind1,' : ',abs(yind1-yind0)
908
             print*,'zind : ',zind0,' -> ',zind1,' : ',abs(zind1-zind0)
909
             print*
910
         endif
911
 
912
      endif
913
 
914
c     Exit point for subroutine
915
 100  continue
916
 
917
      return
918
 
919
      end
920
 
921
c     -------------------------------------------------------------------
922
c     Runge-Kutta (4th order) time-step
923
c     -------------------------------------------------------------------
924
 
925
      subroutine runge(x1,y1,p1,left,x0,y0,p0,reltpos0,reltpos1,
926
     >                 deltat,numit,jump,mdv,wfactor,fbflag,
927
     >		           zuv,zw,zb,uut0,uut1,vvt0,vvt1,wwt0,wwt1,
928
     >                 xmin,ymin,dx,dy,per,hem,nx,ny,nz)
929
 
930
      implicit none
931
 
932
c     Declaration of subroutine parameters
933
      integer      nx,ny,nz
934
      real         x1,y1,p1
935
      integer      left
936
      real	       x0,y0,p0
937
      real         reltpos0,reltpos1
938
      real   	   deltat
939
      integer      numit
940
      integer      jump
941
      real         wfactor
942
      integer      fbflag
943
      real     	   zb(nx*ny)
944
      real         uut0(nx*ny*nz),uut1(nx*ny*nz)
945
      real 	       vvt0(nx*ny*nz),vvt1(nx*ny*nz)
946
      real         wwt0(nx*ny*nz),wwt1(nx*ny*nz)
947
      real         zuv(nx*ny*nz),zw(nx*ny*nz)
948
      real         xmin,ymin,dx,dy
949
      real         per
950
      integer      hem
951
      real         mdv
952
 
953
c     Numerical and physical constants
954
      real         deltay
955
      parameter    (deltay=1.112E5)  ! Distance in m between 2 lat circles
956
      real         pi                       
957
      parameter    (pi=3.1415927)    ! Pi
958
 
959
c     Auxiliary variables
960
      real         xmax,ymax
961
      real	       xind,yind,pind
962
      real	       u0,v0,w0,u1,v1,w1,u,v,w,sp
963
      integer	   icount,n
964
      real         xs,ys,ps,xk(4),yk(4),pk(4)
965
      real         reltpos
966
 
967
c     Externals    
968
      real         int_index4
969
      external     int_index4
970
 
971
c     Reset the flag for domain-leaving
972
      left=0
973
 
974
c     Set the esat-north bounray of the domain
975
      xmax = xmin+real(nx-1)*dx
976
      ymax = ymin+real(ny-1)*dy
977
 
978
c     Apply the Runge Kutta scheme
979
      do n=1,4
980
 
981
c       Get intermediate position and relative time
982
        if (n.eq.1) then
983
          xs=0.
984
          ys=0.
985
          ps=0.
986
          reltpos=reltpos0
987
        else if (n.eq.4) then
988
          xs=xk(3)
989
          ys=yk(3)
990
          ps=pk(3)
991
          reltpos=reltpos1
992
        else
993
          xs=xk(n-1)/2.
994
          ys=yk(n-1)/2.
995
          ps=pk(n-1)/2.
996
          reltpos=(reltpos0+reltpos1)/2.
997
        endif
998
 
999
C       Calculate new winds for advection
1000
        call get_index4 (xind,yind,pind,x0+xs,y0+ys,p0+ps,reltpos,
1001
     >                   zuv,zuv,zb,zb,3,nx,ny,nz,xmin,ymin,dx,dy,mdv)
1002
        u = int_index4(uut0,uut1,nx,ny,nz,xind,yind,pind,reltpos,mdv)
1003
        v = int_index4(vvt0,vvt1,nx,ny,nz,xind,yind,pind,reltpos,mdv)
1004
 
1005
        call get_index4 (xind,yind,pind,x0+xs,y0+ys,p0+ps,reltpos,
1006
     >                   zw,zw,zb,zb,3,nx,ny,nz,xmin,ymin,dx,dy,mdv)
1007
        w = int_index4(wwt0,wwt1,nx,ny,nz,xind,yind,pind,reltpos,mdv)
1008
 
1009
c       Force the near-surface wind to zero
1010
        if (pind.lt.1.) w1=w1*pind
1011
 
1012
c       Update position and keep them
1013
        xk(n)=fbflag*u*deltat/(deltay*cos(y0*pi/180.))
1014
        yk(n)=fbflag*v*deltat/deltay
1015
        pk(n)=fbflag*w*deltat*wfactor
1016
 
1017
      enddo
1018
 
1019
C     Calculate new positions
1020
      x1=x0+(1./6.)*(xk(1)+2.*xk(2)+2.*xk(3)+xk(4))
1021
      y1=y0+(1./6.)*(yk(1)+2.*yk(2)+2.*yk(3)+yk(4))
1022
      p1=p0+(1./6.)*(pk(1)+2.*pk(2)+2.*pk(3)+pk(4))
1023
 
1024
c     Handle pole problems (crossing and near pole trajectory)
1025
      if ((hem.eq.1).and.(y1.gt.90.)) then
1026
         y1=180.-y1
1027
         x1=x1+per/2.
1028
      endif
1029
      if ((hem.eq.1).and.(y1.lt.-90.)) then
1030
         y1=-180.-y1
1031
         x1=x1+per/2.
1032
      endif
1033
      if (y1.gt.89.99) then
1034
         y1=89.99
1035
      endif
1036
 
1037
c     Handle crossings of the dateline
1038
      if ((hem.eq.1).and.(x1.gt.xmin+per-dx)) then
1039
         x1=xmin+amod(x1-xmin,per)
1040
      endif
1041
      if ((hem.eq.1).and.(x1.lt.xmin)) then
1042
         x1=xmin+per+amod(x1-xmin,per)
1043
      endif
1044
 
1045
C     Interpolate surface pressure to actual position
1046
      call get_index4 (xind,yind,pind,x1,y1,0.,reltpos1,
1047
     >                 zuv,zuv,zb,zb,3,
1048
     >                 nx,ny,nz,xmin,ymin,dx,dy,mdv)
1049
      sp = int_index4 (zb,zb,nx,ny,1,xind,yind,1,reltpos,mdv)
1050
 
1051
c     Handle trajectories which cross the lower boundary (jump flag)
1052
      if ((jump.eq.1).and.(p1.lt.sp)) p1=sp+10.
1053
 
1054
C     Check if trajectory leaves data domain
1055
      if ( ( (hem.eq.0).and.(x1.lt.xmin)    ).or.
1056
     >     ( (hem.eq.0).and.(x1.gt.xmax-dx) ).or.
1057
     >       (y1.lt.ymin).or.(y1.gt.ymax).or.(p1.lt.sp) )
1058
     >then
1059
         left=1
1060
         goto 100
1061
      endif
1062
 
1063
c     Exit point fdor subroutine
1064
 100  continue
1065
 
1066
      return
1067
      end