Subversion Repositories lagranto.arpege

Rev

Details | Last modification | View Log | RSS feed

Rev Author Line No. Line
3 michaesp 1
      PROGRAM caltra
2
 
3
C     ********************************************************************
4
C     *                                                                  *
5
C     * Calculates trajectories                                          *
6
C     *                                                                  *
7
C     *	Heini Wernli	   first version:	April 1993               *
8
C     * Michael Sprenger   major upgrade:       2008-2009                *
9
C     *                                                                  *
10
C     ********************************************************************
11
 
12
      implicit none      
13
 
14
c     --------------------------------------------------------------------
15
c     Declaration of parameters
16
c     --------------------------------------------------------------------
17
 
18
c     Maximum number of levels for input files
19
      integer   nlevmax
20
      parameter	(nlevmax=200)
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 for primary and secondary files (typically 'P' and 'S')
47
      character*1 prefix
48
      parameter   (prefix='P')
49
      character*1 srefix
50
      parameter   (srefix='S')
51
 
52
c     Physical constants - needed to compute potential temperature
53
      real      rdcp,tzero
54
      data      rdcp,tzero /0.286,273.15/
55
 
56
c     --------------------------------------------------------------------
57
c     Declaration of variables
58
c     --------------------------------------------------------------------
59
 
60
c     Input parameters
61
      integer                                fbflag          ! Flag for forward/backward mode
62
      integer                                numdat          ! Number of input files
63
      character*11                           dat(ndatmax)    ! Dates of input files
64
      real                                   timeinc         ! Time increment between input files
65
      real                                   per             ! Periodicity (=0 if none)
66
      integer                                ntra            ! Number of trajectories
67
      character*80                           cdfname         ! Name of output files
68
      real                                   ts              ! Time step
69
      real                                   tst,ten         ! Shift of start and end time relative to first data file
70
      integer                                deltout         ! Output time interval (in minutes)
71
      integer                                jflag           ! Jump flag (if =1 ground-touching trajectories reenter atmosphere)
72
      real                                   wfactor         ! Factor for vertical velocity field
73
      character*80                           strname         ! File with start positions
74
      character*80                           timecheck       ! Either 'yes' or 'no'
75
      character*80                           isen            ! Isentropic trajectories ('yes' or 'no')
76
      integer                                thons           ! Isentropic mode: is TH availanle on S
77
      character*80                           modlev          ! 2D (model level) trajectories ('yes' or 'no')
78
 
79
c     Trajectories
80
      integer                                ncol            ! Number of columns for insput trajectories
81
      real,allocatable, dimension (:,:,:) :: trainp          ! Input start coordinates (ntra,1,ncol)
82
      real,allocatable, dimension (:,:,:) :: traout          ! Output trajectories (ntra,ntim,4)
83
      integer                                reftime(6)      ! Reference date
84
      character*80                           vars(200)       ! Field names
85
      real,allocatable, dimension (:)     :: xx0,yy0,pp0     ! Position of air parcels
86
      integer,allocatable, dimension (:)  :: leftflag        ! Flag for domain-leaving
87
      real                                   xx1,yy1,pp1     ! Updated position of air parcel
88
      integer                                leftcount       ! Number of domain leaving trajectories
89
      integer                                ntim            ! Number of output time steps
90
      real,allocatable, dimension (:)     :: theta           ! Potential temperature for isentropic trajectories
91
      real,allocatable, dimension (:)     :: zindex          ! Vertical index for model-level (2D) trajectories
92
 
93
c     Meteorological fields
94
      real,allocatable, dimension (:)     :: spt0,spt1       ! Surface pressure
95
      real,allocatable, dimension (:)     :: uut0,uut1       ! Zonal wind
96
      real,allocatable, dimension (:)     :: vvt0,vvt1       ! Meridional wind
97
      real,allocatable, dimension (:)     :: wwt0,wwt1       ! Vertical wind
98
      real,allocatable, dimension (:)     :: p3t0,p3t1       ! 3d-pressure 
99
      real,allocatable, dimension (:)     :: tht0,tht1       ! 3d potential temperature
100
      real,allocatable, dimension (:)     :: sth0,sth1       ! Surface potential temperature
101
 
102
c     Grid description
103
      real                                   pollon,pollat   ! Longitude/latitude of pole
104
      real                                   ak(nlevmax)     ! Vertical layers and levels
105
      real                                   bk(nlevmax) 
106
      real                                   xmin,xmax       ! Zonal grid extension
107
      real                                   ymin,ymax       ! Meridional grid extension
108
      integer                                nx,ny,nz        ! Grid dimensions
109
      real                                   dx,dy           ! Horizontal grid resolution
110
      integer                                hem             ! Flag for hemispheric domain
111
      real                                   mdv             ! Missing data value
112
 
113
c     Auxiliary variables                 
114
      real                                   delta,rd
115
      integer	                             itm,iloop,i,j,k,filo,lalo
116
      integer                                ierr,stat
117
      integer                                cdfid,fid
118
      real	                                 tstart,time0,time1,time
119
      real                                   reltpos0,reltpos1
120
      real                                   xind,yind,pind,pp,sp,stagz
121
      character*80                           filename,varname
122
      integer                                reftmp(6)
123
      character                              ch
124
      real                                   frac,tload
125
      integer                                itim
126
      integer                                wstep
127
      real                                   x1,y1,p1
128
      real                                   thetamin,thetamax
129
      real                                   zindexmin,zindexmax
130
 
131
c     Externals
132
      real                                   int_index4
133
      external                               int_index4
134
 
135
c     --------------------------------------------------------------------
136
c     Start of program, Read parameters
137
c     --------------------------------------------------------------------
138
 
139
c     Write start message
140
      print*,'========================================================='
141
      print*,'              *** START OF PROGRAM CALTRA ***'
142
      print*
143
 
144
c     Open the parameter file
145
      open(9,file='caltra.param')
146
 
147
c     Read flag for forward/backward mode (fbflag)
148
      read(9,*) fbflag
149
 
150
c     Read number of input files (numdat)
151
      read(9,*) numdat
152
      if (numdat.gt.ndatmax) then
153
        print*,' ERROR: too many input files ',numdat,ndatmax
154
        goto 993
155
      endif
156
 
157
c     Read list of input dates (dat, sort depending on forward/backward mode)
158
      if (fbflag.eq.1) then
159
        do itm=1,numdat
160
          read(9,'(a11)') dat(itm)
161
        enddo
162
      else
163
        do itm=numdat,1,-1
164
          read(9,'(a11)') dat(itm)
165
        enddo
166
      endif
167
 
168
c     Read time increment between input files (timeinc)
169
      read(9,*) timeinc
170
 
171
C     Read if data domain is periodic and its periodicity
172
      read(9,*) per
173
 
174
c     Read the number of trajectories and name of position file
175
      read(9,*) strname
176
      read(9,*) ntra
177
      read(9,*) ncol 
178
      if (ntra.eq.0) goto 991
179
 
180
C     Read the name of the output trajectory file and set the constants file
181
      read(9,*) cdfname
182
 
183
C     Read the timestep for trajectory calculation (convert from minutes to hours)
184
      read(9,*) ts
185
      ts=ts/60.        
186
 
187
C     Read shift of start and end time relative to first data file
188
      read(9,*) tst
189
      read(9,*) ten
190
 
191
C     Read output time interval (in minutes)
192
      read(9,*) deltout
193
 
194
C     Read jumpflag (if =1 ground-touching trajectories reenter the atmosphere)
195
      read(9,*) jflag
196
 
197
C     Read factor for vertical velocity field
198
      read(9,*) wfactor
199
 
200
c     Read the reference time and the time range
201
      read(9,*) reftime(1)                  ! year
202
      read(9,*) reftime(2)                  ! month
203
      read(9,*) reftime(3)                  ! day
204
      read(9,*) reftime(4)                  ! hour
205
      read(9,*) reftime(5)                  ! min
206
      read(9,*) reftime(6)                  ! time range (in min)
207
 
208
c     Read flag for 'no time check'
209
      read(9,*) timecheck
210
 
211
c     Read flag for isentropic trajectories
212
      read(9,*) isen, thons
213
 
214
c     Read flag for model-level trajectories (2D mode)
215
      read(9,*) modlev
216
 
217
c     Close the input file
218
      close(9)
219
 
