Subversion Repositories lagranto.ecmwf

Rev

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