Subversion Repositories lagranto.ecmwf

Rev

Rev 44 | Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
5 michaesp 1
PROGRAM trace
3 michaesp 2
 
3
  ! ********************************************************************
4
  ! *                                                                  *
5
  ! * Trace fields along trajectories                                  *
6
  ! *                                                                  *
5 michaesp 7
  ! * April 1993: First version (Heini Wernli)                         *
8
  ! * 2008-2009 : Major upgrades (Michael Sprenger)                    *
9
  ! * Mar 2012: Clustering option (Bojan Skerlak)                      *  
10
  ! * Nov 2012: Circle options (")                                     *
11
  ! * Jul 2013: user-defined PV,TH @ clustering mode (")               *
3 michaesp 12
  ! *                                                                  *
13
  ! ********************************************************************
14
 
15
  implicit none
16
 
17
  ! --------------------------------------------------------------------
18
  ! Declaration of parameters
19
  ! --------------------------------------------------------------------
20
 
21
  ! Maximum number of levels for input files
22
  integer :: nlevmax
23
  parameter (nlevmax=100)
24
 
25
  ! Maximum number of input files (dates, length of trajectories)
26
  integer :: ndatmax
27
  parameter (ndatmax=500)
28
 
29
  ! Numerical epsilon (for float comparison)
30
  real :: eps
31
  parameter (eps=0.001)
32
 
33
  ! Conversion factors
34
  real :: pi180                                   ! deg -> rad
35
  parameter (pi180=3.14159/180.)
36
  real :: deg2km                                  ! deg -> km (at equator)
37
  parameter (deg2km=111.)
38
  real :: pir
5 michaesp 39
  parameter (pir=255032235.95489)                 ! 2*Pi*R^2 Bojan
3 michaesp 40
 
41
  ! Prefix for primary and secondary fields
42
  character :: charp
43
  character :: chars
44
  parameter (charp='P')
45
  parameter (chars='S')
46
 
47
  ! --------------------------------------------------------------------
48
  ! Declaration of variables
49
  ! --------------------------------------------------------------------
50
 
51
  ! Input and output format for trajectories (see iotra.f)
52
  integer :: inpmode
53
  integer :: outmode
54
 
55
  ! Input parameters
56
  character(len=80) :: inpfile         ! Input trajectory file
57
  character(len=80) :: outfile         ! Output trajectory file
58
  integer :: ntra                      ! Number of trajectories
59
  integer :: ncol                      ! Number of columns (including time, lon, lat, p)
60
  integer :: ntim                      ! Number of times per trajectory
61
  integer :: ntrace0                   ! Number of trace variables
62
  character(len=80) :: tvar0(200)      ! Tracing variable (with mode specification)
63
  character(len=80) :: tvar(200)       ! Tracing variable name (only the variable)
15 michaesp 64
  character(len=80) :: tfil(200)       ! Filename prefix
3 michaesp 65
  real :: fac(200)                     ! Scaling factor
66
  real :: shift_val(200)               ! Shift in space and time relative to trajectory position
67
  character(len=80) :: shift_dir(200)  ! Direction of shift
39 michaesp 68
  character(len=80) :: shift_rel(200)  ! Operator/relator for variable
3 michaesp 69
  integer :: compfl(200)               ! Computation flag (1=compute)
70
  integer :: numdat                    ! Number of input files
71
  character(len=11) :: dat(ndatmax)    ! Dates of input files
72
  real :: timeinc                      ! Time increment between input files
73
  real :: tst                          ! Time shift of start relative to first data file
74
  real :: ten                          ! Time shift of end relatiev to first data file
75
  character(len=20) :: startdate       ! First time/date on trajectory
76
  character(len=20) :: enddate         ! Last time/date on trajectory
77
  integer :: ntrace1                   ! Count trace and additional variables
78
  character(len=80) :: timecheck       ! Either 'yes' or 'no'
79
  character(len=80) :: intmode         ! Interpolation mode ('normal', 'nearest') Bojan ('clustering','circle_avg','circle_max','circle_min')
80
 
81
  ! Trajectories
82
  real,allocatable, dimension (:,:,:) :: trainp          ! Input trajectories (ntra,ntim,ncol)
83
  real,allocatable, dimension (:,:,:) :: traint          ! Internal trajectories (ntra,ntim,ncol+ntrace1)
84
  real,allocatable, dimension (:,:,:) :: traout          ! Output trajectories (ntra,ntim,ncol+ntrace0)
85
  integer :: reftime(6)                                  ! Reference date
11 michaesp 86
  character(len=80) :: varsinp(200)                      ! Field names for input trajectory
87
  character(len=80) :: varsint(200)                      ! Field names for internal trajectory
88
  character(len=80) :: varsout(200)                      ! Field names for output trajectory
3 michaesp 89
  integer :: fid,fod                                     ! File identifier for inp and out trajectories
90
  real :: x0,y0,p0                                       ! Position of air parcel (physical space)
91
  real :: reltpos0                                       ! Relative time of air parcel
92
  real :: xind,yind,pind                                 ! Position of air parcel (grid space)
93
  integer :: fbflag                                      ! Flag for forward (1) or backward (-1) trajectories
11 michaesp 94
  integer :: fok(200)                                    ! Flag whether field is ready
3 michaesp 95
 
96
  ! Meteorological fields
97
  real,allocatable, dimension (:)     :: spt0,spt1       ! Surface pressure
98
  real,allocatable, dimension (:)     :: p3t0,p3t1       ! 3d-pressure
99
  real,allocatable, dimension (:)     :: f3t0,f3t1       ! 3d field for tracing
100
  character(len=80) :: svars(100)                        ! List of variables on S file
101
  character(len=80) :: pvars(100)                        ! List of variables on P file
102
  integer :: n_svars                                     ! Number of variables on S file
103
  integer :: n_pvars                                     ! Number of variables on P file
104
 
105
  ! Grid description
106
  real :: pollon,pollat   ! Longitude/latitude of pole
107
  real :: ak(100)         ! Vertical layers and levels
108
  real :: bk(100)
109
  real :: xmin,xmax       ! Zonal grid extension
110
  real :: ymin,ymax       ! Meridional grid extension
111
  integer :: nx,ny,nz     ! Grid dimensions
112
  real :: dx,dy           ! Horizontal grid resolution
113
  integer :: hem          ! Flag for hemispheric domain
114
  integer :: per          ! Flag for periodic domain
115
  real :: stagz           ! Vertical staggering
116
  real :: mdv             ! Missing data value
117
 
118
  ! Auxiliary variables
119
  integer :: i,j,k,l,m,n
120
  real :: rd
121
  character(len=80) :: filename,varname
122
  real :: time0,time1,reltpos
123
  integer :: itime0,itime1
124
  integer :: stat
125
  real :: tstart
126
  integer :: iloaded0,iloaded1
127
  real :: f0
128
  real :: frac
129
  real :: tload,tfrac
130
  integer :: isok
131
  character :: ch
132
  integer :: ind
133
  integer :: ind1,ind2,ind3,ind4,ind5
134
  integer :: ind6,ind7,ind8,ind9,ind0
135
  integer :: noutside
136
  real :: delta
137
  integer :: itrace0
138
  character(len=80) :: string
5 michaesp 139
  integer err_c1,err_c2,err_c3
3 michaesp 140
 
141
  ! Bojan
142
  real    :: dist,circlesum,circlemax,circlemin,circleavg,radius       ! distance (great circle), sum/max/min/avg in circle, radius of circle
143
  integer :: ist,jst,kst,sp,ml,mr,nd,nu                                ! ijk in stack, sp=stack counter, ml (left), mr (right), nd (down), nu (up)
144
  integer :: lci,lcj,xindb,xindf                                       ! label count i and j, xind back, xind forward
145
  integer :: yindb,yindf,pindb,pindf                                   ! yind back, yind forward, pind back, pind forward
146
  integer :: pvpos,thpos                                               ! position of variables PV and TH in trajectory
5 michaesp 147
  real    :: tropo_pv,tropo_th                                         ! values of PV and TH at dynamical tropopause
3 michaesp 148
  integer, allocatable, dimension (:)   :: stackx,stacky               ! lon/lat of stack
149
  integer, allocatable, dimension (:)   :: lblcount                    ! counter for label
150
  integer, allocatable, dimension (:,:) :: connect                     ! array that keeps track of the visited grid points
151
  real, allocatable, dimension (:)      :: lbl                         ! label
152
  real, allocatable, dimension (:)      :: circlelon,circlelat         ! value of f, lon and lat in circle
153
  real, allocatable, dimension (:)      :: circlef,circlearea          ! value of f, lon and lat in circle
154
  real, allocatable, dimension (:)      :: longrid,latgrid             ! arrays of longitude and latitude of grid points
155
  ! Bojan
156
 
157
  ! Externals
158
  real :: int_index4
5 michaesp 159
  external                               int_index4
160
  real :: sdis ! Bojan: need function sdis (calculates great circle distance)
161
  external                               sdis ! Bojan: need function sdis
3 michaesp 162
 
163
  ! --------------------------------------------------------------------
164
  ! Start of program, Read parameters, get grid parameters
165
  ! --------------------------------------------------------------------
166
 
167
  ! Write start message
168
  print*,'========================================================='
169
  print*,'              *** START OF PROGRAM TRACE ***'
170
  print*
171
 
172
  ! Read parameters
173
  open(10,file='trace.param')
174
   read(10,*) inpfile
175
   read(10,*) outfile
176
   read(10,*) startdate
177
   read(10,*) enddate
178
   read(10,*) fbflag
179
   read(10,*) numdat
180
   if ( fbflag.eq.1) then
181
      do i=1,numdat
182
         read(10,'(a11)') dat(i)
183
      enddo
184
   else
185
      do i=numdat,1,-1
186
         read(10,'(a11)') dat(i)
187
      enddo
188
   endif
189
   read(10,*) timeinc
190
   read(10,*) tst
191
   read(10,*) ten
192
   read(10,*) ntra
193
   read(10,*) ntim
194
   read(10,*) ncol
195
   read(10,*) ntrace0
196
   do i=1,ntrace0
197
      read(10,*) tvar(i), fac(i), compfl(i), tfil(i)
198
   enddo
199
   read(10,*) n_pvars
200
   do i=1,n_pvars
201
      read(10,*) pvars(i)
202
   enddo
203
   read(10,*) n_svars
204
   do i=1,n_svars
205
      read(10,*) svars(i)
206
   enddo
207
   read(10,*) timecheck
208
   read(10,*) intmode
209
   read(10,*) radius
5 michaesp 210
   read(10,*) tropo_pv
211
   read(10,*) tropo_th
3 michaesp 212
  close(10)
213
 
214
  ! Bojan: error if radius < 0
5 michaesp 215
  if (((intmode .eq. "circle_avg") .or. (intmode .eq. "circle_min") .or. (intmode .eq. "circle_max")) .and. (radius .lt. 0)) then
3 michaesp 216
     print*,'ERROR (circle): radius < 0!'
217
     stop
218
  endif
219
 
220
  ! Remove commented tracing fields
221
  itrace0 = 1
222
  do while ( itrace0.le.ntrace0)
223
     string = tvar(itrace0)
224
     if ( string(1:1).eq.'#' ) then
225
        do i=itrace0,ntrace0-1
226
           tvar(i)   = tvar(i+1)
227
           fac(i)    = fac(i+1)
228
           compfl(i) = compfl(i+1)
229
           tfil(i)   = tfil(i+1)
230
        enddo
231
        ntrace0 = ntrace0 - 1
232
     else
233
        itrace0 = itrace0 + 1
234
     endif
235
  enddo
236
 
237
  ! Save the tracing variable (including all mode specifications)
238
  do i=1,ntrace0
239
     tvar0(i) = tvar(i)
240
  enddo
241
 
242
  ! Set the formats of the input and output files
243
  call mode_tra(inpmode,inpfile)
244
  if (inpmode.eq.-1) inpmode=1
245
  call mode_tra(outmode,outfile)
246
  if (outmode.eq.-1) outmode=1
247
 
248
  ! Convert time shifts <tst,ten> from <hh.mm> into fractional time
249
  call hhmm2frac(tst,frac)
250
  tst = frac
251
  call hhmm2frac(ten,frac)
252
  ten = frac
253
 
254
  ! Set the time for the first data file (depending on forward/backward mode)
255
  if (fbflag.eq.1) then
256
    tstart = -tst
257
  else
258
    tstart = tst
259
  endif
260
 
261
  ! Read the constant grid parameters (nx,ny,nz,xmin,xmax,ymin,ymax,pollon,pollat)
262
  ! The negative <-fid> of the file identifier is used as a flag for parameter retrieval
46 michaesp 263
 filename = charp//dat(1)
264
  varname  = 'nil'
265
  do i=1,n_pvars
266
     if ( (pvars(i).eq.'U').and.(varname.eq.'nil') ) varname  = 'U'
267
     if ( (pvars(i).eq.'T').and.(varname.eq.'nil') ) varname  = 'T'
268
  enddo
269
  if ( varname.eq.'nil' ) varname = tvar(1)
3 michaesp 270
  call input_open (fid,filename)
271
  call input_grid (-fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny,tstart,pollon,pollat,rd,rd,nz,rd,rd,rd,timecheck)
272
  call input_close(fid)
46 michaesp 273
  if ( nz.eq.1 ) then
274
     print*,' ERROR: the first tracing variable must be 3D'
275
     stop
276
  endif
3 michaesp 277
 
46 michaesp 278
 
279
 ! filename = charp//dat(1)
280
 ! varname  = 'U'
281
 ! call input_open (fid,filename)
282
 ! call input_grid (-fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny,tstart,pollon,pollat,rd,rd,nz,rd,rd,rd,timecheck)
283
 ! call input_close(fid)
284
 
3 michaesp 285
  ! Allocate memory for some meteorological arrays
286
  allocate(spt0(nx*ny),stat=stat)
287
  if (stat.ne.0) print*,'*** error allocating array spt0 ***'   ! Surface pressure
288
  allocate(spt1(nx*ny),stat=stat)
289
  if (stat.ne.0) print*,'*** error allocating array spt1 ***'
290
  allocate(p3t0(nx*ny*nz),stat=stat)
291
  if (stat.ne.0) print*,'*** error allocating array p3t0 ***'   ! Pressure
292
  allocate(p3t1(nx*ny*nz),stat=stat)
293
  if (stat.ne.0) print*,'*** error allocating array p3t1 ***'
294
  allocate(f3t0(nx*ny*nz),stat=stat)
295
  if (stat.ne.0) print*,'*** error allocating array p3t0 ***'   ! Tracing field
296
  allocate(f3t1(nx*ny*nz),stat=stat)
297
  if (stat.ne.0) print*,'*** error allocating array p3t1 ***'
298
 
299
  ! Get memory for trajectory arrays
300
  allocate(trainp(ntra,ntim,ncol),stat=stat)
301
  if (stat.ne.0) print*,'*** error allocating array tra      ***'
302
 
303
  ! Bojan
304
  ! allocate memory for clustering mode 
305
  if (intmode .eq. 'clustering') then
306
     allocate(lbl(8),stat=stat)
307
     if (stat.ne.0) print*,'*** error allocating array lbl      ***'
308
     allocate(lblcount(5),stat=stat)
309
     if (stat.ne.0) print*,'*** error allocating array lblcount ***'
310
  endif
311
  ! allocate memory for circle mode
312
  if ( (intmode.eq.'circle_avg') .or. (intmode.eq.'circle_min') .or. (intmode.eq.'circle_max') ) then
313
     allocate(connect(nx,ny),stat=stat)
314
     if (stat.ne.0) print*,'*** error allocating connect ***'
315
     allocate(stackx(nx*ny),stat=stat)
316
     if (stat.ne.0) print*,'*** error allocating stackx ***'
317
     allocate(stacky(nx*ny),stat=stat)
318
     if (stat.ne.0) print*,'*** error allocating stacky ***'
319
     allocate(circlelon(nx*ny),stat=stat)
320
     if (stat.ne.0) print*,'*** error allocating circlelon ***'
321
     allocate(circlelat(nx*ny),stat=stat)
322
     if (stat.ne.0) print*,'*** error allocating circlelat ***'
323
     allocate(circlef(nx*ny),stat=stat)
324
     if (stat.ne.0) print*,'*** error allocating circlef ***'
325
     allocate(circlearea(nx*ny),stat=stat)
326
     if (stat.ne.0) print*,'*** error allocating circlearea ***'
327
     allocate(longrid(nx),stat=stat)
328
     if (stat.ne.0) print*,'*** error allocating longrid ***'
329
     allocate(latgrid(ny),stat=stat)
330
     if (stat.ne.0) print*,'*** error allocating latgrid ***'
331
     do m=1,nx
332
        longrid(m)=xmin+dx*(m-1)
333
     enddo
334
     do n=1,ny
335
        latgrid(n)=ymin+dy*(n-1)
336
     enddo
337
  endif
338
  ! Bojan
339
 
340
  ! Set the flags for periodic domains
341
  if ( abs(xmax-xmin-360.).lt.eps ) then
342
     per = 1
343
  elseif ( abs(xmax-xmin-360.+dx).lt.eps ) then
344
     per = 2
345
  else
346
     per = 0
347
  endif
348
 
349
  ! Set logical flag for periodic data set (hemispheric or not)
350
  hem = 0
351
  if (per.eq.0.) then
352
     delta=xmax-xmin-360.
353
     if (abs(delta+dx).lt.eps) then               ! Program aborts: arrays must be closed
354
        print*,' ERROR: arrays must be closed... Stop'
355
     else if (abs(delta).lt.eps) then ! Periodic and hemispheric
356
       hem=1
357
       per=360.
358
    endif
359
  else                                            ! Periodic and hemispheric
360
     hem=1
361
  endif
362
 
363
  ! Write some status information
364
  print*,'---- INPUT PARAMETERS -----------------------------------'
365
  print*
366
  print*,'  Input trajectory file  : ',trim(inpfile)
367
  print*,'  Output trajectory file : ',trim(outfile)
368
  print*,'  Format of input file   : ',inpmode
369
  print*,'  Format of output file  : ',outmode
370
  print*,'  Forward/backward       : ',fbflag
371
  print*,'  #tra                   : ',ntra
372
  print*,'  #col                   : ',ncol
373
  print*,'  #tim                   : ',ntim
374
  print*,'  No time check          : ',trim(timecheck)
375
  print*,'  Interpolation mode     : ',trim(intmode)
5 michaesp 376
  ! Bojan
377
  if (trim(intmode) .eq. "clustering") then
378
     print*,'  Tropopause PV [pvu]    : ',tropo_pv
379
     print*,'  Tropopause TH [K]      : ',tropo_th
380
  endif
3 michaesp 381
  do i=1,ntrace0
382
     if (compfl(i).eq.0) then
383
        print*,'  Tracing field          : ', trim(tvar(i)), fac(i), ' 0 ', tfil(i)
384
     else
385
        print*,'  Tracing field          : ', trim(tvar(i)), fac(i), '  1 ', tfil(i)
386
     endif
387
  enddo
388
  print*
389
  print*,'---- INPUT DATA FILES -----------------------------------'
390
  print*
391
  call frac2hhmm(tstart,tload)
392
  print*,'  Time of 1st data file  : ',tload
393
  print*,'  #input files           : ',numdat
394
  print*,'  time increment         : ',timeinc
395
  call frac2hhmm(tst,tload)
396
  print*,'  Shift of start         : ',tload
397
  call frac2hhmm(ten,tload)
398
  print*,'  Shift of end           : ',tload
399
  print*,'  First/last input file  : ',trim(dat(1)), ' ... ', trim(dat(numdat))
400
  print*,'  Primary variables      : ',trim(pvars(1))
401
  do i=2,n_pvars
402
     print*,'                         : ',trim(pvars(i))
403
  enddo
404
  if ( n_svars.ge.1 ) then
405
     print*,'  Secondary variables    : ',trim(svars(1))
406
     do i=2,n_svars
407
        print*,'                         : ',trim(svars(i))
408
     enddo
409
  endif
410
  print*
411
  print*,'---- CONSTANT GRID PARAMETERS ---------------------------'
412
  print*
413
  print*,'  xmin,xmax     : ',xmin,xmax
414
  print*,'  ymin,ymax     : ',ymin,ymax
415
  print*,'  dx,dy         : ',dx,dy
416
  print*,'  pollon,pollat : ',pollon,pollat
417
  print*,'  nx,ny,nz      : ',nx,ny,nz
418
  print*,'  per, hem      : ',per,hem
419
  print*
420
 
421
  ! --------------------------------------------------------------------
422
  ! Load the input trajectories
423
  ! --------------------------------------------------------------------
424
 
425
  ! Read the input trajectory file
426
  call ropen_tra(fid,inpfile,ntra,ntim,ncol,reftime,varsinp,inpmode)
427
  call read_tra (fid,trainp,ntra,ntim,ncol,inpmode)
428
  call close_tra(fid,inpmode)
429
 
430
  ! Check that first four columns correspond to time,lon,lat,p
431
  if ( (varsinp(1).ne.'time' ).or. &
432
       (varsinp(2).ne.'xpos' ).and.(varsinp(2).ne.'lon' ).or. &
433
       (varsinp(3).ne.'ypos' ).and.(varsinp(3).ne.'lat' ).or. &
434
       (varsinp(4).ne.'ppos' ).and.(varsinp(4).ne.'p'   ) ) then
435
     print*,' ERROR: problem with input trajectories ...'
436
     stop
437
  endif
438
  varsinp(1) = 'time'
439
  varsinp(2) = 'lon'
440
  varsinp(3) = 'lat'
441
  varsinp(4) = 'p'
442
 
443
  ! Write some status information of the input trajectories
444
  print*,'---- INPUT TRAJECTORIES ---------------------------------'
445
  print*
446
  print*,' Start date             : ',trim(startdate)
447
  print*,' End date               : ',trim(enddate)
448
  print*,' Reference time (year)  : ',reftime(1)
449
  print*,'                (month) : ',reftime(2)
450
  print*,'                (day)   : ',reftime(3)
451
  print*,'                (hour)  : ',reftime(4)
452
  print*,'                (min)   : ',reftime(5)
453
  print*,' Time range (min)       : ',reftime(6)
454
  do i=1,ncol
455
     print*,' Var                    :',i,trim(varsinp(i))
456
  enddo
457
  print*
458
 
459
  ! Check that first time is 0 - otherwise the tracing will produce
460
  ! wrong results because later in the code only absolute times are
461
  ! considered: <itime0   = int(abs(tfrac-tstart)/timeinc) + 1>. This
462
  ! will be changed in a future version.
463
  if ( abs( trainp(1,1,1) ).gt.eps ) then
464
     print*,' ERROR: First time of trajectory must be 0, i.e. '
465
     print*,'     correspond to the reference date. Otherwise'
466
     print*,'     the tracing will give wrong results... STOP'
467
     stop
468
  endif
469
 
470
  ! --------------------------------------------------------------------
471
  ! Check dependencies for trace fields which must be calculated
472
  ! --------------------------------------------------------------------
473
 
474
  ! Set the counter for extra fields
475
  ntrace1 = ntrace0
476
 
477
  ! Loop over all tracing variables
478
  i = 1
479
  do while (i.le.ntrace1)
480
 
481
     ! Skip fields which must be available on the input files
482
     if (i.le.ntrace0) then
483
        if (compfl(i).eq.0) goto 100
484
     endif
485
 
486
     ! Get the dependencies for potential temperature (TH)
487
     if ( tvar(i).eq.'TH' ) then
488
        varname='P'                                   ! P
489
        call add2list(varname,tvar,ntrace1)
490
        varname='T'                                   ! T
491
        call add2list(varname,tvar,ntrace1)
492
 
493
     ! Get the dependencies for potential temperature (TH)
494
     elseif ( tvar(i).eq.'RHO' ) then
495
        varname='P'                                   ! P
496
        call add2list(varname,tvar,ntrace1)
497
        varname='T'                                   ! T
498
        call add2list(varname,tvar,ntrace1)
499
 
500
     ! Get the dependencies for relative humidity (RH)
501
     elseif ( tvar(i).eq.'RH' ) then
502
        varname='P'                                   ! P
503
        call add2list(varname,tvar,ntrace1)
504
        varname='T'                                   ! T
505
        call add2list(varname,tvar,ntrace1)
506
        varname='Q'                                   ! Q
507
        call add2list(varname,tvar,ntrace1)
508
 
509
     ! Get the dependencies for equivalent potential temperature (THE)
510
     elseif ( tvar(i).eq.'THE' ) then
511
        varname='P'                                   ! P
512
        call add2list(varname,tvar,ntrace1)
513
        varname='T'                                   ! T
514
        call add2list(varname,tvar,ntrace1)
515
        varname='Q'                                   ! Q
516
        call add2list(varname,tvar,ntrace1)
517
 
518
     ! Get the dependencies for latent heating rate (LHR)
519
     elseif ( tvar(i).eq.'LHR' ) then
520
        varname='P'                                   ! P
521
        call add2list(varname,tvar,ntrace1)
522
        varname='T'                                   ! T
523
        call add2list(varname,tvar,ntrace1)
524
        varname='Q'                                   ! Q
525
        call add2list(varname,tvar,ntrace1)
526
        varname='OMEGA'                               ! OMEGA
527
        call add2list(varname,tvar,ntrace1)
528
        varname='RH'                                  ! RH
529
        call add2list(varname,tvar,ntrace1)
530
 
531
     ! Get the dependencies for wind speed (VEL)
532
     elseif ( tvar(i).eq.'VEL' ) then
533
        varname='U'                                   ! U
534
        call add2list(varname,tvar,ntrace1)
535
        varname='V'                                   ! V
536
        call add2list(varname,tvar,ntrace1)
537
 
538
     ! Get the dependencies for wind direction (DIR)
539
     elseif ( tvar(i).eq.'DIR' ) then
540
        varname='U'                                   ! U
541
        call add2list(varname,tvar,ntrace1)
542
        varname='V'                                   ! V
543
        call add2list(varname,tvar,ntrace1)
544
 
545
     ! Get the dependencies for du/dx (DUDX)
546
     elseif ( tvar(i).eq.'DUDX' ) then
547
        varname='U:+1DLON'                            ! U:+1DLON
548
        call add2list(varname,tvar,ntrace1)
549
        varname='U:-1DLON'                            ! U:-1DLON
550
        call add2list(varname,tvar,ntrace1)
551
 
552
     ! Get the dependencies for dv(dx (DVDX)
553
     elseif ( tvar(i).eq.'DVDX' ) then
554
        varname='V:+1DLON'                            ! V:+1DLON
555
        call add2list(varname,tvar,ntrace1)
556
        varname='V:-1DLON'                            ! V:-1DLON
557
        call add2list(varname,tvar,ntrace1)
558
 
559
     ! Get the dependencies for du/dy (DUDY)
560
     elseif ( tvar(i).eq.'DUDY' ) then
561
        varname='U:+1DLAT'                            ! U:+1DLAT
562
        call add2list(varname,tvar,ntrace1)
563
        varname='U:-1DLAT'                            ! U:-1DLAT
564
        call add2list(varname,tvar,ntrace1)
565
 
566
     ! Get the dependencies for dv/dy (DVDY)
567
     elseif ( tvar(i).eq.'DVDY' ) then
568
        varname='V:+1DLAT'                            ! V:+1DLAT
569
        call add2list(varname,tvar,ntrace1)
570
        varname='V:-1DLAT'                            ! V:-1DLAT
571
        call add2list(varname,tvar,ntrace1)
572
 
573
     ! Get the dependencies for du/dp (DUDP)
574
     elseif ( tvar(i).eq.'DUDP' ) then
575
        varname='U:+1DP'                              ! U:+1DP
576
        call add2list(varname,tvar,ntrace1)
577
        varname='U:-1DP'                              ! U:-1DP
578
        call add2list(varname,tvar,ntrace1)
579
        varname='P:+1DP'                              ! P:+1DP
580
        call add2list(varname,tvar,ntrace1)
581
        varname='P:-1DP'                              ! P:-1DP
582
        call add2list(varname,tvar,ntrace1)
583
 
584
     ! Get the dependencies for dv/dp (DVDP)
585
     elseif ( tvar(i).eq.'DVDP' ) then
586
        varname='V:+1DP'                              ! V:+1DP
587
        call add2list(varname,tvar,ntrace1)
588
        varname='V:-1DP'                              ! V:-1DP
589
        call add2list(varname,tvar,ntrace1)
590
        varname='P:+1DP'                              ! P:+1DP
591
        call add2list(varname,tvar,ntrace1)
592
        varname='P:-1DP'                              ! P:-1DP
593
        call add2list(varname,tvar,ntrace1)
594
 
595
     ! Get the dependencies for dt/dx (DTDX)
596
     elseif ( tvar(i).eq.'DTDX' ) then
597
        varname='T:+1DLON'                            ! T:+1DLON
598
        call add2list(varname,tvar,ntrace1)
599
        varname='T:-1DLON'                            ! T:-1DLON
600
        call add2list(varname,tvar,ntrace1)
601
 
602
     ! Get the dependencies for dth/dy (DTHDY)
603
     elseif ( tvar(i).eq.'DTHDY' ) then
604
        varname='T:+1DLAT'                            ! T:+1DLON
605
        call add2list(varname,tvar,ntrace1)
606
        varname='T:-1DLAT'                            ! T:-1DLON
607
        call add2list(varname,tvar,ntrace1)
608
        varname='P'                                   ! P
609
        call add2list(varname,tvar,ntrace1)
610
 
611
     ! Get the dependencies for dth/dx (DTHDX)
612
     elseif ( tvar(i).eq.'DTHDX' ) then
613
        varname='T:+1DLON'                            ! T:+1DLON
614
        call add2list(varname,tvar,ntrace1)
615
        varname='T:-1DLON'                            ! T:-1DLON
616
        call add2list(varname,tvar,ntrace1)
617
        varname='P'                                   ! P
618
        call add2list(varname,tvar,ntrace1)
619
 
620
     ! Get the dependencies for dt/dy (DTDY)
621
     elseif ( tvar(i).eq.'DTDY' ) then
622
        varname='T:+1DLAT'                            ! T:+1DLON
623
        call add2list(varname,tvar,ntrace1)
624
        varname='T:-1DLAT'                            ! T:-1DLON
625
        call add2list(varname,tvar,ntrace1)
626
 
627
     ! Get the dependencies for dt/dp (DTDP)
628
     elseif ( tvar(i).eq.'DTDP' ) then
629
        varname='T:+1DP'                              ! T:+1DP
630
        call add2list(varname,tvar,ntrace1)
631
        varname='T:-1DP'                              ! T:-1DP
632
        call add2list(varname,tvar,ntrace1)
633
        varname='P:+1DP'                              ! P:+1DP
634
        call add2list(varname,tvar,ntrace1)
635
        varname='P:-1DP'                              ! P:-1DP
636
        call add2list(varname,tvar,ntrace1)
637
 
638
     ! Get the dependencies for dth/dp (DTHDP)
639
     elseif ( tvar(i).eq.'DTHDP' ) then
640
        varname='T:+1DP'                              ! T:+1DP
641
        call add2list(varname,tvar,ntrace1)
642
        varname='T:-1DP'                              ! T:-1DP
643
        call add2list(varname,tvar,ntrace1)
644
        varname='T'                                   ! T
645
        call add2list(varname,tvar,ntrace1)
646
        varname='P:+1DP'                              ! P:+1DP
647
        call add2list(varname,tvar,ntrace1)
648
        varname='P:-1DP'                              ! P:-1DP
649
        call add2list(varname,tvar,ntrace1)
650
        varname='P'                                   ! P
651
        call add2list(varname,tvar,ntrace1)
652
 
653
     ! Get the dependencies for squared Brunt Vaiäläa frequency (NSQ)
654
     elseif ( tvar(i).eq.'NSQ' ) then
655
        varname='DTHDP'                                ! DTHDP
656
        call add2list(varname,tvar,ntrace1)
657
        varname='TH'                                   ! TH
658
        call add2list(varname,tvar,ntrace1)
659
        varname='RHO'                                  ! RHO
660
        call add2list(varname,tvar,ntrace1)
661
 
662
     ! Get the dependencies for relative vorticity (RELVORT)
663
     elseif ( tvar(i).eq.'RELVORT' ) then
664
        varname='U'                                    ! U
665
        call add2list(varname,tvar,ntrace1)
666
        varname='DUDY'                                 ! DUDY
667
        call add2list(varname,tvar,ntrace1)
668
        varname='DVDX'                                 ! DVDX
669
        call add2list(varname,tvar,ntrace1)
670
 
671
     ! Get the dependencies for relative vorticity (ABSVORT)
672
     elseif ( tvar(i).eq.'ABSVORT' ) then
673
        varname='U'                                    ! U
674
        call add2list(varname,tvar,ntrace1)
675
        varname='DUDY'                                 ! DUDY
676
        call add2list(varname,tvar,ntrace1)
677
        varname='DVDX'                                 ! DVDX
678
        call add2list(varname,tvar,ntrace1)
679
 
680
     ! Get the dependencies for divergence (DIV)
681
     elseif ( tvar(i).eq.'DIV' ) then
682
        varname='V'                                    ! U
683
        call add2list(varname,tvar,ntrace1)
684
        varname='DUDX'                                 ! DUDX
685
        call add2list(varname,tvar,ntrace1)
686
        varname='DVDY'                                 ! DVDY
687
        call add2list(varname,tvar,ntrace1)
688
 
689
     ! Get the dependencies for deformation (DEF)
690
     elseif ( tvar(i).eq.'DEF' ) then
691
        call add2list(varname,tvar,ntrace1)
692
        varname='DUDX'                                 ! DUDX
693
        call add2list(varname,tvar,ntrace1)
694
        varname='DVDY'                                 ! DVDY
695
        call add2list(varname,tvar,ntrace1)
696
        varname='DUDY'                                 ! DUDY
697
        call add2list(varname,tvar,ntrace1)
698
        varname='DVDX'                                 ! DVDX
699
        call add2list(varname,tvar,ntrace1)
700
 
701
     ! Get the dependencies for potential vorticity (PV)
702
     elseif ( tvar(i).eq.'PV' ) then
703
        varname='ABSVORT'                              ! ABSVORT
704
        call add2list(varname,tvar,ntrace1)
705
        varname='DTHDP'                                ! DTHDP
706
        call add2list(varname,tvar,ntrace1)
707
        varname='DUDP'                                 ! DUDP
708
        call add2list(varname,tvar,ntrace1)
709
        varname='DVDP'                                 ! DVDP
710
        call add2list(varname,tvar,ntrace1)
711
        varname='DTHDX'                                ! DTHDX
712
        call add2list(varname,tvar,ntrace1)
713
        varname='DTHDY'                                ! DTHDY
714
        call add2list(varname,tvar,ntrace1)
715
 
716
     ! Get the dependencies for Richardson number (RI)
717
     elseif ( tvar(i).eq.'RI' ) then
718
        varname='DUDP'                                 ! DUDP
719
        call add2list(varname,tvar,ntrace1)
720
        varname='DVDP'                                 ! DVDP
721
        call add2list(varname,tvar,ntrace1)
722
        varname='NSQ'                                  ! NSQ
723
        call add2list(varname,tvar,ntrace1)
724
        varname='RHO'                                  ! RHO
725
        call add2list(varname,tvar,ntrace1)
726
 
727
     ! Get the dependencies for Ellrod&Knapp's turbulence index (TI)
728
     elseif ( tvar(i).eq.'TI' ) then
729
        varname='DEF'                                  ! DEF
730
        call add2list(varname,tvar,ntrace1)
731
        varname='DUDP'                                 ! DUDP
732
        call add2list(varname,tvar,ntrace1)
733
        varname='DVDP'                                 ! DVDP
734
        call add2list(varname,tvar,ntrace1)
735
        varname='RHO'                                  ! RHO
736
        call add2list(varname,tvar,ntrace1)
737
 
738
     endif
739
 
740
     ! Exit point for handling additional fields
741
 100   continue
742
     i = i + 1
743
 
744
  enddo
745
 
746
  ! Save the full variable name (including shift specification)
747
  do i=1,ncol
748
     varsint(i)      = varsinp(i)
749
  enddo
750
  do i=1,ntrace1
751
     varsint(i+ncol) = tvar(i)
752
  enddo
753
 
754
  ! Bojan: check that PV and TH are on trajectory
755
  if (intmode .eq. 'clustering') then
756
     pvpos=-1
757
     thpos=-1
758
     do i=1,ncol+ntrace1
759
        if (varsint(i) .eq. 'TH') then
760
        thpos=i
761
        print*,'Clustering: Found TH at position:',thpos
762
        endif
763
        if (varsint(i) .eq. 'PV') then
764
        pvpos=i
765
        print*,'Clustering: Found PV at position:',pvpos
766
        endif
767
     enddo
768
     if (thpos .eq. -1) then
769
        print*,'WARNING (clustering): Did not find TH'
770
        stop
771
     endif
772
     if (pvpos .eq. -1) then
773
        print*,'WARNING (clustering): Did not find PV'
774
        stop
775
     endif
776
  endif
777
  ! Bojan
778
 
5 michaesp 779
  ! Split the tracing variables
780
  do i=1,ntrace0
39 michaesp 781
    call splitvar(tvar(i),shift_val(i),shift_dir(i),shift_rel(i) )
5 michaesp 782
  enddo
783
 
784
 
3 michaesp 785
  ! Split the variable name and set flags
5 michaesp 786
  do i=ntrace0+1,ntrace1
3 michaesp 787
 
788
     ! Set the scaling factor
5 michaesp 789
     fac(i) = 1.
3 michaesp 790
 
791
     ! Set the base variable name, the shift and the direction
39 michaesp 792
     call splitvar(tvar(i),shift_val(i),shift_dir(i),shift_rel(i) )
3 michaesp 793
 
5 michaesp 794
     ! Set the prefix of the file name for additional fields
795
     tfil(i)='*'
796
     do j=1,n_pvars
797
        if ( tvar(i).eq.pvars(j) ) tfil(i)=charp
798
     enddo
799
     do j=1,n_svars
800
        if ( tvar(i).eq.svars(j) ) tfil(i)=chars
801
     enddo
3 michaesp 802
 
803
     ! Set the computational flag
804
     if ( (tvar(i).eq.'P'   ).or. &
805
          (tvar(i).eq.'PLAY').or. &
806
          (tvar(i).eq.'PLEV') ) then
807
        compfl(i) = 0
808
        tfil(i)   = charp
809
     elseif ( ( tfil(i).eq.charp ).or.( tfil(i).eq.chars ) ) then
810
        compfl(i) = 0
811
     else
812
        compfl(i) = 1
813
     endif
814
 
815
  enddo
816
 
817
  ! Check whether the shift modes are supported
818
  do i=1,ntrace1
819
     if ( ( shift_dir(i).ne.'nil'     ).and. &
820
          ( shift_dir(i).ne.'DLON'    ).and. &
821
          ( shift_dir(i).ne.'DLAT'    ).and. &
822
          ( shift_dir(i).ne.'DP'      ).and. &
823
          ( shift_dir(i).ne.'HPA'     ).and. &
15 michaesp 824
          ( shift_dir(i).ne.'HPA(ABS)').and. &
3 michaesp 825
          ( shift_dir(i).ne.'KM(LON)' ).and. &
826
          ( shift_dir(i).ne.'KM(LAT)' ).and. &
827
          ( shift_dir(i).ne.'H'       ).and. &
828
          ( shift_dir(i).ne.'MIN'     ).and. &
39 michaesp 829
          ( shift_dir(i).ne.'INDP'    ).and. & 
830
          ( shift_dir(i).ne.'PMIN'    ).and. &
831
          ( shift_dir(i).ne.'PMAX'    ) ) then
3 michaesp 832
        print*,' ERROR: shift mode ',trim(shift_dir(i)), ' not supported'
833
        stop
834
     endif
835
  enddo
836
 
837
  ! Write status information
838
  print*
839
  print*,'---- COMPLETE TABLE FOR TRACING -------------------------'
840
  print*
841
  do i=1,ntrace1
842
     if ( ( shift_dir(i).ne.'nil' ) ) then
39 michaesp 843
        write(*,'(i4,a4,a8,f10.2,a8,a8,3x,a4,i5)') i,' : ',trim(tvar(i)), &
844
             shift_val(i),trim(shift_dir(i)),trim(shift_rel(i)),tfil(i),compfl(i)
3 michaesp 845
     else
15 michaesp 846
        write(*,'(i4,a4,a8,10x,8x,3x,a4,i5)') &
3 michaesp 847
             i,' : ',trim(tvar(i)),tfil(i),compfl(i)
848
     endif
849
  enddo
850
 
851
  ! --------------------------------------------------------------------
852
  ! Prepare the internal and output trajectories
853
  ! --------------------------------------------------------------------
854
 
855
  ! Allocate memory for internal trajectories
856
  allocate(traint(ntra,ntim,ncol+ntrace1),stat=stat)
857
  if (stat.ne.0) print*,'*** error allocating array traint   ***'
858
 
859
  ! Copy input to output trajectory
860
  do i=1,ntra
861
     do j=1,ntim
862
        do k=1,ncol
863
           traint(i,j,k)=trainp(i,j,k)
864
        enddo
865
     enddo
866
  enddo
867
 
868
  ! Set the flags for ready fields/colums - at begin only read-in fields are ready
869
  do i=1,ncol
870
     fok(i) = 1
871
  enddo
872
  do i=ncol+1,ntrace1
873
     fok(i) = 0
874
  enddo
875
 
876
  ! --------------------------------------------------------------------
877
  ! Trace the fields (fields available on input files)
878
  ! --------------------------------------------------------------------
879
 
880
  print*
881
  print*,'---- TRACING FROM PRIMARY AND SECONDARY DATA FILES ------'
882
 
883
  ! Loop over all tracing fields
884
  do i=1,ntrace1
885
 
886
      ! Skip fields which must be computed (compfl=1), will be handled later
887
      if (compfl(i).ne.0)  goto 110
888
 
889
      ! Write some status information
890
      print*
891
      print*,' Now tracing             : ', trim(tvar(i)),shift_val(i),trim(shift_dir(i)),compfl(i),' ',trim(tfil(i))
892
 
893
      ! Set the flag for ready field/column
894
      fok(ncol+i) = 1
895
 
896
      ! Reset flags for load manager
897
      iloaded0 = -1
898
      iloaded1 = -1
899
 
900
      ! Reset the counter for fields outside domain
901
      noutside = 0
5 michaesp 902
      err_c1   = 0
903
      err_c2   = 0
904
      err_c3   = 0
3 michaesp 905
 
906
      ! Loop over all times
907
      do j=1,ntim
908
 
909
         ! Convert trajectory time from hh.mm to fractional time
910
         call hhmm2frac(trainp(1,j,1),tfrac)
911
 
912
         ! Shift time if requested
913
         if ( shift_dir(i).eq.'H' ) then
914
            tfrac = tfrac + shift_val(i)
915
         elseif ( shift_dir(i).eq.'MIN' ) then
916
            tfrac = tfrac + shift_val(i)/60.
917
         endif
918
 
919
         ! Get the times which are needed
920
         itime0   = int(abs(tfrac-tstart)/timeinc) + 1
921
         time0    = tstart + fbflag * real(itime0-1) * timeinc
922
         itime1   = itime0 + 1
923
         time1    = time0 + fbflag * timeinc
924
         if ( itime1.gt.numdat ) then
925
            itime1 = itime0
926
            time1  = time0
927
         endif
928
 
929
         ! Load manager: Check whether itime0 can be copied from itime1
930
         if ( itime0.eq.iloaded1 ) then
931
            f3t0     = f3t1
932
            p3t0     = p3t1
933
            spt0     = spt1
934
            iloaded0 = itime0
935
         endif
936
 
937
         ! Load manager: Check whether itime1 can be copied from itime0
938
         if ( itime1.eq.iloaded0 ) then
939
            f3t1     = f3t0
940
            p3t1     = p3t0
941
            spt1     = spt0
942
            iloaded1 = itime1
943
         endif
944
 
945
         ! Load manager:  Load first time (tracing variable and grid)
946
         if ( itime0.ne.iloaded0 ) then
947
 
15 michaesp 948
            filename = trim(tfil(i))//trim(dat(itime0))
3 michaesp 949
            call frac2hhmm(time0,tload)
950
            varname  = tvar(i)
951
            write(*,'(a23,a20,a3,a5,f7.2)') '    ->  loading          : ', trim(filename),'  ',trim(varname),tload
952
            call input_open (fid,filename)
39 michaesp 953
            mdv = -999.
3 michaesp 954
            call input_wind &
955
                 (fid,varname,f3t0,tload,stagz,mdv, &
956
                 xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
957
            call input_grid &
958
                 (fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny, &
959
                 tload,pollon,pollat,p3t0,spt0,nz,ak,bk,stagz, &
960
                 timecheck)
961
            call input_close(fid)
962
 
963
            iloaded0 = itime0
964
 
965
         endif
966
 
967
         ! Load manager: Load second time (tracing variable and grid)
968
         if ( itime1.ne.iloaded1 ) then
969
 
15 michaesp 970
            filename = trim(tfil(i))//trim(dat(itime1))
3 michaesp 971
            call frac2hhmm(time1,tload)
972
            varname  = tvar(i)
973
            write(*,'(a23,a20,a3,a5,f7.2)') '    ->  loading          : ', trim(filename),'  ',trim(varname),tload
974
            call input_open (fid,filename)
975
            call input_wind &
976
                 (fid,varname,f3t1,tload,stagz,mdv, &
977
                 xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz,timecheck)
978
            call input_grid &
979
                 (fid,varname,xmin,xmax,ymin,ymax,dx,dy,nx,ny, &
980
                 tload,pollon,pollat,p3t1,spt1,nz,ak,bk,stagz, &
981
                 timecheck)
982
            call input_close(fid)
983
 
984
            iloaded1 = itime1
985
 
986
         endif
987
 
988
         ! Loop over all trajectories
989
         do k=1,ntra
990
 
991
            ! Set the horizontal position where to interpolate to
992
            x0       = traint(k,j,2)                          ! Longitude
993
            y0       = traint(k,j,3)                          ! Latitude
994
 
995
            ! Set the vertical position where to interpolate to
996
            if ( nz.gt.1 ) then
997
               p0 = traint(k,j,4)                             ! Pressure (3D tracing)
998
            else
999
               p0 = 1050.                                     ! Lowest level (2D tracing)
1000
            endif
1001
 
1002
            ! Set negative pressures to mdv
5 michaesp 1003
            if (p0.lt.0.) then
1004
	        f0 = mdv
1005
                goto 109
1006
            endif
3 michaesp 1007
 
1008
            ! Set the relative time
1009
            call hhmm2frac(traint(k,j,1),tfrac)
1010
            reltpos0 = fbflag * (tfrac-time0)/timeinc
1011
 
1012
            ! Make adjustments depending on the shift flag
1013
            if ( shift_dir(i).eq.'DLON' ) then                         ! DLON
1014
               x0 = x0 + shift_val(i)
1015
 
1016
            elseif  ( shift_dir(i).eq.'DLAT' ) then                    ! DLAT
1017
               y0 = y0 + shift_val(i)
1018
 
1019
            elseif ( shift_dir(i).eq.'KM(LON)' ) then                  ! KM(LON)
1020
               x0 = x0 + shift_val(i)/deg2km * 1./cos(y0*pi180)
1021
 
1022
            elseif ( shift_dir(i).eq.'KM(LAT)' ) then                  ! KM(LAT)
1023
               y0 = y0 + shift_val(i)/deg2km
1024
 
1025
            elseif ( shift_dir(i).eq.'HPA' ) then                      ! HPA
1026
               p0 = p0 + shift_val(i)
1027
 
15 michaesp 1028
            elseif ( shift_dir(i).eq.'HPA(ABS)' ) then                 ! HPA(ABS)
1029
               p0 = shift_val(i)
1030
 
3 michaesp 1031
            elseif ( shift_dir(i).eq.'DP' ) then                       ! DP
1032
               call get_index4 (xind,yind,pind,x0,y0,p0,reltpos0, &
1033
                    p3t0,p3t1,spt0,spt1,3,nx,ny,nz,xmin,ymin,dx,dy,mdv)
1034
               pind = pind - shift_val(i)
1035
               p0   = int_index4(p3t0,p3t1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
1036
 
1037
            elseif ( shift_dir(i).eq.'INDP' ) then
1038
               p0   = int_index4(p3t0,p3t1,nx,ny,nz,xind,yind,shift_val(i),reltpos0,mdv)
1039
 
1040
            endif
1041
 
1042
            ! Handle periodic boundaries in zonal direction
1043
            if ( (x0.gt.xmax).and.(per.ne.0) ) x0 = x0 - 360.
1044
            if ( (x0.lt.xmin).and.(per.ne.0) ) x0 = x0 + 360.
1045
 
1046
            ! Handle pole problems for hemispheric data (taken from caltra.f)
1047
            if ((hem.eq.1).and.(y0.gt.90.)) then
1048
               print*,'WARNING: y0>90 ',y0,' => setting to 180-y0 ',180.-y0
1049
               y0=180.-y0
1050
               x0=x0+per/2.
1051
            endif
1052
            if ((hem.eq.1).and.(y0.lt.-90.)) then
5 michaesp 1053
               print*,'WARNING: y0<-90 ',y0,' => setting to -180-y0 ',-180.-y0
3 michaesp 1054
               y0=-180.-y0
1055
               x0=x0+per/2.
1056
            endif
1057
 
1058
            ! Get the index where to interpolate (x0,y0,p0)
1059
            if ( (abs(x0-mdv).gt.eps).and. (abs(y0-mdv).gt.eps) ) then
1060
               call get_index4 (xind,yind,pind,x0,y0,p0,reltpos0, &
1061
                    p3t0,p3t1,spt0,spt1,3,nx,ny,nz,xmin,ymin,dx,dy,mdv)
1062
            else
1063
               xind = mdv
1064
               yind = mdv
1065
               pind = mdv
1066
            endif
1067
 
5 michaesp 1068
           ! Check if point is within grid (keep indices if ok)
3 michaesp 1069
            if ( (xind.ge.1.).and.(xind.le.real(nx)).and. &
1070
                 (yind.ge.1.).and.(yind.le.real(ny)).and. &
1071
                 (pind.ge.1.).and.(pind.le.real(nz)) ) then
1072
		 xind = xind
1073
		 yind = yind
1074
		 pind = pind
1075
 
1076
            ! Check if pressure is outside, but rest okay => adjust to lowest or highest level
1077
            elseif ( (xind.ge.1.).and.(xind.le.real(nx)).and. (yind.ge.1.).and.(yind.le.real(ny)) ) then ! only vertical problem
5 michaesp 1078
 
3 michaesp 1079
               if ( pind.gt.nz ) then ! pressure too low, index too high
5 michaesp 1080
                 err_c1 = err_c1 + 1
1081
                 if ( err_c1.lt.10 ) then
39 michaesp 1082
                    write(*,'(a,f7.3,a)') ' WARNING: pressure too low (pind = ',pind,') => adjusted to highest level (pind=nz.)'
5 michaesp 1083
                    print*,'(x0,y0,p0)=',x0,y0,p0
1084
                    pind = real(nz)
1085
                 elseif ( err_c1.eq.10 ) then
1086
                    print*,' WARNING: more pressures too low -> adjusted to highest level '
1087
                    pind = real(nz)
1088
                 else
1089
                    pind = real(nz)
1090
                 endif
1091
 
3 michaesp 1092
               elseif (pind.lt.1.) then ! pressure too high, index too low
5 michaesp 1093
                 err_c2 = err_c2 + 1
1094
                 if ( err_c2.lt.10 ) then
39 michaesp 1095
                    write(*,'(a,f7.3,a)') ' WARNING: pressure too high (pind = ',pind,') => adjusted to lowest level, (pind=1.)'
5 michaesp 1096
                    print*,'(x0,y0,p0)=',x0,y0,p0
1097
                    pind = 1.
1098
                 elseif ( err_c2.eq.10 ) then
1099
                    print*,' WARNING: more pressures too high -> adjusted to lowest level '
1100
                    pind = 1.
1101
                 else
1102
                    pind = 1.
1103
                 endif
1104
 
3 michaesp 1105
              endif
1106
 
5 michaesp 1107
            ! Grid point is outside!
3 michaesp 1108
            else
5 michaesp 1109
 
1110
               err_c3 = err_c3 + 1
1111
               if ( err_c3.lt.10 ) then
1112
                  print*,'ERROR: point is outside grid (horizontally)'
1113
                  print*,'   Trajectory # ',k
1114
                  print*,'   Position     ',x0,y0,p0
1115
                  print*,'  (xind,yind):  ',xind,yind
1116
                  xind          = mdv
1117
                  yind          = mdv
1118
                  pind          = mdv
1119
                  traint(k,j,2) = mdv  
1120
                  traint(k,j,3) = mdv  
1121
                  traint(k,j,4) = mdv  
1122
               elseif ( err_c3.eq.10 ) then
1123
                  print*,'ERROR: more points outside grid (horizontally)'
1124
                  xind          = mdv
1125
                  yind          = mdv
1126
                  pind          = mdv
1127
                  traint(k,j,2) = mdv  
1128
                  traint(k,j,3) = mdv  
1129
                  traint(k,j,4) = mdv  
1130
               else
1131
                  xind          = mdv
1132
                  yind          = mdv
1133
                  pind          = mdv
1134
                  traint(k,j,2) = mdv  
1135
                  traint(k,j,3) = mdv  
1136
                  traint(k,j,4) = mdv  
1137
               endif
1138
 
3 michaesp 1139
            endif
1140
 
1141
            ! ------------------------ NEAREST mode ------------------------------- 
5 michaesp 1142
            ! Interpolate to nearest grid point
3 michaesp 1143
            if ( intmode.eq.'nearest') then
1144
 
1145
               xind = real( nint(xind) )
1146
               yind = real( nint(yind) )
1147
               pind = real( nint(pind) )
1148
 
1149
               if ( xind.lt.1.  ) xind = 1.
1150
               if ( xind.gt.nx  ) xind = real(nx)
1151
               if ( yind.lt.1.  ) yind = 1.
1152
               if ( yind.gt.ny  ) yind = real(ny)
1153
 
1154
               if ( pind.lt.1.  ) pind = 1.
1155
               if ( pind.gt.nz  ) pind = real(nz)
1156
 
1157
               if ( abs(reltpos0).ge.eps ) then
1158
                  print*,'ERROR (nearest): reltpos != 0',reltpos0
1159
                  stop
1160
               endif
1161
 
1162
               ! interpolate
1163
               f0 = int_index4(f3t0,f3t1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
1164
 
5 michaesp 1165
            ! ------------------------ end NEAREST mode ------------------------------- 
1166
 
3 michaesp 1167
            ! ------------------------ CLUSTERING mode ------------------------------- Bojan
1168
            elseif (intmode.eq.'clustering') then
5 michaesp 1169
               if (varname.ne.'LABEL' ) then
1170
                  print*,'ERROR (clustering): varname is not LABEL'
1171
                  stop
1172
               endif
3 michaesp 1173
 
1174
            ! Get indices of box around the point
5 michaesp 1175
            xindb=floor(xind)
1176
            xindf=ceiling(xind)
1177
            yindb=floor(yind)
1178
            yindf=ceiling(yind)
1179
            pindb=floor(pind)
1180
            pindf=ceiling(pind)
3 michaesp 1181
 
1182
            ! Make sure all points are within grid
5 michaesp 1183
            if ( xindb.lt.1 ) xindb = 1
1184
            if ( xindf.lt.1 ) xindf = 1
3 michaesp 1185
            if ( xindb.gt.nx ) xindb = nx
1186
            if ( xindf.gt.nx ) xindf = nx
5 michaesp 1187
            if ( yindb.lt.1 ) yindb = 1
1188
            if ( yindf.lt.1 ) yindf = 1
3 michaesp 1189
            if ( yindb.gt.ny ) yindb = ny
1190
            if ( yindf.gt.ny ) yindf = ny
5 michaesp 1191
            if ( pindb.lt.1 ) pindb = 1
1192
            if ( pindf.lt.1 ) pindf = 1
3 michaesp 1193
            if ( pindb.gt.nz ) pindb = nz
1194
            if ( pindf.gt.nz ) pindf = nz
1195
 
1196
            ! Shift one point if they are equal
1197
            if ( xindb.eq.xindf ) then
1198
               if ( xindf.eq.nx ) then
1199
                  xindb=nx-1
1200
               else
1201
                  xindf=xindb+1
1202
               endif
1203
            endif
1204
            if ( yindb.eq.yindf ) then
1205
               if ( yindf.eq.ny ) then
1206
                  yindb=ny-1
1207
               else
1208
                  yindf=yindb+1
1209
               endif
1210
            endif
1211
            if ( pindb.eq.pindf ) then
1212
               if ( pindf.eq.nz ) then
1213
                  pindb=nz-1
1214
               else
1215
                  pindf=pindb+1
1216
               endif
1217
            endif
1218
            ! Give warnings and stop if problems occur
1219
            if ( xindb.eq.xindf ) then
1220
               print*,'ERROR (clustering): xindb=xindf'
1221
               print*,xind,xindb,xindf
1222
               stop
1223
            endif
1224
            if ( yindb.eq.yindf ) then
1225
               print*,'ERROR (clustering): yindb=yindf'
1226
               print*,yind,yindb,yindf
1227
               stop
1228
            endif
1229
            if ( pindb.eq.pindf ) then
1230
               print*,'ERROR (clustering): pindb=pindf'
1231
               print*,pind,pindb,pindf
1232
               stop
1233
            endif
1234
            if ( ( xindb.lt.1 ).or.( xindf.gt.nx ) ) then
1235
               print*,'ERROR (clustering): xindb/f outside'
1236
               print*,xind,xindb,xindf
1237
               stop
1238
            endif
1239
            if ( ( yindb.lt.1 ).or.( yindf.gt.ny ) ) then
1240
               print*,'ERROR (clustering): yindb/f outside'
1241
               print*,yind,yindb,yindf
1242
               stop
1243
            endif
1244
            if ( ( pindb.lt.1 ).or.( pindf.gt.nz ) ) then
1245
               print*,'ERROR (clustering): pindb/f outside'
1246
               print*,pind,pindb,pindf
1247
               stop
1248
            endif
1249
            if ( abs(reltpos0).ge.eps ) then
1250
               print*,'ERROR (clustering): reltpos != 0',reltpos0
1251
               stop
1252
            endif
1253
 
1254
            ! Get Value in Box
1255
            lblcount=(/0,0,0,0,0/)
1256
 
1257
            lbl(1) = f3t0( xindb + nx*(yindb-1) + nx*ny*(pindb-1) )
1258
            lbl(2) = f3t0( xindf + nx*(yindb-1) + nx*ny*(pindb-1) )
1259
            lbl(3) = f3t0( xindb + nx*(yindf-1) + nx*ny*(pindb-1) )
1260
            lbl(4) = f3t0( xindf + nx*(yindf-1) + nx*ny*(pindb-1) )
1261
            lbl(5) = f3t0( xindb + nx*(yindb-1) + nx*ny*(pindf-1) )
1262
            lbl(6) = f3t0( xindf + nx*(yindb-1) + nx*ny*(pindf-1) )
1263
            lbl(7) = f3t0( xindb + nx*(yindf-1) + nx*ny*(pindf-1) )
1264
            lbl(8) = f3t0( xindf + nx*(yindf-1) + nx*ny*(pindf-1) )
1265
 
1266
            ! Count the number of times every label appears
1267
            do lci=1,5
1268
               do lcj=1,8
1269
                  if ( abs(lbl(lcj)-lci).lt.eps ) then
1270
                     lblcount(lci)=lblcount(lci)+1
1271
                  endif
1272
               enddo
1273
            enddo
1274
 
1275
            ! Set to -9 to detect if no label was assigned in the end
1276
            f0=-9
1277
 
1278
            ! Stratosphere (PV)
5 michaesp 1279
            if ( abs(traint(k,j,pvpos)) .ge. tropo_pv ) then
3 michaesp 1280
               if ( (lblcount(2).ge.lblcount(3)).and. (lblcount(2).ge.lblcount(5)) ) then
1281
                  f0=2
1282
               elseif ( lblcount(3).ge.lblcount(5) ) then
1283
                  f0=3
1284
               elseif ( lblcount(5).gt.lblcount(3) ) then
1285
                  f0=5
1286
               endif
1287
            endif
1288
 
1289
            ! Troposphere (PV)
5 michaesp 1290
            if ( abs(traint(k,j,pvpos)) .lt. tropo_pv ) then
3 michaesp 1291
               if ( lblcount(1).ge.lblcount(4) ) then
1292
               f0=1
1293
               elseif ( lblcount(4).gt.lblcount(1) ) then
1294
               f0=4
1295
               endif
1296
            endif
1297
 
1298
            ! Stratosphere (TH)
5 michaesp 1299
            if ( traint(k,j,thpos) .ge. tropo_th ) then
3 michaesp 1300
               f0=2
1301
            endif
1302
 
1303
            if (f0.eq.-9) then
1304
               print*,'ERROR (Clustering): No label assigned!'
1305
               stop
1306
            endif
5 michaesp 1307
            ! ------------------------ end CLUSTERING mode -------------------------------
3 michaesp 1308
 
1309
            ! ------------------------ CIRCLE modes ------------------------------- Bojan
5 michaesp 1310
            ! elseif (not clustering but one of the possible circle modes)
3 michaesp 1311
            elseif ( (intmode.eq.'circle_avg') .or. (intmode.eq.'circle_min') .or. (intmode.eq.'circle_max') ) then
1312
 
1313
            ! reset arrays for this point
5 michaesp 1314
            connect=0
1315
            stackx=0
1316
            stacky=0
1317
            circlelon=0
1318
            circlelat=0
1319
            circlef=0
1320
            circlearea=0
3 michaesp 1321
 
1322
            ! Get indices of one coarse grid point within search radius (nint=round to next integer)
1323
            if ( sdis(x0,y0,longrid(nint(xind)),latgrid(nint(yind))) .gt. radius) then
1324
               print*,'ERROR (circle): Search radius is too small... (1). r =',radius
1325
               print*,'Distance to nint grid point (minimum search radius)=',sdis(x0,y0,longrid(nint(xind)),latgrid(nint(yind)))
1326
               stop
1327
            endif
1328
 
1329
            ! Initialize stack with nint(xind),nint(yind)
5 michaesp 1330
            kst=0 ! counts the number of points in circle
1331
            stackx(1)=nint(xind)
1332
            stacky(1)=nint(yind)
1333
            sp=1 ! stack counter
3 michaesp 1334
            do while (sp.ne.0)
1335
 
1336
            ! Get an element from stack
1337
             ist=stackx(sp)
1338
             jst=stacky(sp)
1339
             sp=sp-1
1340
 
1341
            ! Get distance from reference point
1342
             dist=sdis(x0,y0,longrid(ist),latgrid(jst))
1343
 
1344
            ! Check whether distance is smaller than search radius: connected
1345
             if (dist.lt.radius) then
1346
 
1347
            ! Increase total stack index
1348
              kst=kst+1
1349
              circlelon(kst)=longrid(ist)
1350
              circlelat(kst)=latgrid(jst)
1351
 
1352
            ! Interpolate field to position of point (interpolation in time!) 
1353
              circlef(kst) = int_index4(f3t0,f3t1,nx,ny,nz,real(ist),real(jst),pind,reltpos0,mdv)
1354
 
1355
            ! Calculate area of point (for circle_avg mode only)
1356
              if ( intmode .eq. 'circle_avg' ) then
1357
                 circlearea(kst) = pir/(nx-1)*(sin(pi180*abs(circlelat(kst)))-sin(pi180*(abs(circlelat(kst))-dy)))
1358
              endif
1359
 
1360
            ! Mark this point as visited
1361
              connect(ist,jst)=1
1362
 
1363
            ! Get coordinates of neighbouring points and implement periodicity
1364
              mr=ist+1
1365
              if (mr.gt.nx) mr=1
1366
              ml=ist-1
1367
              if (ml.lt.1) ml=nx
1368
              nu=jst+1
1369
              if (nu.gt.ny) nu=ny
1370
              nd=jst-1
1371
              if (nd.lt.1) nd=1
1372
 
1373
            ! Update stack with neighbouring points
1374
              if (connect(mr,jst).ne. 1) then
1375
                 connect(mr,jst)=1
1376
                 sp=sp+1
1377
                 stackx(sp)=mr
1378
                 stacky(sp)=jst
1379
              endif
1380
              if (connect(ml,jst).ne. 1) then
1381
                 connect(ml,jst)=1
1382
                 sp=sp+1
1383
                 stackx(sp)=ml
1384
                 stacky(sp)=jst
1385
              endif
1386
              if (connect(ist,nd).ne. 1) then
1387
                 connect(ist,nd)=1
1388
                 sp=sp+1
1389
                 stackx(sp)=ist
1390
                 stacky(sp)=nd
1391
              endif
1392
              if (connect(ist,nu).ne. 1) then
1393
                 connect(ist,nu)=1
1394
                 sp=sp+1
1395
                 stackx(sp)=ist
1396
                 stacky(sp)=nu
1397
              endif
1398
 
1399
             endif ! endif radius is smaller => end of updating stack
1400
 
1401
            end do ! end working on stack 
1402
 
1403
            if (kst.ge.1) then
5 michaesp 1404
               ! Choose output depending on intmode
3 michaesp 1405
               if ( intmode .eq. 'circle_avg' ) then
5 michaesp 1406
                  ! calculate area-weighted average of f in circle
3 michaesp 1407
                  circlesum=0.
1408
                  do l=1,kst
1409
                     circlesum=circlesum+circlef(l)*circlearea(l)
1410
                  enddo
1411
                  circleavg=circlesum/sum(circlearea(1:kst))
5 michaesp 1412
                  !print*,'area-weighted average of f in circle=',circleavg
3 michaesp 1413
                  f0=circleavg
1414
               elseif ( intmode .eq. 'circle_min' ) then
5 michaesp 1415
                  ! calculate minimum in circle
3 michaesp 1416
                  circlemin=circlef(1)
1417
                  do l=1,kst
1418
                     if (circlef(l) .lt. circlemin) then
1419
                        circlemin=circlef(l)
1420
                     endif
1421
                  enddo
5 michaesp 1422
                  !print*,'minimum of f in circle=',circlemin       
3 michaesp 1423
                  f0=circlemin
1424
               elseif ( intmode .eq. 'circle_max' ) then             
5 michaesp 1425
                  ! calculate maximum in circle
3 michaesp 1426
                  circlemax=circlef(1)
1427
                  do l=1,kst
1428
                     if (circlef(l) .gt. circlemax) then
1429
                        circlemax=circlef(l)
1430
                     endif
1431
                  enddo
5 michaesp 1432
                  !print*,'maximum of f in circle=',circlemax
3 michaesp 1433
                  f0=circlemax
1434
               else 
1435
                  print*,'ERROR (circle): intmode not valid!'
1436
                  stop
1437
               endif
1438
            else
1439
               print*,'ERROR (circle): Search radius is too small... (2). r =',radius
1440
               stop
1441
            endif
1442
 
5 michaesp 1443
            ! ------------------------ end CIRCLE modes -------------------------------
3 michaesp 1444
 
1445
            ! ------------------------ NORMAL mode -------------------------------
5 michaesp 1446
            else ! not clustering nor circle: NORMAL mode
3 michaesp 1447
 
5 michaesp 1448
            ! Check if point is within grid
1449
!            if ( (xind.ge.1.).and.(xind.le.real(nx)).and. &
1450
!                 (yind.ge.1.).and.(yind.le.real(ny)).and. &
1451
!                 (pind.ge.1.).and.(pind.le.real(nz)) ) then
1452
!
1453
            ! Do the interpolation: everthing is ok
39 michaesp 1454
            if ( ( shift_dir(i).ne.'PMIN' ).and.( shift_dir(i).ne.'PMAX' ) ) then
3 michaesp 1455
               f0 = int_index4(f3t0,f3t1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
39 michaesp 1456
 
1457
            elseif  ( shift_dir(i).eq.'PMIN' ) then
1458
               f0 = int_index4(f3t0,f3t1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
1459
 
1460
               ! Handle > and < operator
1461
               if ( (shift_rel(i).eq.'>').and.(f0.gt.shift_val(i)) ) then
1462
                 do while ( (f0.gt.shift_val(i)).and.(pind.lt.nz) )
1463
                     pind = pind + 0.1
1464
                     f0 = int_index4(f3t0,f3t1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
1465
                 enddo
1466
                 f0 = int_index4(p3t0,p3t1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
1467
 
1468
               elseif ( (shift_rel(i).eq.'<').and.(f0.lt.shift_val(i)) ) then
1469
                 do while ( (f0.lt.shift_val(i)).and.(pind.lt.nz) )
1470
                     pind = pind + 0.1
1471
                     f0 = int_index4(f3t0,f3t1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
1472
                 enddo
1473
                 f0 = int_index4(p3t0,p3t1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
1474
 
1475
              else
1476
                 f0 = mdv
1477
              endif
1478
 
1479
 
1480
           elseif  ( shift_dir(i).eq.'PMAX' ) then
1481
               f0 = int_index4(f3t0,f3t1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
1482
 
1483
               ! Handle > and < operator
1484
               if ( (shift_rel(i).eq.'>').and.(f0.gt.shift_val(i)) ) then
1485
                 do while ( (f0.gt.shift_val(i)).and.(pind.gt.1) )
1486
                     pind = pind - 0.1
1487
                     f0 = int_index4(f3t0,f3t1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
1488
                 enddo
1489
                 f0 = int_index4(p3t0,p3t1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
1490
 
1491
               elseif ( (shift_rel(i).eq.'<').and.(f0.lt.shift_val(i)) ) then
1492
                 do while ( (f0.lt.shift_val(i)).and.(pind.gt.1) )
1493
                     pind = pind - 0.1
1494
                     f0 = int_index4(f3t0,f3t1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
1495
                 enddo
1496
                 f0 = int_index4(p3t0,p3t1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
1497
 
1498
              else
1499
                 f0 = mdv
1500
              endif
1501
 
1502
          endif  
3 michaesp 1503
 
5 michaesp 1504
!            ! Check if pressure is outside, but rest okay: adjust to lowest or highest level
1505
!            elseif ( (xind.ge.1.).and.(xind.le.real(nx)).and. (yind.ge.1.).and.(yind.le.real(ny)) ) then ! only vertical problem
1506
!               if ( pind.gt.nz ) then ! pressure too low, index too high
1507
!                 pind = real(nz)
1508
!                 print*,' Warning: pressure too low -> adjusted to highest level, pind=nz.'
1509
!                 print*,'(x0,y0,p0)=',x0,y0,p0
1510
!               elseif (pind.lt.1.) then ! pressure too high, index too low
1511
!                 pind = 1.
1512
!                 print*,' Warning: pressure too high -> adjusted to lowest level, pind=1.'
1513
!                 print*,'(x0,y0,p0)=',x0,y0,p0
1514
!              endif
1515
!              f0 = int_index4(f3t0,f3t1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
3 michaesp 1516
 
5 michaesp 1517
!            ! Less than 10 outside
1518
!            elseif (noutside.lt.10) then
1519
!               print*,' ',trim(tvar(i)),' @ ',x0,y0,p0,'outside'
1520
!               f0       = mdv
1521
!               noutside = noutside + 1
1522
!
1523
!            ! More than 10 outside
1524
!            elseif (noutside.eq.10) then
1525
!               print*,' ...more than 10 outside...'
1526
!               f0       = mdv
1527
!               noutside = noutside + 1
1528
 
1529
!            ! Else (not everything okay and also not 'tolerated cases') set to missing data
1530
!            else
1531
!               f0       = mdv
1532
!            endif
1533
 
1534
            ! ------------------------ end NORMAL mode -------------------------------
1535
            endif ! end if nearest case
1536
 
39 michaesp 1537
           ! Exit for loop over all times -save interpolated value
5 michaesp 1538
 109        continue
1539
 
3 michaesp 1540
            ! Save the new field
1541
            if ( abs(f0-mdv).gt.eps) then
1542
               traint(k,j,ncol+i) = f0
1543
            else
1544
               traint(k,j,ncol+i) = mdv
1545
            endif
1546
 
1547
         enddo ! end loop over all trajectories
1548
 
1549
      enddo ! end loop over all times
1550
 
1551
      ! Exit point for loop over all tracing variables
1552
 110   continue
1553
 
1554
   enddo ! end loop over all variables
1555
 
1556
  ! --------------------------------------------------------------------
1557
  ! Calculate additional fields along the trajectories
1558
  ! --------------------------------------------------------------------
1559
 
1560
  print*
1561
  print*,'---- CALCULATE ADDITIONAL FIELDS FROM TRAJECTORY TABLE --'
1562
 
1563
  ! Loop over all tracing fields
1564
  do i=ntrace1,1,-1
1565
 
1566
      ! Skip fields which must not be computed (compfl=0)
1567
      if (compfl(i).eq.0) goto 120
1568
 
1569
      ! Write some status information
1570
      print*
1571
      write(*,'(a10,f10.2,a5,i3,3x,a2)') &
1572
           trim(tvar(i)),shift_val(i),trim(shift_dir(i)),compfl(i),trim(tfil(i))
1573
 
1574
      ! Loop over trajectories and times
1575
      do j=1,ntra
1576
      do k=1,ntim
1577
 
1578
         ! Potential temperature (TH)
1579
         if  ( varsint(i+ncol).eq.'TH' ) then
1580
 
1581
            varname='T'
1582
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1583
            varname='p'
1584
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1585
 
1586
            call calc_TH  (traint(j,k,ncol+i), traint(j,k,ind1), &
1587
                 traint(j,k,ind2) )
1588
 
1589
         ! Density (RHO)
1590
         elseif  ( varsint(i+ncol).eq.'RHO' ) then
1591
 
1592
            varname='T'
1593
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1594
            varname='p'
1595
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1596
 
1597
            call calc_RHO (traint(j,k,ncol+i), traint(j,k,ind1), &
1598
                 traint(j,k,ind2) )
1599
 
1600
         ! Relative humidity (RH)
1601
         elseif  ( varsint(i+ncol).eq.'RH' ) then
1602
 
1603
            varname='T'
1604
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1605
            varname='p'
1606
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1607
            varname='Q'
1608
            call list2ind (ind3,varname,varsint,fok,ncol+ntrace1)
1609
 
1610
            call calc_RH (traint(j,k,ncol+i), traint(j,k,ind1), &
1611
                 traint(j,k,ind2),traint(j,k,ind3) )
1612
 
1613
         ! Equivalent potential temperature (THE)
1614
         elseif  ( varsint(i+ncol).eq.'THE' ) then
1615
 
1616
            varname='T'
1617
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1618
            varname='p'
1619
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1620
            varname='Q'
1621
            call list2ind (ind3,varname,varsint,fok,ncol+ntrace1)
1622
 
1623
            call calc_THE (traint(j,k,ncol+i), traint(j,k,ind1), &
1624
                 traint(j,k,ind2),traint(j,k,ind3) )
1625
 
1626
         ! Latent heating rate (LHR)
1627
         elseif  ( varsint(i+ncol).eq.'LHR' ) then
1628
 
1629
            varname='T'
1630
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1631
            varname='p'
1632
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1633
            varname='Q'
1634
            call list2ind (ind3,varname,varsint,fok,ncol+ntrace1)
1635
            varname='OMEGA'
1636
            call list2ind (ind4,varname,varsint,fok,ncol+ntrace1)
1637
            varname='RH'
1638
            call list2ind (ind5,varname,varsint,fok,ncol+ntrace1)
1639
 
1640
            call calc_LHR (traint(j,k,ncol+i), traint(j,k,ind1), &
1641
                 traint(j,k,ind2),traint(j,k,ind3),traint(j,k,ind4),traint(j,k,ind5) )
1642
 
1643
         ! Wind speed (VEL)
1644
         elseif  ( varsint(i+ncol).eq.'VEL' ) then
1645
 
1646
            varname='U'
1647
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1648
            varname='V'
1649
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1650
 
1651
            call calc_VEL (traint(j,k,ncol+i), traint(j,k,ind1), &
1652
                 traint(j,k,ind2) )
1653
 
1654
         ! Wind direction (DIR)
1655
         elseif  ( varsint(i+ncol).eq.'DIR' ) then
1656
 
1657
            varname='U'
1658
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1659
            varname='V'
1660
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1661
 
1662
            call calc_DIR (traint(j,k,ncol+i), traint(j,k,ind1), &
1663
                 traint(j,k,ind2) )
1664
 
1665
         ! Zonal gradient of U (DUDX)
1666
         elseif  ( varsint(i+ncol).eq.'DUDX' ) then
1667
 
1668
            varname='U:+1DLON'
1669
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1670
            varname='U:-1DLON'
1671
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1672
            varname='lat'
1673
            call list2ind (ind3,varname,varsint,fok,ncol+ntrace1)
1674
 
1675
            call calc_DUDX (traint(j,k,ncol+i), traint(j,k,ind1), &
1676
                 traint(j,k,ind2),traint(j,k,ind3) )
1677
 
1678
         ! Zonal gradient of V (DVDX)
1679
         elseif  ( varsint(i+ncol).eq.'DVDX' ) then
1680
 
1681
            varname='V:+1DLON'
1682
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1683
            varname='V:-1DLON'
1684
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1685
            varname='lat'
1686
            call list2ind (ind3,varname,varsint,fok,ncol+ntrace1)
1687
 
1688
            call calc_DVDX (traint(j,k,ncol+i), traint(j,k,ind1), &
1689
                 traint(j,k,ind2),traint(j,k,ind3) )
1690
 
1691
         ! Zonal gradient of T (DTDX)
1692
         elseif  ( varsint(i+ncol).eq.'DVDX' ) then
1693
 
1694
            varname='T:+1DLON'
1695
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1696
            varname='T:-1DLON'
1697
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1698
            varname='lat'
1699
            call list2ind (ind3,varname,varsint,fok,ncol+ntrace1)
1700
 
1701
            call calc_DTDX (traint(j,k,ncol+i), traint(j,k,ind1), &
1702
                 traint(j,k,ind2),traint(j,k,ind3) )
1703
 
1704
         ! Zonal gradient of TH (DTHDX)
1705
         elseif  ( varsint(i+ncol).eq.'DTHDX' ) then
1706
 
1707
            varname='T:+1DLON'
1708
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1709
            varname='T:-1DLON'
1710
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1711
            varname='P'
1712
            call list2ind (ind3,varname,varsint,fok,ncol+ntrace1)
1713
            varname='lat'
1714
            call list2ind (ind4,varname,varsint,fok,ncol+ntrace1)
1715
 
1716
            call calc_DTHDX (traint(j,k,ncol+i), traint(j,k,ind1), &
1717
                 traint(j,k,ind2),traint(j,k,ind3),traint(j,k,ind4) )
1718
 
1719
         ! Meridional gradient of U (DUDY)
1720
         elseif  ( varsint(i+ncol).eq.'DUDY' ) then
1721
 
1722
            varname='U:+1DLAT'
1723
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1724
            varname='U:-1DLAT'
1725
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1726
 
1727
            call calc_DUDY (traint(j,k,ncol+i), traint(j,k,ind1), &
1728
                 traint(j,k,ind2) )
1729
 
1730
         ! Meridional gradient of V (DVDY)
1731
         elseif  ( varsint(i+ncol).eq.'DVDY' ) then
1732
 
1733
            varname='V:+1DLAT'
1734
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1735
            varname='V:-1DLAT'
1736
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1737
 
1738
            call calc_DVDY (traint(j,k,ncol+i), traint(j,k,ind1), &
1739
                 traint(j,k,ind2) )
1740
 
1741
         ! Meridional gradient of T (DTDY)
1742
         elseif  ( varsint(i+ncol).eq.'DTDY' ) then
1743
 
1744
            varname='T:+1DLAT'
1745
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1746
            varname='T:-1DLAT'
1747
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1748
 
1749
            call calc_DTDY (traint(j,k,ncol+i), traint(j,k,ind1), &
1750
                 traint(j,k,ind2) )
1751
 
1752
         ! Meridional gradient of TH (DTHDY)
1753
         elseif  ( varsint(i+ncol).eq.'DTHDY' ) then
1754
 
1755
            varname='T:+1DLAT'
1756
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1757
            varname='T:-1DLAT'
1758
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1759
            varname='P'
1760
            call list2ind (ind3,varname,varsint,fok,ncol+ntrace1)
1761
 
1762
            call calc_DTDY (traint(j,k,ncol+i), traint(j,k,ind1), &
1763
                 traint(j,k,ind2),traint(j,k,ind3) )
1764
 
1765
 
1766
         ! Vertical wind shear DU/DP (DUDP)
1767
         elseif  ( varsint(i+ncol).eq.'DUDP' ) then
1768
 
1769
            varname='U:+1DP'
1770
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1771
            varname='U:-1DP'
1772
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1773
            varname='P:+1DP'
1774
            call list2ind (ind3,varname,varsint,fok,ncol+ntrace1)
1775
            varname='P:-1DP'
1776
            call list2ind (ind4,varname,varsint,fok,ncol+ntrace1)
1777
 
1778
            call calc_DUDP (traint(j,k,ncol+i), traint(j,k,ind1), &
1779
                 traint(j,k,ind2),traint(j,k,ind3),traint(j,k,ind4) )
1780
 
1781
         ! Vertical wind shear DV/DP (DVDP)
1782
         elseif  ( varsint(i+ncol).eq.'DVDP' ) then
1783
 
1784
            varname='V:+1DP'
1785
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1786
            varname='V:-1DP'
1787
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1788
            varname='P:+1DP'
1789
            call list2ind (ind3,varname,varsint,fok,ncol+ntrace1)
1790
            varname='P:-1DP'
1791
            call list2ind (ind4,varname,varsint,fok,ncol+ntrace1)
1792
 
1793
            call calc_DVDP (traint(j,k,ncol+i), traint(j,k,ind1), &
1794
                 traint(j,k,ind2),traint(j,k,ind3),traint(j,k,ind4) )
1795
 
1796
         ! Vertical derivative of T (DTDP)
1797
         elseif  ( varsint(i+ncol).eq.'DTDP' ) then
1798
 
1799
            varname='T:+1DP'
1800
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1801
            varname='T:-1DP'
1802
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1803
            varname='P:+1DP'
1804
            call list2ind (ind3,varname,varsint,fok,ncol+ntrace1)
1805
            varname='P:-1DP'
1806
            call list2ind (ind4,varname,varsint,fok,ncol+ntrace1)
1807
 
1808
            call calc_DTDP (traint(j,k,ncol+i), traint(j,k,ind1), &
1809
                 traint(j,k,ind2),traint(j,k,ind3),traint(j,k,ind4) )
1810
 
1811
         ! Vertical derivative of TH (DTHDP)
1812
         elseif  ( varsint(i+ncol).eq.'DTHDP' ) then
1813
 
1814
            varname='T:+1DP'
1815
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1816
            varname='T:-1DP'
1817
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1818
            varname='P:+1DP'
1819
            call list2ind (ind3,varname,varsint,fok,ncol+ntrace1)
1820
            varname='P:-1DP'
1821
            call list2ind (ind4,varname,varsint,fok,ncol+ntrace1)
1822
            varname='P'
1823
            call list2ind (ind5,varname,varsint,fok,ncol+ntrace1)
1824
            varname='T'
1825
            call list2ind (ind6,varname,varsint,fok,ncol+ntrace1)
1826
 
1827
            call calc_DTHDP (traint(j,k,ncol+i), traint(j,k,ind1), &
1828
                 traint(j,k,ind2),traint(j,k,ind3),traint(j,k,ind4),traint(j,k,ind5),traint(j,k,ind6) )
1829
 
1830
         ! Squared Brunt-Vaisäla frequency (NSQ)
1831
         elseif  ( varsint(i+ncol).eq.'NSQ' ) then
1832
 
1833
            varname='DTHDP'
1834
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1835
            varname='TH'
1836
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1837
            varname='RHO'
1838
            call list2ind (ind3,varname,varsint,fok,ncol+ntrace1)
1839
 
1840
            call calc_NSQ (traint(j,k,ncol+i), traint(j,k,ind1), &
1841
                 traint(j,k,ind2),traint(j,k,ind3))
1842
 
1843
         ! Relative vorticity (RELVORT)
1844
         elseif  ( varsint(i+ncol).eq.'RELVORT' ) then
1845
 
1846
            varname='DUDY'
1847
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1848
            varname='DVDX'
1849
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1850
            varname='U'
1851
            call list2ind (ind3,varname,varsint,fok,ncol+ntrace1)
1852
            varname='lat'
1853
            call list2ind (ind4,varname,varsint,fok,ncol+ntrace1)
1854
 
1855
            call calc_RELVORT (traint(j,k,ncol+i), traint(j,k,ind1), &
1856
                 traint(j,k,ind2),traint(j,k,ind3),traint(j,k,ind4))
1857
 
1858
         ! Absolute vorticity (ABSVORT)
1859
         elseif  ( varsint(i+ncol).eq.'ABSVORT' ) then
1860
 
1861
            varname='DUDY'
1862
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1863
            varname='DVDX'
1864
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1865
            varname='U'
1866
            call list2ind (ind3,varname,varsint,fok,ncol+ntrace1)
1867
            varname='lat'
1868
            call list2ind (ind4,varname,varsint,fok,ncol+ntrace1)
1869
 
1870
            call calc_ABSVORT (traint(j,k,ncol+i), traint(j,k,ind1), &
1871
                 traint(j,k,ind2),traint(j,k,ind3),traint(j,k,ind4))
1872
 
1873
         ! Divergence (DIV)
1874
         elseif  ( varsint(i+ncol).eq.'DIV' ) then
1875
 
1876
            varname='DUDX'
1877
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1878
            varname='DVDY'
1879
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1880
            varname='V'
1881
            call list2ind (ind3,varname,varsint,fok,ncol+ntrace1)
1882
            varname='lat'
1883
            call list2ind (ind4,varname,varsint,fok,ncol+ntrace1)
1884
 
1885
            call calc_DIV (traint(j,k,ncol+i), traint(j,k,ind1), &
1886
                 traint(j,k,ind2),traint(j,k,ind3),traint(j,k,ind4))
1887
 
1888
         ! Deformation (DEF)
1889
         elseif  ( varsint(i+ncol).eq.'DEF' ) then
1890
 
1891
            varname='DUDX'
1892
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1893
            varname='DVDX'
1894
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1895
            varname='DUDY'
1896
            call list2ind (ind3,varname,varsint,fok,ncol+ntrace1)
1897
            varname='DVDY'
1898
            call list2ind (ind4,varname,varsint,fok,ncol+ntrace1)
1899
 
1900
            call calc_DEF (traint(j,k,ncol+i), traint(j,k,ind1), &
1901
                 traint(j,k,ind2),traint(j,k,ind3),traint(j,k,ind4))
1902
 
1903
         ! Potential Vorticity (PV)
1904
         elseif  ( varsint(i+ncol).eq.'PV' ) then
1905
 
1906
            varname='ABSVORT'
1907
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1908
            varname='DTHDP'
1909
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1910
            varname='DUDP'
1911
            call list2ind (ind3,varname,varsint,fok,ncol+ntrace1)
1912
            varname='DVDP'
1913
            call list2ind (ind4,varname,varsint,fok,ncol+ntrace1)
1914
            varname='DTHDX'
1915
            call list2ind (ind5,varname,varsint,fok,ncol+ntrace1)
1916
            varname='DTHDY'
1917
            call list2ind (ind6,varname,varsint,fok,ncol+ntrace1)
1918
 
1919
            call calc_PV (traint(j,k,ncol+i), traint(j,k,ind1), &
1920
                 traint(j,k,ind2),traint(j,k,ind3),traint(j,k,ind4),traint(j,k,ind5),traint(j,k,ind6) )
1921
 
1922
         ! Richardson number (RI)
1923
         elseif  ( varsint(i+ncol).eq.'RI' ) then
1924
 
1925
            varname='DUDP'
1926
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1927
            varname='DVDP'
1928
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1929
            varname='NSQ'
1930
            call list2ind (ind3,varname,varsint,fok,ncol+ntrace1)
1931
            varname='RHO'
1932
            call list2ind (ind4,varname,varsint,fok,ncol+ntrace1)
1933
 
1934
            call calc_RI (traint(j,k,ncol+i), traint(j,k,ind1), &
1935
                 traint(j,k,ind2),traint(j,k,ind3),traint(j,k,ind4) )
1936
 
1937
         ! Ellrod and Knapp's turbulence idicator (TI)
1938
         elseif  ( varsint(i+ncol).eq.'TI' ) then
1939
 
1940
            varname='DEF'
1941
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1942
            varname='DUDP'
1943
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1944
            varname='DVDP'
1945
            call list2ind (ind3,varname,varsint,fok,ncol+ntrace1)
1946
            varname='RHO'
1947
            call list2ind (ind4,varname,varsint,fok,ncol+ntrace1)
1948
 
1949
            call calc_TI (traint(j,k,ncol+i), traint(j,k,ind1), &
1950
                 traint(j,k,ind2),traint(j,k,ind3),traint(j,k,ind4) )
1951
 
1952
         ! Spherical distance from starting position (DIST0)
1953
         elseif  ( varsint(i+ncol).eq.'DIST0' ) then
1954
 
1955
            varname='lon'
1956
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1957
            varname='lat'
1958
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1959
 
1960
            call calc_DIST0 (traint(j,k,ncol+i), traint(j,k,ind1), &
1961
                 traint(j,k,ind2),traint(j,1,ind1),traint(j,1,ind2) )
1962
 
1963
         ! Spherical distance length of trajectory (DIST)
1964
         elseif  ( varsint(i+ncol).eq.'DIST' ) then
1965
 
1966
            varname='lon'
1967
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1968
            varname='lat'
1969
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1970
 
1971
            if ( k.eq.1 ) then
1972
               traint(j,k,ncol+i) = 0.
1973
            else
1974
               call calc_DIST0 (delta, traint(j,k  ,ind1), &
1975
                   traint(j,k  ,ind2),traint(j,k-1,ind1),traint(j,k-1,ind2) )
1976
               traint(j,k,ncol+i) = traint(j,k-1,ncol+i) + delta
1977
            endif
1978
 
1979
         ! Heading of the trajectory (HEAD)
1980
         elseif  ( varsint(i+ncol).eq.'HEAD' ) then
1981
 
1982
            varname='lon'
1983
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1984
            varname='lat'
1985
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
1986
 
1987
            if (k.eq.ntim) then
1988
               traint(j,k,ncol+i) = mdv
1989
            else
1990
               call calc_HEAD (traint(j,k,ncol+i), &
1991
                    traint(j,k  ,ind1),traint(j,k  ,ind2),traint(j,k+1,ind1),traint(j,k+1,ind2) )
1992
            endif
44 michaesp 1993
 
1994
        ! Directional change (DANGLE)
1995
         elseif  ( varsint(i+ncol).eq.'DANGLE' ) then
3 michaesp 1996
 
44 michaesp 1997
            varname='lon'
1998
            call list2ind (ind1,varname,varsint,fok,ncol+ntrace1)
1999
            varname='lat'
2000
            call list2ind (ind2,varname,varsint,fok,ncol+ntrace1)
3 michaesp 2001
 
44 michaesp 2002
            if (k.eq.ntim) then
2003
               traint(j,k,ncol+i) = mdv
2004
            else
2005
               if ( k.eq.1 ) then
2006
                   traint(j,k,ncol+i) = mdv
2007
               elseif ( k.eq.ntim ) then
2008
                   traint(j,k,ncol+i) = mdv
2009
               else
2010
                   call calc_DANGLE (traint(j,k,ncol+i),      &
2011
                       traint(j,k-1,ind1),traint(j,k-1,ind2), &
2012
                       traint(j,k  ,ind1),traint(j,k  ,ind2), &
2013
                       traint(j,k+1,ind1),traint(j,k+1,ind2) )
2014
               endif
2015
            endif
2016
 
2017
 
3 michaesp 2018
         ! Invalid tracing variable
2019
         else
2020
 
2021
            print*,' ERROR: invalid tracing variable ', trim(varsint(i+ncol))
2022
            stop
2023
 
2024
 
2025
         endif
2026
 
2027
      ! End loop over all trajectories and times
2028
      enddo
2029
      enddo
2030
 
2031
      ! Set the flag for a ready field/column
2032
      fok(ncol+i) = 1
2033
 
2034
 
2035
      ! Exit point for loop over all tracing fields
2036
 120   continue
2037
 
2038
   enddo
2039
 
2040
  ! --------------------------------------------------------------------
2041
  ! Write output to output trajectory file
2042
  ! --------------------------------------------------------------------
2043
 
2044
  ! Write status information
2045
  print*
2046
  print*,'---- WRITE OUTPUT TRAJECTORIES --------------------------'
2047
  print*
2048
 
2049
  ! Allocate memory for internal trajectories
2050
  allocate(traout(ntra,ntim,ncol+ntrace0),stat=stat)
2051
  if (stat.ne.0) print*,'*** error allocating array traout   ***'
2052
 
2053
  ! Copy input to output trajectory (apply scaling of output)
2054
  do i=1,ntra
2055
     do j=1,ntim
2056
        do k=1,ncol+ntrace0
2057
           if ( k.le.ncol ) then
2058
              traout(i,j,k) = traint(i,j,k)
2059
           elseif ( abs(traint(i,j,k)-mdv).gt.eps ) then
2060
              traout(i,j,k) = fac(k-ncol) * traint(i,j,k)
2061
           else
2062
              traout(i,j,k) = mdv
2063
           endif
2064
        enddo
2065
     enddo
2066
  enddo
2067
 
2068
  ! Set the variable names for output trajectory
2069
  do i=1,ncol+ntrace0
39 michaesp 2070
        varsout(i)      = varsint(i)
3 michaesp 2071
  enddo
39 michaesp 2072
  do i=1,ntrace0
2073
     if ( (shift_dir(i).eq.'PMIN').or.(shift_dir(i).eq.'PMAX') ) then
2074
        varsout(ncol+i) = trim(tvar(i))//':'//trim(shift_dir(i))
2075
     endif
2076
  enddo
2077
 
3 michaesp 2078
  ! Write trajectories
2079
  call wopen_tra(fod,outfile,ntra,ntim,ncol+ntrace0, &
2080
       reftime,varsout,outmode)
2081
  call write_tra(fod,traout ,ntra,ntim,ncol+ntrace0,outmode)
2082
  call close_tra(fod,outmode)
2083
 
2084
  ! Write some status information, and end of program message
2085
  print*
2086
  print*,'---- STATUS INFORMATION --------------------------------'
2087
  print*
2088
  print*,' ok'
2089
  print*
2090
  print*,'              *** END OF PROGRAM TRACE ***'
2091
  print*,'========================================================='
2092
 
2093
 
5 michaesp 2094
end program trace
3 michaesp 2095
 
2096
 
2097
 
2098
! ******************************************************************
2099
! * SUBROUTINE SECTION                                             *
2100
! ******************************************************************
2101
 
2102
! ------------------------------------------------------------------
2103
! Add a variable to the list if not yet included in this list
2104
! ------------------------------------------------------------------
2105
 
2106
subroutine add2list (varname,list,nlist)
2107
 
2108
  implicit none
2109
 
2110
  ! Declaration of subroutine parameters
2111
  character(len=80) :: varname
2112
  character(len=80) :: list(200)
2113
  integer :: nlist
2114
 
2115
  ! Auxiliray variables
2116
  integer :: i,j
2117
  integer :: isok
2118
 
2119
  ! Expand the list, if necessary
2120
  isok = 0
2121
  do i=1,nlist
2122
     if ( list(i).eq.varname ) isok = 1
2123
  enddo
2124
  if ( isok.eq.0 ) then
2125
     nlist       = nlist + 1
2126
     list(nlist) = varname
2127
  endif
2128
 
2129
  ! Check for too large number of fields
2130
  if ( nlist.ge.200) then
2131
     print*,' ERROR: too many additional fields for tracing ...'
2132
     stop
2133
  endif
2134
 
2135
end subroutine add2list
2136
 
2137
 
2138
! ------------------------------------------------------------------
2139
! Get the index of a variable in the list
2140
! ------------------------------------------------------------------
2141
 
2142
subroutine list2ind (ind,varname,list,fok,nlist)
2143
 
2144
  implicit none
2145
 
2146
  ! Declaration of subroutine parameters
2147
  integer :: ind
2148
  character(len=80) :: varname
2149
  character(len=80) :: list(200)
2150
  integer :: fok(200)
2151
  integer :: nlist
2152
 
2153
  ! Auxiliray variables
2154
  integer :: i,j
2155
  integer :: isok
2156
 
2157
  ! Get the index - error message if not found
2158
  ind = 0
2159
  do i=1,nlist
2160
     if ( varname.eq.list(i) ) then
2161
        ind = i
2162
        goto 100
2163
     endif
2164
  enddo
2165
 
2166
  if ( ind.eq.0) then
2167
     print*
2168
     print*,' ERROR: cannot find ',trim(varname),' in list ...'
2169
     do i=1,nlist
2170
        print*,i,trim(list(i))
2171
     enddo
2172
     print*
2173
     stop
2174
  endif
2175
 
2176
  ! Exit point
2177
 100   continue
2178
 
2179
  ! Check whether the field/column is ready
2180
  if ( fok(ind).eq.0 ) then
2181
     print*
2182
     print*,' ERROR: unresolved dependence : ',trim(list(ind))
2183
     print*
2184
     stop
2185
  endif
2186
 
2187
end subroutine list2ind
2188
 
2189
 
2190
! ------------------------------------------------------------------
2191
! Split the variable name into name, shift and direction
2192
! ------------------------------------------------------------------
2193
 
39 michaesp 2194
subroutine splitvar (tvar,shift_val,shift_dir,shift_rel)
3 michaesp 2195
 
2196
  implicit none
2197
 
2198
  ! Declaration of subroutine parameters
2199
  character(len=80) :: tvar
2200
  real :: shift_val
2201
  character(len=80) :: shift_dir
39 michaesp 2202
  character(len=80) :: shift_rel
3 michaesp 2203
 
2204
  ! Auxiliary variables
2205
  integer :: i,j
39 michaesp 2206
  integer :: icolon,inumber,irelator
3 michaesp 2207
  character(len=80) :: name
2208
  character :: ch
15 michaesp 2209
  integer      isabsval
3 michaesp 2210
 
2211
  ! Save variable name
2212
  name = tvar
39 michaesp 2213
  shift_rel = 'nil'
2214
  shift_dir = 'nil'
3 michaesp 2215
 
2216
  ! Search for colon
2217
  icolon=0
2218
  do i=1,80
2219
     if ( (name(i:i).eq.':').and.(icolon.ne.0) ) goto 100
2220
     if ( (name(i:i).eq.':').and.(icolon.eq.0) ) icolon=i
2221
  enddo
2222
 
2223
  ! If there is a colon, split the variable name
2224
  if ( icolon.ne.0 ) then
2225
 
2226
     tvar = name(1:(icolon-1))
2227
 
15 michaesp 2228
     ! Get the index for number
3 michaesp 2229
     do i=icolon+1,80
2230
        ch = name(i:i)
2231
        if ( ( ch.ne.'0' ).and. ( ch.ne.'1' ).and.( ch.ne.'2' ).and. &
2232
             ( ch.ne.'3' ).and. ( ch.ne.'4' ).and.( ch.ne.'5' ).and. &
2233
             ( ch.ne.'6' ).and. ( ch.ne.'7' ).and.( ch.ne.'8' ).and. &
2234
             ( ch.ne.'9' ).and. ( ch.ne.'+' ).and.( ch.ne.'-' ).and. &
2235
             ( ch.ne.'.' ).and. ( ch.ne.' ' )  ) then
2236
           inumber = i
2237
           exit
2238
        endif
2239
     enddo
2240
 
39 michaesp 2241
     ! If the variable name is e.g. PMIN:UMF>0, the variable to be read
2242
     ! is UMF, the value is 0, and the direction is 'PMIN'
2243
     shift_rel = 'nil'
2244
     if ( (tvar.eq.'PMIN').or.(tvar.eq.'PMAX') ) then
2245
         shift_dir = tvar
2246
         irelator = 0
2247
         do i=icolon+1,80
2248
            ch = name(i:i)
2249
            if ( (ch.eq.'>').or.(ch.eq.'<') ) then
2250
               irelator = i
2251
            endif
2252
         enddo
2253
         if ( irelator.eq.0 ) then
2254
             print*,' ERROR: dont know how to interpret ',trim(name)
2255
             stop
2256
         endif
2257
         tvar = name(icolon+1:irelator-1)
2258
         read(name( (irelator+1):80 ),*) shift_val
2259
         shift_rel = name(irelator:irelator)
2260
 
2261
         goto 90
2262
 
2263
     endif
2264
 
15 michaesp 2265
     ! Get the number
3 michaesp 2266
     read(name( (icolon+1):(inumber-1) ),*) shift_val
2267
 
15 michaesp 2268
     ! Decide whether it is a shift relatiev to trajectory or absolute value
2269
     ! If the number starts with + or -, it is relative to the trajectory
2270
     isabsval = 1
2271
     do i=icolon+1,inumber-1
2272
       ch = name(i:i)
2273
       if ( (ch.eq.'+').or.(ch.eq.'-') ) isabsval = 0
2274
     enddo
2275
 
2276
     ! Get the unit/shift axis
3 michaesp 2277
     shift_dir = name(inumber:80)
15 michaesp 2278
     if ( isabsval.eq.1 ) then
2279
       shift_dir=trim(shift_dir)//'(ABS)'
2280
     endif
3 michaesp 2281
 
2282
  else
2283
 
2284
     shift_dir = 'nil'
2285
     shift_val = 0.
2286
 
2287
  endif
2288
 
39 michaesp 2289
 90 continue
2290
 
2291
    return
3 michaesp 2292
 
2293
  ! Error handling
2294
 100   continue
2295
 
2296
  print*,' ERROR: cannot split variable name ',trim(tvar)
2297
  stop
2298
 
2299
end subroutine splitvar