220
c     Calculate the number of output time steps
221
      ntim = abs(reftime(6)/deltout) + 1
222
 
223
c     Set the formats of the input and output files
224
      call mode_tra(inpmode,strname)
225
      call mode_tra(outmode,cdfname)
226
      if (outmode.eq.-1) outmode=1
227
 
228
c     Write some status information
229
      print*,'---- INPUT PARAMETERS -----------------------------------'
230
      print* 
231
      print*,'  Forward/Backward       : ',fbflag
232
      print*,'  #input files           : ',numdat
233
      print*,'  First/last input file  : ',trim(dat(1)),' ... ',
234
     >                                     trim(dat(numdat))
235
      print*,'  time increment         : ',timeinc
236
      print*,'  Output file            : ',trim(cdfname)
237
      print*,'  Time step (min)        : ',60.*ts
238
      write(*,'(a27,f7.2,f7.2)') '   Time shift (start,end) : ',tst,ten
239
      print*,'  Output time interval   : ',deltout
240
      print*,'  Jump flag              : ',jflag
241
      print*,'  Vertical wind (scale)  : ',wfactor
242
      print*,'  Trajectory pos file    : ',trim(strname)
243
      print*,'  # of trajectories      : ',ntra
244
      print*,'  # of output timesteps  : ',ntim
245
      if ( inpmode.eq.-1) then
246
         print*,'  Input format           : (lon,lat,p)-list'
247
      else
248
         print*,'  Input format           : ',inpmode
249
      endif
250
      print*,'  Output format          : ',outmode
251
      print*,'  Periodicity            : ',per
252
      print*,'  Time check             : ',trim(timecheck)
253
      print*,'  Isentropic trajectories: ',trim(isen),thons
254
      print*,'  Model-level trajs (2D) : ',trim(modlev)
255
      print*
256
 
257
      if ( (isen.eq.'yes').and.(modlev.eq.'yes') ) then
258
         print*,
259
     >    '  WARNING: isentropic and 2D mode chosen -> 2D accepted'
260
         print*
261
         isen = 'no'
262
      endif
263
 
264
c     Init missing data value
265
      mdv = -999.
266
 
267
      print*,'---- FIXED NUMERICAL PARAMETERS -------------------------'
268
      print*
269
      print*,'  Numerical scheme       : ',imethod
270
      print*,'  Number of iterations   : ',numit
271
      print*,'  Filename prefix        : ',prefix
272
      print*,'  Missing data value     : ',mdv
273
      print*
274
 
275
c     --------------------------------------------------------------------
276
c     Read grid parameters, checks and allocate memory
277
c     --------------------------------------------------------------------
278
 
279
c     Read the constant grid parameters (nx,ny,nz,xmin,xmax,ymin,ymax,pollon,pollat)
280
c     The negative <-fid> of the file identifier is used as a flag for parameter retrieval  
281
      filename = prefix//dat(1)
282
      varname  = 'U'
283
      nx       = 1
284
      ny       = 1
285
      nz       = 1
286
      tload    = -tst
287
 
288
      call input_open (fid,filename)
289
c     read grid information on input file, subroutine in ${LAGRANTO}/lib/ioip_mf.f
290
      call input_grid (-fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny,
291
     >                 tload,pollon,pollat,rd,rd,nz,rd,rd,rd,timecheck)
292
 
293
      call input_close(fid)
294
 
295
C     Check if the number of levels is too large
296
      if (nz.gt.nlevmax) goto 993
297
 
298
C     Set logical flag for periodic data set (hemispheric or not)
299
      hem = 0
300
      if (per.eq.0.) then
301
         delta=xmax-xmin-360.
302
         if (abs(delta+dx).lt.eps) then               ! Program aborts: arrays must be closed
303
            goto 992
304
        else if (abs(delta).lt.eps) then              ! Periodic and hemispheric
305
           hem=1
306
           per=360.
307
        endif
308
      else                                            ! Periodic and hemispheric
309
         hem=1
310
      endif
311
 
312
C     Allocate memory for some meteorological arrays
313
      allocate(spt0(nx*ny),stat=stat)
314
      if (stat.ne.0) print*,'*** error allocating array spt0 ***'   ! Surface pressure
315
      allocate(spt1(nx*ny),stat=stat)
316
      if (stat.ne.0) print*,'*** error allocating array spt1 ***'
317
      allocate(uut0(nx*ny*nz),stat=stat)
318
      if (stat.ne.0) print*,'*** error allocating array uut0 ***'   ! Zonal wind
319
      allocate(uut1(nx*ny*nz),stat=stat)
320
      if (stat.ne.0) print*,'*** error allocating array uut1 ***'
321
      allocate(vvt0(nx*ny*nz),stat=stat)
322
      if (stat.ne.0) print*,'*** error allocating array vvt0 ***'   ! Meridional wind
323
      allocate(vvt1(nx*ny*nz),stat=stat)
324
      if (stat.ne.0) print*,'*** error allocating array vvt1 ***'
325
      allocate(wwt0(nx*ny*nz),stat=stat)
326
      if (stat.ne.0) print*,'*** error allocating array wwt0 ***'   ! Vertical wind
327
      allocate(wwt1(nx*ny*nz),stat=stat)
328
      if (stat.ne.0) print*,'*** error allocating array wwt1 ***'
329
      allocate(p3t0(nx*ny*nz),stat=stat)
330
      if (stat.ne.0) print*,'*** error allocating array p3t0 ***'   ! Pressure
331
      allocate(p3t1(nx*ny*nz),stat=stat)
332
      if (stat.ne.0) print*,'*** error allocating array p3t1 ***'
333
      allocate(tht0(nx*ny*nz),stat=stat)
334
      if (stat.ne.0) print*,'*** error allocating array tht0 ***'   ! Potential temperature
335
      allocate(tht1(nx*ny*nz),stat=stat)
336
      if (stat.ne.0) print*,'*** error allocating array tht1 ***'
337
      allocate(sth0(nx*ny),stat=stat)
338
      if (stat.ne.0) print*,'*** error allocating array spt0 ***'   ! Surface potential temperature
339
      allocate(sth1(nx*ny),stat=stat)
340
      if (stat.ne.0) print*,'*** error allocating array sth1 ***'
341
 
342
C     Get memory for trajectory arrays
343
      allocate(trainp(ntra,1,ncol),stat=stat)
344
      if (stat.ne.0) print*,'*** error allocating array trainp   ***' ! Input start coordinates
345
      allocate(traout(ntra,ntim,4),stat=stat)
346
      if (stat.ne.0) print*,'*** error allocating array traout   ***' ! Output trajectories
347
      allocate(xx0(ntra),stat=stat)
348
      if (stat.ne.0) print*,'*** error allocating array xx0      ***' ! X position (longitude)
349
      allocate(yy0(ntra),stat=stat)
350
      if (stat.ne.0) print*,'*** error allocating array yy0      ***' ! Y position (latitude)
351
      allocate(pp0(ntra),stat=stat)
352
      if (stat.ne.0) print*,'*** error allocating array pp0      ***' ! Pressure
353
      allocate(leftflag(ntra),stat=stat)
354
      if (stat.ne.0) print*,'*** error allocating array leftflag ***' ! Leaving-domain flag
355
      allocate(theta(ntra),stat=stat)
356
      if (stat.ne.0) print*,'*** error allocating array theta ***'    ! Potential temperature for isentropic trajectories
357
      allocate(zindex(ntra),stat=stat)
358
      if (stat.ne.0) print*,'*** error allocating array kind ***'     ! Vertical index for model-level trajectories
359
 
360
c     Write some status information
361
      print*,'---- CONSTANT GRID PARAMETERS ---------------------------'
362
      print*
363
      print*,'  xmin,xmax     : ',xmin,xmax
364
      print*,'  ymin,ymax     : ',ymin,ymax
365
      print*,'  dx,dy         : ',dx,dy
366
      print*,'  pollon,pollat : ',pollon,pollat
367
      print*,'  nx,ny,nz      : ',nx,ny,nz
368
      print*,'  per, hem      : ',per,hem
369
      print*
370
 
371
c     --------------------------------------------------------------------
372
c     Initialize the trajectory calculation
373
c     --------------------------------------------------------------------
374
 
375
c     Read start coordinates from file - Format (lon,lat,lev)
376
      if (inpmode.eq.-1) then
377
         open(fid,file=strname)
