Subversion Repositories lagranto.ecmwf

Rev

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