Subversion Repositories lagranto.20cr

Rev

Details | Last modification | View Log | RSS feed

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