378
          do i=1,ntra
379
             read(fid,*) xx0(i),yy0(i),pp0(i)
380
          enddo
381
         close(fid)
382
 
383
c     Read start coordinates from trajectory file - check consistency of ref time
384
      else
385
         call ropen_tra(cdfid,strname,ntra,1,ncol,reftmp,vars,inpmode)
386
         call read_tra (cdfid,trainp,ntra,1,ncol,inpmode)
387
         do i=1,ntra
388
            time   = trainp(i,1,1)
389
            xx0(i) = trainp(i,1,2) 
390
            yy0(i) = trainp(i,1,3) 
391
            pp0(i) = trainp(i,1,4) 
392
         enddo
393
         call close_tra(cdfid,inpmode)
394
 
395
         if ( ( reftime(1).ne.reftmp(1) ).or.
396
     >        ( reftime(2).ne.reftmp(2) ).or.
397
     >        ( reftime(3).ne.reftmp(3) ).or.
398
     >        ( reftime(4).ne.reftmp(4) ).or.
399
     >        ( reftime(5).ne.reftmp(5) ) )
400
     >   then
401
            print*,' WARNING: Inconsistent reference times'
402
            write(*,'(5i8)') (reftime(i),i=1,5)
403
            write(*,'(5i8)') (reftmp (i),i=1,5)
404
            print*,'Enter a key to proceed...'
405
            stop
406
         endif
407
      endif
408
 
409
c     Set sign of time range
410
      reftime(6) = fbflag * reftime(6)
411
 
412
c     Write some status information
413
      print*,'---- REFERENCE DATE---------- ---------------------------'
414
      print*
415
      print*,' Reference time (year)  :',reftime(1)
416
      print*,'                (month) :',reftime(2)
417
      print*,'                (day)   :',reftime(3)
418
      print*,'                (hour)  :',reftime(4)
419
      print*,'                (min)   :',reftime(5)
420
      print*,' Time range             :',reftime(6),' min'
421
      print*
422
 
423
C     Save starting positions 
424
      itim = 1
425
      do i=1,ntra
426
         traout(i,itim,1) = 0.
427
         traout(i,itim,2) = xx0(i)
428
         traout(i,itim,3) = yy0(i)
429
         traout(i,itim,4) = pp0(i)
430
      enddo
431
 
432
c     Init the flag and the counter for trajectories leaving the domain
433
      leftcount=0
434
      do i=1,ntra
435
         leftflag(i)=0
436
      enddo
437
 
438
C     Convert time shifts <tst,ten> from <hh.mm> into fractional time
439
      call hhmm2frac(tst,frac)
440
      tst = frac
441
      call hhmm2frac(ten,frac)
442
      ten = frac
443
 
444
c     Get 3D and surface pressure from first data file 
445
      filename = prefix//dat(1)
446
      call input_open (fid,filename)
447
      if ( modlev.eq.'no' ) then 
448
         varname = 'P'
449
         call input_grid
450
     >       (fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny,
451
     >        tload,pollon,pollat,p3t1,spt1,nz,ak,bk,stagz,timecheck)
452
      else  
453
         varname = 'P.ML'
454
         call input_grid
455
     >       (fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny,
456
     >        tload,pollon,pollat,p3t1,spt1,nz,ak,bk,stagz,timecheck)
457
      endif
458
      call input_close(fid)
459
 
460
 
461
c     Check that all starting positions are above topography    
462
      do i=1,ntra
463
 
464
C       Interpolate surface pressure to actual position (from first input file)
465
        x1 = xx0(i)
466
        y1 = yy0(i)
467
        call get_index4 (xind,yind,pind,x1,y1,1050.,0.,
468
     >                   p3t1,p3t1,spt1,spt1,3,
469
     >                   nx,ny,nz,xmin,ymin,dx,dy,mdv)
470
        sp = int_index4 (spt1,spt1,nx,ny,1,xind,yind,1.,0.,mdv)
471
 
472
cc       Decide whether to keep the trajectory
473
        if ( pp0(i).gt.sp ) then
474
            write(*,'(a30,3f10.2)')
475
     >               'WARNING: starting point below topography ',
476
     >               xx0(i),yy0(i),pp0(i)
477
            leftflag(i) = 1
478
        endif
479
 
480
      enddo
481
 
482
 
483
c     Special handling for isentropic trajectories - read potential
484
c     temperature from S file or calculate it based on temperature and
485
c     pressure from P file; then, calculate for each trajectory its
486
c     potential temperature - which will stay fixed over time
487
      if ( isen.eq.'yes' ) then
488
 
489
c         Get potential temperature from S file
490
          if ( thons.eq.1 ) then
491
              filename = srefix//dat(1)
492
              print*,' TH <- ',trim(filename)
493
              call input_open (fid,filename)
494
              varname='TH'
495
              call input_wind
