Subversion Repositories lagranto.ecmwf

Rev

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