Subversion Repositories lagranto.ecmwf

Rev

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