496
     >            (fid,varname,tht1,tload,stagz,mdv,
497
     >              xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
498
              call input_close(fid)
499
 
500
c         Calculate potential temperature from P file
501
          else
502
              filename = prefix//dat(1)
503
              print*,' TH = T * (1000/P)^RDCP <- ',trim(filename)
504
              call input_open (fid,filename)
505
              varname='T'
506
              call input_wind
507
     >            (fid,varname,tht1,tload,stagz,mdv,
508
     >              xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
509
              call input_close(fid)
510
              do i=1,nx*ny*nz
511
                if (tht1(i).lt.100.) then
512
                   tht1(i)=(tht1(i)+tzero)*( (1000./p3t1(i))**rdcp )
513
                else
514
                   tht1(i)=tht1(i)*( (1000./p3t1(i))**rdcp )
515
                endif
516
              enddo
517
          endif
518
 
519
c         Take surface potential temperature from lowest level
520
          do i=1,nx*ny
521
            sth1(i) = tht1(i)
522
          enddo
523
 
524
c         Calculate now the potential temperature of all trajectories
525
          do i=1,ntra
526
 
527
             x1 = xx0(i)
528
             y1 = yy0(i)
529
             p1 = pp0(i)
530
             call get_index4 (xind,yind,pind,x1,y1,p1,0.,
531
     >                 p3t1,p3t1,spt1,spt1,3,
532
     >                 nx,ny,nz,xmin,ymin,dx,dy,mdv)
533
             theta(i) =
534
     >        int_index4(tht1,tht1,nx,ny,nz,xind,yind,pind,0.,mdv)
535
 
536
          enddo
537
 
538
c         Write info about theta range of starting positions
539
          thetamin = theta(1)
540
          thetamax = theta(1)
541
          do i=2,ntra
542
            if ( theta(i).gt.thetamax ) thetamax = theta(i)
543
            if ( theta(i).lt.thetamin ) thetamin = theta(i)
544
          enddo
545
 
546
c         Write some status information
547
          print*
548
          print*,
549
     >     '---- THETA RANGE OF ISENTROPIC TRAJECTORIES -------------'
550
          print*
551
          print*,' Theta(min)             :',thetamin
552
          print*,' Theta(max)             :',thetamax
553
 
554
      endif
555
 
556
c     Special handling for model-level (2D) trajectories - get the
557
c     vertical index for each trajectory - which will remain fixed
558
      if ( modlev.eq.'yes' ) then
559
          do i=1,ntra
560
             x1 = xx0(i)
561
             y1 = yy0(i)
562
             p1 = pp0(i)
563
             call get_index4 (xind,yind,pind,x1,y1,p1,0.,
564
     >                 p3t1,p3t1,spt1,spt1,3,
565
     >                 nx,ny,nz,xmin,ymin,dx,dy,mdv)
566
             zindex(i) = pind
567
          enddo
568
 
569
          do i=1,nz
570
            print*,i,p3t1(189+(141-1)*nx+(i-1)*nx*ny)
571
          enddo
572
             print*,x1,y1,p1
573
             print*,xind,yind,pind
574
 
575
c         Write info about zindex range of starting positions
576
          zindexmin = zindex(1)
577
          zindexmax = zindex(1)
578
          do i=2,ntra
579
            if ( zindex(i).gt.zindexmax ) zindexmax = zindex(i)
580
            if ( zindex(i).lt.zindexmin ) zindexmin = zindex(i)
581
          enddo
582
 
583
c         Write some status information
584
          print*
585
          print*,
586
     >     '---- INDEX RANGE OF MODEL-LEVEL TRAJECTORIES ------------'
587
          print*
588
          print*,' Zindex(min)            :',zindexmin
589
          print*,' Zindex(max)            :',zindexmax
590
 
591
 
592
      endif
593
 
594
c     -----------------------------------------------------------------------
595
c     Loop to calculate trajectories
596
c     -----------------------------------------------------------------------
597
 
598
c     Write some status information
599
      print*
600
      print*,'---- TRAJECTORIES ----------- ---------------------------'
601
      print*    
602
 
603
C     Set the time for the first data file (depending on forward/backward mode)
604
      if (fbflag.eq.1) then
605
        tstart = -tst
606
      else
607
        tstart = tst
608
      endif
609
 
610
c     Set the minute counter for output
611
      wstep = 0
612
 
613
c     Read wind fields and vertical grid from first file
614
      filename = prefix//dat(1)
615
 
616
      call frac2hhmm(tstart,tload)
617
 
618
      write(*,'(a16,a20,f9.2)') '  (file,time) : ',
619
     >                       trim(filename),tload
620
 
621
      call input_open (fid,filename)
622
 
623
      varname='U'                                      ! U
624
      call input_wind 
625
     >    (fid,varname,uut1,tload,stagz,mdv,
626
     >     xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
627
 
628
      varname='V'                                      ! V
629
      call input_wind 
630
     >    (fid,varname,vvt1,tload,stagz,mdv,
631
     >     xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
632
 
633
      if ( (modlev.eq.'no').and.(isen.eq.'no') ) then
634
         varname='OMEGA'                                  ! OMEGA
635
         call input_wind
636
     >       (fid,varname,wwt1,tload,stagz,mdv,
637
     >        xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
638
      endif
639
 
640
      if ( modlev.eq.'no' ) then
641
         call input_grid                                  ! GRID - AK,BK -> P
642
     >       (fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny,
643
     >        tload,pollon,pollat,p3t1,spt1,nz,ak,bk,stagz,timecheck)
644
      else
645
         varname='P.ML'                                   ! GRID - P,PS
646
         call input_grid                                  !
647
     >       (fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny,
648
     >        tload,pollon,pollat,p3t1,spt1,nz,ak,bk,stagz,timecheck)
649
      endif
650
 
651
      call input_close(fid)
652
 
653
c     Special handling for isentropic trajectories - read potential
654
c     temperature from S file or calculate it based on temperature and
655
c     pressure from P file
656
      if ( isen.eq.'yes' ) then
657
 
658
c         Get potential temperature from S file
659
          if ( thons.eq.1 ) then
660
              filename = srefix//dat(1)
661
              print*,' TH <- ',trim(filename)
662
              call input_open (fid,filename)
663
              varname='TH'
664
              call input_wind
665
     >            (fid,varname,tht1,tload,stagz,mdv,
666
     >              xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
667
              call input_close(fid)
668
 
669
c         Calculate potential temperature from P file
670
          else
671
              filename = prefix//dat(1)
672
              print*,' TH = T * (1000/P)^RDCP <- ',trim(filename)
673
              call input_open (fid,filename)
674
              varname='T'
675
              call input_wind
676
     >            (fid,varname,tht1,tload,stagz,mdv,
677
     >              xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
678
              call input_close(fid)
679
              do i=1,nx*ny*nz
680
                if (tht1(i).lt.100.) then
681
                   tht1(i)=(tht1(i)+tzero)*( (1000./p3t1(i))**rdcp )
682
                else
683
                   tht1(i)=tht1(i)*( (1000./p3t1(i))**rdcp )
684
                endif
685
              enddo
686
          endif
687
 
688
c         Take surface potential temperature from lowest level
689
          do i=1,nx*ny
690
            sth1(i) = tht1(i)
691
          enddo
692
      endif
693
 
694
c     Loop over all input files (time step is <timeinc>)
695
      do itm=1,numdat-1
696
 
697
c       Calculate actual and next time
698
        time0 = tstart+real(itm-1)*timeinc*fbflag
699
        time1 = time0+timeinc*fbflag
700
 
701
c       Copy old velocities and pressure fields to new ones       
702
        do i=1,nx*ny*nz
703
           uut0(i)=uut1(i)
704
           vvt0(i)=vvt1(i)
705
           wwt0(i)=wwt1(i)
706
           p3t0(i)=p3t1(i)
707
           tht0(i)=tht1(i)
708
c           print*,'i',uut1(i),vvt1(i),p3t1(i)
709
        enddo
710
        do i=1,nx*ny
711
           spt0(i)=spt1(i)
712
           sth0(i)=sth1(i)
713
        enddo
714
 
715
c       Read wind fields and surface pressure at next time
716
        filename = prefix//dat(itm+1)
717
 
718
        call frac2hhmm(time1,tload)
719
        write(*,'(a16,a20,f9.2)') '  (file,time) : ',
720
     >                          trim(filename),tload
721
 
722
        call input_open (fid,filename)
723
 
724
        varname='U'                                     ! U
725
        call input_wind 
726
     >       (fid,varname,uut1,tload,stagz,mdv,
727
     >        xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
728
 
729
        varname='V'                                     ! V
730
        call input_wind 
731
     >       (fid,varname,vvt1,tload,stagz,mdv,
732
     >        xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
733
 
734
        if ( (modlev.eq.'no').and.(isen.eq.'no') ) then
735
           varname='OMEGA'                              ! OMEGA
736
           call input_wind
737
     >          (fid,varname,wwt1,tload,stagz,mdv,
738
     >           xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
739
        endif
740
 
741
        if ( modlev.eq.'no' ) then
742
           call input_grid                                  ! GRID - AK,NK -> P
743
     >          (fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny,
744
     >           tload,pollon,pollat,p3t1,spt1,nz,ak,bk,stagz,timecheck)
745
        else
746
           varname='P.ML'                                   ! GRID - P,PS
747
           call input_grid                                  !
748
     >       (fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny,
749
     >        tload,pollon,pollat,p3t1,spt1,nz,ak,bk,stagz,timecheck)
750
        endif
751
 
752
        call input_close(fid)
753
 
754
c     Special handling for isentropic trajectories - read potential
755
c     temperature from S file or calculate it based on temperature and
756
c     pressure from P file
757
      if ( isen.eq.'yes' ) then
758
 
759
c         Get TH from S file
760
          if ( thons.eq.1 ) then
761
              filename = srefix//dat(itm+1)
762
              print*,' TH <- ',trim(filename)
763
              call input_open (fid,filename)
764
              varname='TH'
765
              call input_wind
766
     >            (fid,varname,tht1,tload,stagz,mdv,
767
     >              xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
768
              call input_close(fid)
769
 
770
c         Calculate potential temperature from P file
771
          else
772
              filename = prefix//dat(itm+1)
773
              print*,' TH = T * (1000/P)^RDCP <- ',trim(filename)
774
              call input_open (fid,filename)
775
              varname='T'
776
              call input_wind
777
     >            (fid,varname,tht1,tload,stagz,mdv,
778
     >              xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
779
              call input_close(fid)
780
              do i=1,nx*ny*nz
781
                if (tht1(i).lt.100.) then
782
                   tht1(i)=(tht1(i)+tzero)*( (1000./p3t1(i))**rdcp )
783
                else
784
                   tht1(i)=tht1(i)*( (1000./p3t1(i))**rdcp )
785
                endif
786
              enddo
787
          endif
788
 
789
c         Take surface potential temperature from lowest level
790
          do i=1,nx*ny
791
            sth1(i) = tht1(i)
792
          enddo
793
      endif
794
 
795
C       Determine the first and last loop indices
796
        if (numdat.eq.2) then
797
          filo = nint(tst/ts)+1
798
          lalo = nint((timeinc-ten)/ts) 
799
        elseif ( itm.eq.1 ) then
800
          filo = nint(tst/ts)+1
801
          lalo = nint(timeinc/ts)
802
        else if (itm.eq.numdat-1) then
803
          filo = 1
804
          lalo = nint((timeinc-ten)/ts)
805
        else
806
          filo = 1
807
          lalo = nint(timeinc/ts)
808
        endif
809
 
810
c       Split the interval <timeinc> into computational time steps <ts>
811
        do iloop=filo,lalo
812
 
813
C         Calculate relative time position in the interval timeinc (0=beginning, 1=end)
814
          reltpos0 = ((real(iloop)-1.)*ts)/timeinc
815
          reltpos1 = real(iloop)*ts/timeinc
816
 
817
c         Timestep for all trajectories
818
          do i=1,ntra
819
 
820
C           Check if trajectory has already left the data domain
821
            if (leftflag(i).ne.1) then	
822
 
823
c             3D: Iterative Euler timestep (x0,y0,p0 -> x1,y1,p1)
824
              if ( (imethod.eq.1   ).and.
825
     >             (isen   .eq.'no').and.
826
     >             (modlev .eq.'no') )
827
     >        then
828
                 call euler_3d(
829
     >               xx1,yy1,pp1,leftflag(i),
830
     >               xx0(i),yy0(i),pp0(i),reltpos0,reltpos1,
831
     >               ts*3600,numit,jflag,mdv,wfactor,fbflag,
832
     >               spt0,spt1,p3t0,p3t1,uut0,uut1,vvt0,vvt1,wwt0,wwt1,
833
     >               xmin,ymin,dx,dy,per,hem,nx,ny,nz)
834
 
835
c             3D: Runge-Kutta timestep (x0,y0,p0 -> x1,y1,p1)
836
              else if ( (imethod.eq.2   ).and.
837
     >                  (isen   .eq.'no').and.
838
     >                  (modlev .eq.'no') )
839
     >        then
840
                 call runge(
841
     >               xx1,yy1,pp1,leftflag(i),
842
     >               xx0(i),yy0(i),pp0(i),reltpos0,reltpos1,
843
     >               ts*3600,numit,jflag,mdv,wfactor,fbflag,
844
     >               spt0,spt1,p3t0,p3t1,uut0,uut1,vvt0,vvt1,wwt0,wwt1,
845
     >               xmin,ymin,dx,dy,per,hem,nx,ny,nz)
846
 
847
c             ISENTROPIC: Iterative Euler timestep (x0,y0,p0 -> x1,y1,p1)
848
              else if ( (imethod.eq.1    ).and.
849
     >                  (isen   .eq.'yes').and.
850
     >                  (modlev .eq.'no' ) )
851
     >        then
852
                 call euler_isen(
853
     >               xx1,yy1,pp1,leftflag(i),
854
     >               xx0(i),yy0(i),pp0(i),theta(i),reltpos0,reltpos1,
855
     >               ts*3600,numit,jflag,mdv,wfactor,fbflag,
856
     >               spt0,spt1,p3t0,p3t1,uut0,uut1,vvt0,vvt1,
857
     >               sth0,sth1,tht0,tht1,
858
     >               xmin,ymin,dx,dy,per,hem,nx,ny,nz)
859
 
860
c             MODEL-LEVEL (2D): Iterative Euler timestep (x0,y0,p0 -> x1,y1,p1)
861
              else if ( (imethod.eq.1    ).and.
862
     >                  (isen   .eq.'no' ).and.
863
     >                  (modlev .eq.'yes') )
864
     >        then
865
                 call euler_2d(
866
     >               xx1,yy1,pp1,leftflag(i),
867
     >               xx0(i),yy0(i),pp0(i),zindex(i),reltpos0,reltpos1,
868
     >               ts*3600,numit,jflag,mdv,wfactor,fbflag,
869
     >               spt0,spt1,p3t0,p3t1,uut0,uut1,vvt0,vvt1,
870
     >               xmin,ymin,dx,dy,per,hem,nx,ny,nz)
871
 
872
              endif
873
 
874
c             Update trajectory position, or increase number of trajectories leaving domain
875
              if (leftflag(i).eq.1) then
876
                leftcount=leftcount+1
877
                if ( leftcount.lt.10 ) then
878
                   print*,'     -> Trajectory ',i,' leaves domain'
879
                elseif ( leftcount.eq.10 ) then
880
                   print*,'     -> N>=10 trajectories leave domain'
881
                endif
882
              else
883
                xx0(i)=xx1      
884
                yy0(i)=yy1
885
                pp0(i)=pp1
886
              endif
887
 
888
c          Trajectory has already left data domain (mark as <mdv>)
889
           else
890
              xx0(i)=mdv
891
              yy0(i)=mdv
892
              pp0(i)=mdv
893
 
894
           endif
895
 
896
          enddo
897
 
898
C         Save positions only every deltout minutes
899
          delta = aint(iloop*60*ts/deltout)-iloop*60*ts/deltout
900
          if (abs(delta).lt.eps) then
901
c          wstep = wstep + abs(ts)
902
c          if ( mod(wstep,deltout).eq.0 ) then
903
            time = time0+reltpos1*timeinc*fbflag
904
            itim = itim + 1
905
            if ( itim.le.ntim ) then
906
              do i=1,ntra
907
                 call frac2hhmm(time,tload)
908
                 traout(i,itim,1) = tload
909
                 traout(i,itim,2) = xx0(i)
910
                 traout(i,itim,3) = yy0(i)
911
                 traout(i,itim,4) = pp0(i)
912
              enddo
913
            endif
914
          endif
915
 
916
        enddo
917
 
918
      enddo
919
 
920
c     Write trajectory file
921
      vars(1)  ='time'
922
      vars(2)  ='lon'
923
      vars(3)  ='lat'
924
      vars(4)  ='p'
925
      call wopen_tra(cdfid,cdfname,ntra,ntim,4,reftime,vars,outmode)
926
      call write_tra(cdfid,traout,ntra,ntim,4,outmode)
927
      call close_tra(cdfid,outmode)   
928
 
929
c     Write some status information, and end of program message
930
      print*  
931
      print*,'---- STATUS INFORMATION --------------------------------'
932
      print*
933
      print*,'  #leaving domain    ', leftcount
934
      print*,'  #staying in domain ', ntra-leftcount
935
      print*
936
      print*,'              *** END OF PROGRAM CALTRA ***'
937
      print*,'========================================================='
938
 
939
      stop
940
 
941
c     ------------------------------------------------------------------
942
c     Exception handling
943
c     ------------------------------------------------------------------
944
 
945
 991  write(*,*) '*** ERROR: all start points outside the data domain'
946
      call exit(1)
947
 
948
 992  write(*,*) '*** ERROR: close arrays on files (prog. closear)'
949
      call exit(1)
950
 
951
 993  write(*,*) '*** ERROR: problems with array size'
952
      call exit(1)
953
 
954
      end 
955
 
956
 
957
c     *******************************************************************
958
c     * Time step : either Euler or Runge-Kutta                         *
959
c     *******************************************************************
960
 
961
C     Time-step from (x0,y0,p0) to (x1,y1,p1)
962
C
963
C     (x0,y0,p0) input	coordinates (long,lat,p) for starting point
964
C     (x1,y1,p1) output	coordinates (long,lat,p) for end point
965
C     deltat	 input	timestep in seconds
966
C     numit	 input	number of iterations
967
C     jump	 input  flag (=1 trajectories don't enter the ground)
968
C     left	 output	flag (=1 if trajectory leaves data domain)
969
 
970
c     -------------------------------------------------------------------
971
c     Iterative Euler time step (KINEMATIC 3D TRAJECTORIES)
972
c     -------------------------------------------------------------------
973
 
974
      subroutine euler_3d(x1,y1,p1,left,x0,y0,p0,reltpos0,reltpos1,
975
     >                deltat,numit,jump,mdv,wfactor,fbflag,
976
     >		          spt0,spt1,p3d0,p3d1,uut0,uut1,vvt0,vvt1,wwt0,wwt1,
977
     >                xmin,ymin,dx,dy,per,hem,nx,ny,nz)
978
 
979
      implicit none
980
 
981
c     Declaration of subroutine parameters
982
      integer      nx,ny,nz
983
      real         x1,y1,p1
984
      integer      left
985
      real	   x0,y0,p0
986
      real         reltpos0,reltpos1
987
      real   	   deltat
988
      integer      numit
989
      integer      jump
990
      real         wfactor
991
      integer      fbflag
992
      real     	   spt0(nx*ny)   ,spt1(nx*ny)
993
      real         uut0(nx*ny*nz),uut1(nx*ny*nz)
994
      real 	   vvt0(nx*ny*nz),vvt1(nx*ny*nz)
995
      real         wwt0(nx*ny*nz),wwt1(nx*ny*nz)
996
      real         p3d0(nx*ny*nz),p3d1(nx*ny*nz)
997
      real         xmin,ymin,dx,dy
998
      real         per
999
      integer      hem
1000
      real         mdv
1001
 
1002
c     Numerical and physical constants
1003
      real         deltay
1004
      parameter    (deltay=1.112E5)  ! Distance in m between 2 lat circles
1005
      real         pi                       
1006
      parameter    (pi=3.1415927)    ! Pi
1007
 
1008
c     Auxiliary variables
1009
      real         xmax,ymax
1010
      real	   xind,yind,pind
1011
      real	   u0,v0,w0,u1,v1,w1,u,v,w,sp
1012
      integer	   icount
1013
      character    ch
1014
 
1015
c     Externals    
1016
      real         int_index4
1017
      external     int_index4
1018
 
1019
c     Reset the flag for domain-leaving
1020
      left=0
1021
 
1022
c     Set the esat-north bounray of the domain
1023
      xmax = xmin+real(nx-1)*dx
1024
      ymax = ymin+real(ny-1)*dy
1025
 
1026
C     Interpolate wind fields to starting position (x0,y0,p0)
1027
      call get_index4 (xind,yind,pind,x0,y0,p0,reltpos0,
1028
     >                 p3d0,p3d1,spt0,spt1,3,
1029
     >                 nx,ny,nz,xmin,ymin,dx,dy,mdv)
1030
      u0 = int_index4(uut0,uut1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
1031
      v0 = int_index4(vvt0,vvt1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
1032
      w0 = int_index4(wwt0,wwt1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
1033
 
1034
c     Force the near-surface wind to zero
1035
      if (pind.lt.1.) w0=w0*pind
1036
 
1037
C     For first iteration take ending position equal to starting position
1038
      x1=x0
1039
      y1=y0
1040
      p1=p0
1041
 
1042
C     Iterative calculation of new position
1043
      do icount=1,numit
1044
 
1045
C        Calculate new winds for advection
1046
         call get_index4 (xind,yind,pind,x1,y1,p1,reltpos1,
1047
     >                    p3d0,p3d1,spt0,spt1,3,
1048
     >                    nx,ny,nz,xmin,ymin,dx,dy,mdv)
1049
         u1 = int_index4(uut0,uut1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
1050
         v1 = int_index4(vvt0,vvt1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
1051
         w1 = int_index4(wwt0,wwt1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
1052
 
1053
c        Force the near-surface wind to zero
1054
         if (pind.lt.1.) w1=w1*pind
1055
 
1056
c        Get the new velocity in between
1057
         u=(u0+u1)/2.
1058
         v=(v0+v1)/2.
1059
         w=(w0+w1)/2.
1060
 
1061
C        Calculate new positions
1062
         x1 = x0 + fbflag*u*deltat/(deltay*cos(y0*pi/180.))
1063
         y1 = y0 + fbflag*v*deltat/deltay
1064
         p1 = p0 + fbflag*wfactor*w*deltat/100.
1065
 
1066
c       Handle pole problems (crossing and near pole trajectory)
1067
        if ((hem.eq.1).and.(y1.gt.90.)) then
1068
          y1=180.-y1
1069
          x1=x1+per/2.
1070
        endif
1071
        if ((hem.eq.1).and.(y1.lt.-90.)) then
1072
          y1=-180.-y1
1073
          x1=x1+per/2.
1074
        endif
1075
        if (y1.gt.89.99) then
1076
           y1=89.99
1077
        endif
1078
 
1079
c       Handle crossings of the dateline
1080
        if ((hem.eq.1).and.(x1.gt.xmin+per-dx)) then
1081
           x1=xmin+amod(x1-xmin,per)
1082
        endif
1083
        if ((hem.eq.1).and.(x1.lt.xmin)) then
1084
           x1=xmin+per+amod(x1-xmin,per)
1085
        endif
1086
 
1087
C       Interpolate surface pressure to actual position
1088
        call get_index4 (xind,yind,pind,x1,y1,1050.,reltpos1,
1089
     >                   p3d0,p3d1,spt0,spt1,3,
1090
     >                   nx,ny,nz,xmin,ymin,dx,dy,mdv)
1091
        sp = int_index4 (spt0,spt1,nx,ny,1,xind,yind,1.,reltpos1,mdv)
1092
 
1093
c       Handle trajectories which cross the lower boundary (jump flag)
1094
        if ((jump.eq.1).and.(p1.gt.sp)) p1=sp-10.
1095
 
1096
C       Check if trajectory leaves data domain
1097
        if ( ( (hem.eq.0).and.(x1.lt.xmin)    ).or.
1098
     >       ( (hem.eq.0).and.(x1.gt.xmax-dx) ).or.
1099
     >         (y1.lt.ymin).or.(y1.gt.ymax).or.(p1.gt.sp) )
1100
     >  then
1101
          left=1
1102
          goto 100
1103
        endif
1104
 
1105
      enddo
1106
 
1107
c     Exit point for subroutine
1108
 100  continue
1109
 
1110
      return
1111
 
1112
      end
1113
 
1114
c     -------------------------------------------------------------------
1115
c     Iterative Euler time step (ISENTROPIC)
1116
c     -------------------------------------------------------------------
1117
 
1118
      subroutine euler_isen(x1,y1,p1,left,x0,y0,p0,theta,
1119
     >                reltpos0,reltpos1,
1120
     >                deltat,numit,jump,mdv,wfactor,fbflag,
1121
     >                spt0,spt1,p3d0,p3d1,uut0,uut1,vvt0,vvt1,
1122
     >                sth0,sth1,tht0,tht1,
1123
     >                xmin,ymin,dx,dy,per,hem,nx,ny,nz)
1124
 
1125
      implicit none
1126
 
1127
c     Declaration of subroutine parameters
1128
      integer      nx,ny,nz
1129
      real         x1,y1,p1
1130
      integer      left
1131
      real         x0,y0,p0
1132
      real         reltpos0,reltpos1
1133
      real         deltat
1134
      integer      numit
1135
      integer      jump
1136
      real         wfactor
1137
      integer      fbflag
1138
      real         spt0(nx*ny)   ,spt1(nx*ny)
1139
      real         sth0(nx*ny)   ,sth1(nx*ny)
1140
      real         uut0(nx*ny*nz),uut1(nx*ny*nz)
1141
      real         vvt0(nx*ny*nz),vvt1(nx*ny*nz)
1142
      real         p3d0(nx*ny*nz),p3d1(nx*ny*nz)
1143
      real         tht0(nx*ny*nz),tht1(nx*ny*nz)
1144
      real         xmin,ymin,dx,dy
1145
      real         per
1146
      integer      hem
1147
      real         mdv
1148
      real         theta
1149
 
1150
c     Numerical and physical constants
1151
      real         deltay
1152
      parameter    (deltay=1.112E5)  ! Distance in m between 2 lat circles
1153
      real         pi
1154
      parameter    (pi=3.1415927)    ! Pi
1155
 
1156
c     Auxiliary variables
1157
      real         xmax,ymax
1158
      real         xind,yind,pind
1159
      real         u0,v0,w0,u1,v1,w1,u,v,w,sp
1160
      integer      icount
1161
      character    ch
1162
 
1163
c     Externals
1164
      real         int_index4
1165
      external     int_index4
1166
 
1167
c     Reset the flag for domain-leaving
1168
      left=0
1169
 
1170
c     Set the esat-north bounray of the domain
1171
      xmax = xmin+real(nx-1)*dx
1172
      ymax = ymin+real(ny-1)*dy
1173
 
1174
C     Interpolate wind fields to starting position (x0,y0,p0)
1175
      call get_index4 (xind,yind,pind,x0,y0,p0,reltpos0,
1176
     >                 p3d0,p3d1,spt0,spt1,3,
1177
     >                 nx,ny,nz,xmin,ymin,dx,dy,mdv)
1178
      u0 = int_index4(uut0,uut1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
1179
      v0 = int_index4(vvt0,vvt1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
1180
 
1181
C     For first iteration take ending position equal to starting position
1182
      x1=x0
1183
      y1=y0
1184
      p1=p0
1185
 
1186
C     Iterative calculation of new position
1187
      do icount=1,numit
1188
 
1189
C        Calculate new winds for advection
1190
         call get_index4 (xind,yind,pind,x1,y1,p1,reltpos1,
1191
     >                    p3d0,p3d1,spt0,spt1,3,
1192
     >                    nx,ny,nz,xmin,ymin,dx,dy,mdv)
1193
         u1 = int_index4(uut0,uut1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
1194
         v1 = int_index4(vvt0,vvt1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
1195
 
1196
c        Get the new velocity in between
1197
         u=(u0+u1)/2.
1198
         v=(v0+v1)/2.
1199
 
1200
C        Calculate new positions
1201
         x1 = x0 + fbflag*u*deltat/(deltay*cos(y0*pi/180.))
1202
         y1 = y0 + fbflag*v*deltat/deltay
1203
 
1204
c        Get the pressure on the isentropic surface at the new position
1205
         call get_index4 (xind,yind,pind,x1,y1,theta,reltpos1,
1206
     >                    tht0,tht1,sth0,sth1,1,
1207
     >                    nx,ny,nz,xmin,ymin,dx,dy,mdv)
1208
         p1 = int_index4(p3d0,p3d1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
1209
 
1210
c       Handle pole problems (crossing and near pole trajectory)
1211
        if ((hem.eq.1).and.(y1.gt.90.)) then
1212
          y1=180.-y1
1213
          x1=x1+per/2.
1214
        endif
1215
        if ((hem.eq.1).and.(y1.lt.-90.)) then
1216
          y1=-180.-y1
1217
          x1=x1+per/2.
1218
        endif
1219
        if (y1.gt.89.99) then
1220
           y1=89.99
1221
        endif
1222
 
1223
c       Handle crossings of the dateline
1224
        if ((hem.eq.1).and.(x1.gt.xmin+per-dx)) then
1225
           x1=xmin+amod(x1-xmin,per)
1226
        endif
1227
        if ((hem.eq.1).and.(x1.lt.xmin)) then
1228
           x1=xmin+per+amod(x1-xmin,per)
1229
        endif
1230
 
1231
C       Interpolate surface pressure to actual position
1232
        call get_index4 (xind,yind,pind,x1,y1,1050.,reltpos1,
1233
     >                   p3d0,p3d1,spt0,spt1,3,
1234
     >                   nx,ny,nz,xmin,ymin,dx,dy,mdv)
1235
        sp = int_index4 (spt0,spt1,nx,ny,1,xind,yind,1.,reltpos1,mdv)
1236
 
1237
c       Handle trajectories which cross the lower boundary (jump flag)
1238
        if ((jump.eq.1).and.(p1.gt.sp)) p1=sp-10.
1239
 
1240
 
1241
 
1242
C       Check if trajectory leaves data domain
1243
        if ( ( (hem.eq.0).and.(x1.lt.xmin)    ).or.
1244
     >       ( (hem.eq.0).and.(x1.gt.xmax-dx) ).or.
1245
     >         (y1.lt.ymin).or.(y1.gt.ymax).or.(p1.gt.sp) )
1246
     >  then
1247
          left=1
1248
          goto 100
1249
        endif
1250
 
1251
      enddo
1252
 
1253
c     Exit point for subroutine
1254
 100  continue
1255
 
1256
      return
1257
 
1258
      end
1259
 
1260
c     -------------------------------------------------------------------
1261
c     Iterative Euler time step (MODEL-LEVEL, 2D)
1262
c     -------------------------------------------------------------------
1263
 
1264
      subroutine euler_2d(x1,y1,p1,left,x0,y0,p0,zindex,
1265
     >                reltpos0,reltpos1,
1266
     >                deltat,numit,jump,mdv,wfactor,fbflag,
1267
     >                spt0,spt1,p3d0,p3d1,uut0,uut1,vvt0,vvt1,
1268
     >                xmin,ymin,dx,dy,per,hem,nx,ny,nz)
1269
 
1270
      implicit none
1271
 
1272
c     Declaration of subroutine parameters
1273
      integer      nx,ny,nz
1274
      real         x1,y1,p1
1275
      integer      left
1276
      real         x0,y0,p0
1277
      real         reltpos0,reltpos1
1278
      real         deltat
1279
      integer      numit
1280
      integer      jump
1281
      real         wfactor
1282
      integer      fbflag
1283
      real         spt0(nx*ny)   ,spt1(nx*ny)
1284
      real         uut0(nx*ny*nz),uut1(nx*ny*nz)
1285
      real         vvt0(nx*ny*nz),vvt1(nx*ny*nz)
1286
      real         p3d0(nx*ny*nz),p3d1(nx*ny*nz)
1287
      real         xmin,ymin,dx,dy
1288
      real         per
1289
      integer      hem
1290
      real         mdv
1291
      real         zindex
1292
 
1293
c     Numerical and physical constants
1294
      real         deltay
1295
      parameter    (deltay=1.112E5)  ! Distance in m between 2 lat circles
1296
      real         pi
1297
      parameter    (pi=3.1415927)    ! Pi
1298
      real         eps
1299
      parameter    (eps=0.001)
1300
 
1301
c     Auxiliary variables
1302
      real         xmax,ymax
1303
      real         xind,yind,pind
1304
      real         u0,v0,w0,u1,v1,w1,u,v,w,sp
1305
      integer      icount
1306
      character    ch
1307
 
1308
c     Externals
1309
      real         int_index4
1310
      external     int_index4
1311
 
1312
c     Reset the flag for domain-leaving
1313
      left=0
1314
 
1315
c     Set the esat-north bounray of the domain
1316
      xmax = xmin+real(nx-1)*dx
1317
      ymax = ymin+real(ny-1)*dy
1318
 
1319
C     Interpolate wind fields to starting position (x0,y0,p0)
1320
      call get_index4 (xind,yind,pind,x0,y0,p0,reltpos0,
1321
     >                 p3d0,p3d1,spt0,spt1,3,
1322
     >                 nx,ny,nz,xmin,ymin,dx,dy,mdv)
1323
      u0 = int_index4(uut0,uut1,nx,ny,nz,xind,yind,zindex,reltpos0,mdv)
1324
      v0 = int_index4(vvt0,vvt1,nx,ny,nz,xind,yind,zindex,reltpos0,mdv)
1325
 
1326
C     For first iteration take ending position equal to starting position
1327
      x1=x0
1328
      y1=y0
1329
      p1=p0
1330
 
1331
C     Iterative calculation of new position
1332
      do icount=1,numit
1333
 
1334
C        Calculate new winds for advection
1335
         call get_index4 (xind,yind,pind,x1,y1,p1,reltpos1,
1336
     >                    p3d0,p3d1,spt0,spt1,3,
1337
     >                    nx,ny,nz,xmin,ymin,dx,dy,mdv)
1338
         u1=int_index4(uut0,uut1,nx,ny,nz,xind,yind,zindex,reltpos1,mdv)
1339
         v1=int_index4(vvt0,vvt1,nx,ny,nz,xind,yind,zindex,reltpos1,mdv)
1340
         if ( abs(u1-mdv).lt.eps ) then
1341
             left = 1
1342
             goto 100
1343
         endif
1344
         if ( abs(v1-mdv).lt.eps ) then
1345
             left = 1
1346
             goto 100
1347
         endif
1348
 
1349
c        Get the new velocity in between
1350
         u=(u0+u1)/2.
1351
         v=(v0+v1)/2.
1352
 
1353
C        Calculate new positions
1354
         x1 = x0 + fbflag*u*deltat/(deltay*cos(y0*pi/180.))
1355
         y1 = y0 + fbflag*v*deltat/deltay
1356
 
1357
c        Get the pressure on the model surface at the new position
1358
         xind = (x1 - xmin ) / dx + 1.
1359
         yind = (y1 - ymin ) / dy + 1.
1360
         p1 =
1361
     >     int_index4(p3d0,p3d1,nx,ny,nz,xind,yind,zindex,reltpos1,mdv)
1362
         if ( abs(p1-mdv).lt.eps ) then
1363
             left = 1
1364
             goto 100
1365
         endif
1366
 
1367
c       Handle pole problems (crossing and near pole trajectory)
1368
        if ((hem.eq.1).and.(y1.gt.90.)) then
1369
          y1=180.-y1
1370
          x1=x1+per/2.
1371
        endif
1372
        if ((hem.eq.1).and.(y1.lt.-90.)) then
1373
          y1=-180.-y1
1374
          x1=x1+per/2.
1375
        endif
1376
        if (y1.gt.89.99) then
1377
           y1=89.99
1378
        endif
1379
 
1380
c       Handle crossings of the dateline
1381
        if ((hem.eq.1).and.(x1.gt.xmin+per-dx)) then
1382
           x1=xmin+amod(x1-xmin,per)
1383
        endif
1384
        if ((hem.eq.1).and.(x1.lt.xmin)) then
1385
           x1=xmin+per+amod(x1-xmin,per)
1386
        endif
1387
 
1388
C       Interpolate surface pressure to actual position
1389
        call get_index4 (xind,yind,pind,x1,y1,1050.,reltpos1,
1390
     >                   p3d0,p3d1,spt0,spt1,3,
1391
     >                   nx,ny,nz,xmin,ymin,dx,dy,mdv)
1392
        sp = int_index4 (spt0,spt1,nx,ny,1,xind,yind,1.,reltpos1,mdv)
1393
 
1394
c       Handle trajectories which cross the lower boundary (jump flag)
1395
        if ((jump.eq.1).and.(p1.gt.sp)) p1=sp-10.
1396
 
1397
C       Check if trajectory leaves data domain
1398
        if ( ( (hem.eq.0).and.(x1.lt.xmin)    ).or.
1399
     >       ( (hem.eq.0).and.(x1.gt.xmax-dx) ).or.
1400
     >         (y1.lt.ymin).or.(y1.gt.ymax).or.(p1.gt.sp) )
1401
     >  then
1402
          left=1
1403
          goto 100
1404
        endif
1405
 
1406
      enddo
1407
 
1408
c     Exit point for subroutine
1409
 100  continue
1410
 
1411
      return
1412
 
1413
      end
1414
 
1415
c     -------------------------------------------------------------------
1416
c     Runge-Kutta (4th order) time-step
1417
c     -------------------------------------------------------------------
1418
 
1419
      subroutine runge(x1,y1,p1,left,x0,y0,p0,reltpos0,reltpos1,
1420
     >                 deltat,numit,jump,mdv,wfactor,fbflag,
1421
     >		       spt0,spt1,p3d0,p3d1,uut0,uut1,vvt0,vvt1,wwt0,wwt1,
1422
     >                 xmin,ymin,dx,dy,per,hem,nx,ny,nz)
1423
 
1424
      implicit none
1425
 
1426
c     Declaration of subroutine parameters
1427
      integer      nx,ny,nz
1428
      real         x1,y1,p1
1429
      integer      left
1430
      real	   x0,y0,p0
1431
      real         reltpos0,reltpos1
1432
      real   	   deltat
1433
      integer      numit
1434
      integer      jump
1435
      real         wfactor
1436
      integer      fbflag
1437
      real     	   spt0(nx*ny)   ,spt1(nx*ny)
1438
      real         uut0(nx*ny*nz),uut1(nx*ny*nz)
1439
      real 	   vvt0(nx*ny*nz),vvt1(nx*ny*nz)
1440
      real         wwt0(nx*ny*nz),wwt1(nx*ny*nz)
1441
      real         p3d0(nx*ny*nz),p3d1(nx*ny*nz)
1442
      real         xmin,ymin,dx,dy
1443
      real         per
1444
      integer      hem
1445
      real         mdv
1446
 
1447
c     Numerical and physical constants
1448
      real         deltay
1449
      parameter    (deltay=1.112E5)  ! Distance in m between 2 lat circles
1450
      real         pi                       
1451
      parameter    (pi=3.1415927)    ! Pi
1452
 
1453
c     Auxiliary variables
1454
      real         xmax,ymax
1455
      real	   xind,yind,pind
1456
      real	   u0,v0,w0,u1,v1,w1,u,v,w,sp
1457
      integer	   icount,n
1458
      real         xs,ys,ps,xk(4),yk(4),pk(4)
1459
      real         reltpos
1460
 
1461
c     Externals    
1462
      real         int_index4
1463
      external     int_index4
1464
 
1465
c     Reset the flag for domain-leaving
1466
      left=0
1467
 
1468
c     Set the esat-north bounray of the domain
1469
      xmax = xmin+real(nx-1)*dx
1470
      ymax = ymin+real(ny-1)*dy
1471
 
1472
c     Apply the Runge Kutta scheme
1473
      do n=1,4
1474
 
1475
c       Get intermediate position and relative time
1476
        if (n.eq.1) then
1477
          xs=0.
1478
          ys=0.
1479
          ps=0.
1480
          reltpos=reltpos0
1481
        else if (n.eq.4) then
1482
          xs=xk(3)
1483
          ys=yk(3)
1484
          ps=pk(3)
1485
          reltpos=reltpos1
1486
        else
1487
          xs=xk(n-1)/2.
1488
          ys=yk(n-1)/2.
1489
          ps=pk(n-1)/2.
1490
          reltpos=(reltpos0+reltpos1)/2.
1491
        endif
1492
 
1493
C       Calculate new winds for advection
1494
        call get_index4 (xind,yind,pind,x0+xs,y0+ys,p0+ps,reltpos,
1495
     >                   p3d0,p3d1,spt0,spt1,3,
1496
     >                   nx,ny,nz,xmin,ymin,dx,dy,mdv)
1497
        u = int_index4 (uut0,uut1,nx,ny,nz,xind,yind,pind,reltpos,mdv)
1498
        v = int_index4 (vvt0,vvt1,nx,ny,nz,xind,yind,pind,reltpos,mdv)
1499
        w = int_index4 (wwt0,wwt1,nx,ny,nz,xind,yind,pind,reltpos,mdv)
1500
 
1501
c       Force the near-surface wind to zero
1502
        if (pind.lt.1.) w1=w1*pind
1503
 
1504
c       Update position and keep them
1505
        xk(n)=fbflag*u*deltat/(deltay*cos(y0*pi/180.))
1506
        yk(n)=fbflag*v*deltat/deltay
1507
        pk(n)=fbflag*w*deltat*wfactor/100.
1508
 
1509
      enddo
1510
 
1511
C     Calculate new positions
1512
      x1=x0+(1./6.)*(xk(1)+2.*xk(2)+2.*xk(3)+xk(4))
1513
      y1=y0+(1./6.)*(yk(1)+2.*yk(2)+2.*yk(3)+yk(4))
1514
      p1=p0+(1./6.)*(pk(1)+2.*pk(2)+2.*pk(3)+pk(4))
1515
 
1516
c     Handle pole problems (crossing and near pole trajectory)
1517
      if ((hem.eq.1).and.(y1.gt.90.)) then
1518
         y1=180.-y1
1519
         x1=x1+per/2.
1520
      endif
1521
      if ((hem.eq.1).and.(y1.lt.-90.)) then
1522
         y1=-180.-y1
1523
         x1=x1+per/2.
1524
      endif
1525
      if (y1.gt.89.99) then
1526
         y1=89.99
1527
      endif
1528
 
1529
c     Handle crossings of the dateline
1530
      if ((hem.eq.1).and.(x1.gt.xmin+per-dx)) then
1531
         x1=xmin+amod(x1-xmin,per)
1532
      endif
1533
      if ((hem.eq.1).and.(x1.lt.xmin)) then
1534
         x1=xmin+per+amod(x1-xmin,per)
1535
      endif
1536
 
1537
 
1538
C     Interpolate surface pressure to actual position
1539
      call get_index4 (xind,yind,pind,x1,y1,1050.,reltpos1,
1540
     >                 p3d0,p3d1,spt0,spt1,3,
1541
     >                 nx,ny,nz,xmin,ymin,dx,dy,mdv)
1542
      sp = int_index4 (spt0,spt1,nx,ny,1,xind,yind,1,reltpos,mdv)
1543
 
1544
c     Handle trajectories which cross the lower boundary (jump flag)
1545
      if ((jump.eq.1).and.(p1.gt.sp)) p1=sp-10.
1546
 
1547
C     Check if trajectory leaves data domain
1548
      if ( ( (hem.eq.0).and.(x1.lt.xmin)    ).or.
1549
     >     ( (hem.eq.0).and.(x1.gt.xmax-dx) ).or.
1550
     >       (y1.lt.ymin).or.(y1.gt.ymax).or.(p1.gt.sp) )
1551
     >then
1552
         left=1
1553
         goto 100
1554
      endif
1555
 
1556
 
1557
c     Exit point fdor subroutine
1558
 100  continue
1559
 
1560
      return
1561
      end