Subversion Repositories lagranto.um

Rev

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