Subversion Repositories lagranto.um

Rev

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