Subversion Repositories lagranto.wrf

Rev

Rev 20 | Go to most recent revision | Details | Last modification | View Log | RSS feed

Rev Author Line No. Line
2 michaesp 1
 
2
PROGRAM trace_label_circle
3
 
4
  ! ********************************************************************
5
  ! *                                                                  *
6
  ! * Trace fields along trajectories                                  *
7
  ! *                                                                  *
8
  ! * Heini Wernli       first version:       April 1993               *
9
  ! * Michael Sprenger   major upgrade:       2008-2009                *
10
  ! * Sebastian Schemm WRF                                             *
11
  ! ********************************************************************
12
 
13
  implicit none
14
 
15
  ! --------------------------------------------------------------------
16
  ! Declaration of parameters
17
  ! --------------------------------------------------------------------
18
 
19
  ! Maximum number of levels for input files
20
  integer :: nlevmax
21
  parameter (nlevmax=100)
22
 
23
  ! Maximum number of input files (dates, length of trajectories)
24
  integer :: ndatmax
25
  parameter (ndatmax=500)
26
 
27
  ! Numerical epsilon (for float comparison)
28
  real :: eps
29
  parameter (eps=0.001)
30
 
31
  ! Prefix for primary and secondary fields
32
  character :: charp
33
  character :: chars
34
  parameter (charp='P')
35
  parameter (chars='S')
36
 
37
  ! --------------------------------------------------------------------
38
  ! Declaration of variables
39
  ! --------------------------------------------------------------------
40
 
41
  ! Input and output format for trajectories (see iotra.f)
42
  integer :: inpmode
43
  integer :: outmode
44
 
45
  ! Input parameters
46
  character(len=80) :: inpfile         ! Input trajectory file
47
  character(len=80) :: outfile         ! Output trajectory file
48
  integer :: ntra                      ! Number of trajectories
49
  integer :: ncol                      ! Number of columns (including time, lon, lat, p)
50
  integer :: ntim                      ! Number of times per trajectory
51
  integer :: ntrace0                   ! Number of trace variables
52
  character(len=80) :: tvar0(200)      ! Tracing variable (with mode specification)
53
  character(len=80) :: tvar(200)       ! Tracing variable name (only the variable)
54
  character(len=1)  :: tfil(200)       ! Filename prefix
55
  real :: fac(200)                     ! Scaling factor
56
  integer :: compfl(200)               ! Computation flag (1=compute)
57
  integer :: numdat                    ! Number of input files
58
  character(len=11) :: dat(ndatmax)    ! Dates of input files
59
  real :: timeinc                      ! Time increment between input files
60
  real :: tst                          ! Time shift of start relative to first data file
61
  real :: ten                          ! Time shift of end relatiev to first data file
62
  character(len=20) :: startdate       ! First time/date on trajectory
63
  character(len=20) :: enddate         ! Last time/date on trajectory
64
  character(len=80) :: timecheck       ! Either 'yes' or 'no'
65
  character(len=80) :: intmode         ! Interpolation mode ('normal', 'nearest')
66
 
67
  ! Trajectories
68
  real,allocatable, dimension (:,:,:) :: trainp          ! Input trajectories (ntra,ntim,ncol)
69
  real,allocatable, dimension (:,:,:) :: traint          ! Internal trajectories (ntra,ntim,ncol+ntrace0)
70
  real,allocatable, dimension (:,:,:) :: traout          ! Output trajectories (ntra,ntim,ncol+ntrace0)
71
  integer :: reftime(6)                                  ! Reference date
72
  character(len=80) :: varsinp(100)                      ! Field names for input trajectory
73
  character(len=80) :: varsint(100)                      ! Field names for internal trajectory
74
  character(len=80) :: varsout(100)                      ! Field names for output trajectory
75
  integer :: fid,fod                                     ! File identifier for inp and out trajectories
76
  real :: x0,y0,p0                                       ! Position of air parcel (physical space)
77
  real :: reltpos0                                       ! Relative time of air parcel
78
  real :: xind,yind,pind                                 ! Position of air parcel (grid space)
79
  integer :: fbflag                                      ! Flag for forward (1) or backward (-1) trajectories
80
  integer :: fok(100)                                    ! Flag whether field is ready
81
 
82
  ! Meteorological fields
83
  real,allocatable, dimension (:)     :: spt0,spt1       ! Surface pressure
84
  real,allocatable, dimension (:)     :: p3t0,p3t1       ! 3d-pressure
85
  real,allocatable, dimension (:)     :: f3t0,f3t1       ! 3d field for tracing
86
  character(len=80) :: svars(100)                        ! List of variables on S file
87
  character(len=80) :: pvars(100)                        ! List of variables on P file
88
  integer :: n_svars                                     ! Number of variables on S file
89
  integer :: n_pvars                                     ! Number of variables on P file
90
 
91
  ! Grid description
92
  real :: pollon,pollat   ! Longitude/latitude of pole
93
  real :: xmin,xmax       ! Zonal grid extension
94
  real :: ymin,ymax       ! Meridional grid extension
95
  integer :: nx,ny,nz     ! Grid dimensions
96
  real :: dx,dy           ! Horizontal grid resolution
97
  integer :: hem          ! Flag for hemispheric domain
98
  integer :: per          ! Flag for periodic domain
99
  real :: stagz           ! Vertical staggering
100
  real :: mdv             ! Missing data value
101
 
102
  ! Auxiliary variables
103
  integer :: i,j,k,l,m,n
104
  real :: rd
105
  character(len=80) :: filename,varname
106
  real :: time0,time1,reltpos
107
  integer :: itime0,itime1
108
  integer :: stat
109
  real :: tstart
110
  integer :: iloaded0,iloaded1
111
  real :: f0
112
  real :: frac
113
  real :: tload,tfrac
114
  integer :: isok
115
  character :: ch
116
  integer :: ind
117
  integer :: ind1,ind2,ind3,ind4,ind5
118
  integer :: ind6,ind7,ind8,ind9,ind0
119
  integer :: noutside
120
  real :: delta
121
  integer :: itrace0
122
  character(len=80) :: string
123
 
124
  ! Externals
125
  real :: int_index4
126
  external int_index4
127
  real :: sdis   
128
  external sdis 
129
 
130
  ! --------------------------------------------------------------------
131
  ! Start of program, Read parameters, get grid parameters
132
  ! --------------------------------------------------------------------
133
 
134
  ! Write start message
135
  print*,'========================================================='
136
  print*,'              *** START OF PROGRAM TRACE ***'
137
  print*
138
 
139
  ! Read parameters
140
  open(10,file='trace.param')
141
   read(10,*) inpfile
142
   read(10,*) outfile
143
   read(10,*) startdate
144
   read(10,*) enddate
145
   read(10,*) fbflag
146
   read(10,*) numdat
147
   if ( fbflag.eq.1) then
148
      do i=1,numdat
149
         read(10,'(a11)') dat(i)
150
      enddo
151
   else
152
      do i=numdat,1,-1
153
         read(10,'(a11)') dat(i)
154
      enddo
155
   endif
156
   read(10,*) timeinc
157
   read(10,*) tst
158
   read(10,*) ten
159
   read(10,*) ntra
160
   read(10,*) ntim
161
   read(10,*) ncol
162
   read(10,*) ntrace0
163
   do i=1,ntrace0
164
      read(10,*) tvar(i), fac(i), compfl(i), tfil(i)
165
   enddo
166
   read(10,*) n_pvars
167
   do i=1,n_pvars
168
      read(10,*) pvars(i)
169
   enddo
170
   read(10,*) n_svars
171
   do i=1,n_svars
172
      read(10,*) svars(i)
173
   enddo
174
   read(10,*) timecheck
175
   read(10,*) intmode
176
  close(10)
177
 
178
  ! Remove commented tracing fields
179
  itrace0 = 1
180
  do while ( itrace0.le.ntrace0)
181
     string = tvar(itrace0)
182
     if ( string(1:1).eq.'#' ) then
183
        do i=itrace0,ntrace0-1
184
           tvar(i)   = tvar(i+1)
185
           fac(i)    = fac(i+1)
186
           compfl(i) = compfl(i+1)
187
           tfil(i)   = tfil(i+1)
188
        enddo
189
        ntrace0 = ntrace0 - 1
190
     else
191
        itrace0 = itrace0 + 1
192
     endif
193
  enddo
194
 
195
  ! Save the tracing variable (including all mode specifications)
196
  do i=1,ntrace0
197
     tvar0(i) = tvar(i)
198
  enddo
199
 
200
  ! Set the formats of the input and output files
201
  call mode_tra(inpmode,inpfile)
202
  if (inpmode.eq.-1) inpmode=1
203
  call mode_tra(outmode,outfile)
204
  if (outmode.eq.-1) outmode=1
205
 
206
  ! Convert time shifts <tst,ten> from <hh.mm> into fractional time
207
  call hhmm2frac(tst,frac)
208
  tst = frac
209
  call hhmm2frac(ten,frac)
210
  ten = frac
211
 
212
  ! Set the time for the first data file (depending on forward/backward mode)
213
  if (fbflag.eq.1) then
214
    tstart = -tst
215
  else
216
    tstart = tst
217
  endif
218
 
219
  ! Read the constant grid parameters (nx,ny,nz,xmin,xmax,ymin,ymax,pollon,pollat)
220
  ! The negative <-fid> of the file identifier is used as a flag for parameter retrieval
221
  filename = charp//dat(1)
222
  varname  = 'U'
223
  call input_open (fid,filename)
224
  call input_grid_wrf(fid,xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz)
225
  call input_close(fid)
226
 
227
  ! Allocate memory for some meteorological arrays
228
  allocate(spt0(nx*ny),stat=stat)
229
  if (stat.ne.0) print*,'*** error allocating array spt0 ***'   ! Surface pressure
230
  allocate(spt1(nx*ny),stat=stat)
231
  if (stat.ne.0) print*,'*** error allocating array spt1 ***'
232
  allocate(p3t0(nx*ny*nz),stat=stat)
233
  if (stat.ne.0) print*,'*** error allocating array p3t0 ***'   ! Pressure
234
  allocate(p3t1(nx*ny*nz),stat=stat)
235
  if (stat.ne.0) print*,'*** error allocating array p3t1 ***'
236
  allocate(f3t0(nx*ny*nz),stat=stat)
237
  if (stat.ne.0) print*,'*** error allocating array p3t0 ***'   ! Tracing field
238
  allocate(f3t1(nx*ny*nz),stat=stat)
239
  if (stat.ne.0) print*,'*** error allocating array p3t1 ***'
240
 
241
  ! Get memory for trajectory arrays
242
  allocate(trainp(ntra,ntim,ncol),stat=stat)
243
  if (stat.ne.0) print*,'*** error allocating array tra      ***'
244
 
245
  ! Write some status information
246
  print*,'---- INPUT PARAMETERS -----------------------------------'
247
  print*
248
  print*,'  Input trajectory file  : ',trim(inpfile)
249
  print*,'  Output trajectory file : ',trim(outfile)
250
  print*,'  Format of input file   : ',inpmode
251
  print*,'  Format of output file  : ',outmode
252
  print*,'  Forward/backward       : ',fbflag
253
  print*,'  #tra                   : ',ntra
254
  print*,'  #col                   : ',ncol
255
  print*,'  #tim                   : ',ntim
256
  print*,'  No time check          : ',trim(timecheck)
257
  print*,'  Interpolation mode     : ',trim(intmode)
258
  do i=1,ntrace0
259
     if (compfl(i).eq.0) then
260
        print*,'  Tracing field          : ', trim(tvar(i)), fac(i), ' 0 ', tfil(i)
261
     else
262
        print*,'  Tracing field          : ', trim(tvar(i)), fac(i), '  1 ', tfil(i)
263
     endif
264
  enddo
265
  print*
266
  print*,'---- INPUT DATA FILES -----------------------------------'
267
  print*
268
  call frac2hhmm(tstart,tload)
269
  print*,'  Time of 1st data file  : ',tload
270
  print*,'  #input files           : ',numdat
271
  print*,'  time increment         : ',timeinc
272
  call frac2hhmm(tst,tload)
273
  print*,'  Shift of start         : ',tload
274
  call frac2hhmm(ten,tload)
275
  print*,'  Shift of end           : ',tload
276
  print*,'  First/last input file  : ',trim(dat(1)), ' ... ', trim(dat(numdat))
277
  print*,'  Primary variables      : ',trim(pvars(1))
278
  do i=2,n_pvars
279
     print*,'                         : ',trim(pvars(i))
280
  enddo
281
  if ( n_svars.ge.1 ) then
282
     print*,'  Secondary variables    : ',trim(svars(1))
283
     do i=2,n_svars
284
        print*,'                         : ',trim(svars(i))
285
     enddo
286
  endif
287
  print*
288
  print*,'---- CONSTANT GRID PARAMETERS ---------------------------'
289
  print*
290
  print*,'  xmin,xmax     : ',xmin,xmax
291
  print*,'  ymin,ymax     : ',ymin,ymax
292
  print*,'  dx,dy         : ',dx,dy
293
  print*,'  pollon,pollat : ',pollon,pollat
294
  print*,'  nx,ny,nz      : ',nx,ny,nz
295
  print*,'  per, hem      : ',per,hem
296
  print*
297
 
298
  ! --------------------------------------------------------------------
299
  ! Load the input trajectories
300
  ! --------------------------------------------------------------------
301
 
302
  ! Read the input trajectory file
303
  call ropen_tra(fid,inpfile,ntra,ntim,ncol,reftime,varsinp,inpmode)
304
  call read_tra (fid,trainp,ntra,ntim,ncol,inpmode)
305
  call close_tra(fid,inpmode)
306
 
307
  ! Check that first four columns correspond to time,x,y,z
308
  if ( (varsinp(1).ne.'time' ).or. &
309
       (varsinp(2).ne.'x'    ).or. &
310
       (varsinp(3).ne.'y'    ).or. &
311
       (varsinp(4).ne.'z'   ) )    &
312
  then
313
     print*,' ERROR: problem with input trajectories ...'
314
     stop
315
  endif
316
  varsinp(1) = 'time'
317
  varsinp(2) = 'x'
318
  varsinp(3) = 'y'
319
  varsinp(4) = 'z'
320
 
321
  ! Transform start coordinates from index space to WRF grid
322
  do i=1,ntra
323
     do j=1,ntim
324
         trainp(i,j,2) = xmin + ( trainp(i,j,2) - 1. ) * dx
325
         trainp(i,j,3) = ymin + ( trainp(i,j,3) - 1. ) * dy
326
     enddo
327
  enddo
328
 
329
  ! Write some status information of the input trajectories
330
  print*,'---- INPUT TRAJECTORIES ---------------------------------'
331
  print*
332
  print*,' Start date             : ',trim(startdate)
333
  print*,' End date               : ',trim(enddate)
334
  print*,' Reference time (year)  : ',reftime(1)
335
  print*,'                (month) : ',reftime(2)
336
  print*,'                (day)   : ',reftime(3)
337
  print*,'                (hour)  : ',reftime(4)
338
  print*,'                (min)   : ',reftime(5)
339
  print*,' Time range (min)       : ',reftime(6)
340
  do i=1,ncol
341
     print*,' Var                    :',i,trim(varsinp(i))
342
  enddo
343
  print*
344
 
345
  ! Check that first time is 0 - otherwise the tracing will produce
346
  ! wrong results because later in the code only absolute times are
347
  ! considered: <itime0   = int(abs(tfrac-tstart)/timeinc) + 1>. This
348
  ! will be changed in a future version.
349
  if ( abs( trainp(1,1,1) ).gt.eps ) then
350
     print*,' ERROR: First time of trajectory must be 0, i.e. '
351
     print*,'     correspond to the reference date. Otherwise'
352
     print*,'     the tracing will give wrong results... STOP'
353
     stop
354
  endif
355
 
356
  ! Save the full variable name (including shift specification)
357
  do i=1,ncol
358
     varsint(i)      = varsinp(i)
359
  enddo
360
  do i=1,ntrace0
361
     varsint(i+ncol) = tvar(i)
362
  enddo
363
 
364
  ! Split the variable name and set flags
365
  do i=1,ntrace0
366
 
367
     ! Set the prefix of the file name
368
     if ( compfl(i).eq.0  ) then
369
        tfil(i)=' '
370
        do j=1,n_pvars
371
           if ( tvar(i).eq.pvars(j) ) tfil(i)=charp
372
        enddo
373
        do j=1,n_svars
374
           if ( tvar(i).eq.svars(j) ) tfil(i)=chars
375
        enddo
376
     elseif ( tvar(i).eq.'Z' ) then
377
        tfil(i) = 'P'
378
     elseif ( tvar(i).eq.'ZB' ) then
379
        tfil(i) = 'P'
380
     else
381
        tfil(i) = '*' 
382
     endif
383
 
384
  enddo
385
 
386
  ! Write status information
387
  print*
388
  print*,'---- COMPLETE TABLE FOR TRACING -------------------------'
389
  print*
390
  do i=1,ntrace0
391
      write(*,'(i4,a4,a8,10x,8x,3x,a1,i5)') &
392
             i,' : ',trim(tvar(i)),tfil(i),compfl(i)
393
  enddo
394
 
395
  ! Allocate memory for internal trajectories
396
  allocate(traint(ntra,ntim,ncol+ntrace0),stat=stat)
397
  if (stat.ne.0) print*,'*** error allocating array traint   ***'
398
 
399
  ! Copy input to output trajectory
400
  do i=1,ntra
401
     do j=1,ntim
402
        do k=1,ncol
403
           traint(i,j,k)=trainp(i,j,k)
404
        enddo
405
     enddo
406
  enddo
407
 
408
  ! --------------------------------------------------------------------
409
  ! Trace the fields (fields available on input files)
410
  ! --------------------------------------------------------------------
411
 
412
  print*
413
  print*,'---- TRACING FROM PRIMARY AND SECONDARY DATA FILES ------'
414
 
415
  ! Loop over all tracing fields
416
  do i=1,ntrace0
417
 
418
      ! Write some status information
419
      print*
420
      print*,' Now tracing             : ', trim(tvar(i)),' ',trim(tfil(i))
421
 
422
      ! Reset flags for load manager
423
      iloaded0 = -1
424
      iloaded1 = -1
425
 
426
      ! Reset the counter for fields outside domain
427
      noutside = 0
428
 
429
      ! Loop over all times
430
      do j=1,ntim
431
 
432
         ! Convert trajectory time from hh.mm to fractional time
433
         call hhmm2frac(trainp(1,j,1),tfrac)
434
 
435
         ! Get the times which are needed
436
         itime0   = int(abs(tfrac-tstart)/timeinc) + 1
437
         time0    = tstart + fbflag * real(itime0-1) * timeinc
438
         itime1   = itime0 + 1
439
         time1    = time0 + fbflag * timeinc
440
         if ( itime1.gt.numdat ) then
441
            itime1 = itime0
442
            time1  = time0
443
         endif
444
 
445
         ! Load manager: Check whether itime0 can be copied from itime1
446
         if ( itime0.eq.iloaded1 ) then
447
            f3t0     = f3t1
448
            p3t0     = p3t1
449
            spt0     = spt1
450
            iloaded0 = itime0
451
         endif
452
 
453
         ! Load manager: Check whether itime1 can be copied from itime0
454
         if ( itime1.eq.iloaded0 ) then
455
            f3t1     = f3t0
456
            p3t1     = p3t0
457
            spt1     = spt0
458
            iloaded1 = itime1
459
         endif
460
 
461
 
462
         ! Load manager:  Load first time (tracing variable and grid)
463
         if ( itime0.ne.iloaded0 ) then
464
 
465
            filename = tfil(i)//dat(itime0)
466
            call frac2hhmm(time0,tload)
467
            varname  = tvar(i)
468
            write(*,'(a23,a20,a3,a5,f7.2)') '    ->  loading          : ', trim(filename),'  ',trim(varname),tload
469
            call input_open (fid,filename)
470
 
471
	        call input_var_wrf(fid,varname,f3t0,nx,ny,nz)	! field to trace
472
 
473
	        varname='geopot'					            ! geopot.height
474
	        call input_var_wrf(fid,varname,p3t0,nx,ny,nz)
475
 
476
	        varname='geopots'					            ! surface geopot.height
477
	        call input_var_wrf(fid,varname,spt0,nx,ny,1 )
478
 
479
            call input_close(fid)
480
 
481
            iloaded0 = itime0
482
 
483
         endif
484
 
485
         ! Load manager: Load second time (tracing variable and grid)
486
         if ( itime1.ne.iloaded1 ) then
487
 
488
            filename = tfil(i)//dat(itime1)
489
            call frac2hhmm(time1,tload)
490
            varname  = tvar(i)
491
            write(*,'(a23,a20,a3,a5,f7.2)') '    ->  loading          : ', trim(filename),'  ',trim(varname),tload
492
            call input_open (fid,filename)
493
 
494
            call input_var_wrf(fid,varname,f3t1,nx,ny,nz)	        ! field to trace
495
 
496
	        varname='geopot'					                    ! geopot.height
497
	        call input_var_wrf(fid,varname,p3t1,nx,ny,nz)
498
 
499
	        varname='geopots'					                    ! geopot.height
500
	        call input_var_wrf(fid,varname,spt1,nx,ny,1)
501
 
502
            call input_close(fid)
503
 
504
            iloaded1 = itime1
505
 
506
         endif
507
 
508
         ! Loop over all trajectories
509
         do k=1,ntra
510
 
511
            ! Set the horizontal position where to interpolate to
512
            x0       = traint(k,j,2)                          ! Longitude
513
            y0       = traint(k,j,3)                          ! Latitude
514
 
515
            ! Set the vertical position where to interpolate to: 
516
            if ( nz.gt.1 ) then
517
               p0 = traint(k,j,4)                            ! Pressure (3D tracing)
518
            else
519
               p0 = 1                                        ! Lowest level (2D tracing)
520
            endif
521
 
522
            ! Get the index where to interpolate (x0,y0,p0)
523
            if ( (abs(x0-mdv).gt.eps).and. (abs(y0-mdv).gt.eps) ) then
524
               call get_index4 (xind,yind,pind,x0,y0,p0,reltpos0, &
525
                    p3t0,p3t1,spt0,spt1,3,nx,ny,nz,xmin,ymin,dx,dy,mdv)
526
            else
527
               xind = mdv
528
               yind = mdv
529
               pind = mdv
530
            endif
531
 
532
 
533
            ! Check if point is within grid (keep indices if ok)
534
            if ( (xind.ge.1.).and.(xind.le.real(nx)).and. &
535
                 (yind.ge.1.).and.(yind.le.real(ny)).and. &
536
                 (pind.ge.1.).and.(pind.le.real(nz)) )    &
537
            then
538
		           xind = xind
539
		           yind = yind
540
		           pind = pind
541
 
542
            endif
543
 
544
            ! ------------------------ NEAREST mode ------------------------------- 
545
 
546
            if ( intmode.eq.'nearest') then
547
 
548
               xind = real( nint(xind) )
549
               yind = real( nint(yind) )
550
               pind = real( nint(pind) )
551
 
552
               if ( xind.lt.1.  ) xind = 1.
553
               if ( xind.gt.nx  ) xind = real(nx)
554
               if ( yind.lt.1.  ) yind = 1.
555
               if ( yind.gt.ny  ) yind = real(ny)
556
 
557
               if ( pind.lt.1.  ) pind = 1.
558
               if ( pind.gt.nz  ) pind = real(nz)
559
 
560
               if ( abs(reltpos0).ge.eps ) then
561
                  print*,'ERROR (nearest): reltpos != 0',reltpos0
562
                  stop
563
               endif
564
 
565
               ! interpolate
566
               f0 = int_index4(f3t0,f3t1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
567
 
568
            ! ------------------------ NORMAL mode -------------------------------
569
            else 
570
 
571
               f0 = int_index4(f3t0,f3t1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
572
 
573
            endif
574
 
575
            ! Save the new field
576
            if ( abs(f0-mdv).gt.eps) then
577
               traint(k,j,ncol+i) = f0
578
            else
579
               traint(k,j,ncol+i) = mdv
580
            endif
581
 
582
         enddo ! end loop over all trajectories
583
 
584
      enddo ! end loop over all times
585
 
586
      ! Exit point for loop over all tracing variables
587
 110   continue
588
 
589
  enddo ! end loop over all tracing varoables
590
 
591
  ! --------------------------------------------------------------------
592
  ! Write output to output trajectory file
593
  ! --------------------------------------------------------------------
594
 
595
  ! Write status information
596
  print*
597
  print*,'---- WRITE OUTPUT TRAJECTORIES --------------------------'
598
  print*
599
 
600
  ! Allocate memory for internal trajectories
601
  allocate(traout(ntra,ntim,ncol+ntrace0),stat=stat)
602
  if (stat.ne.0) print*,'*** error allocating array traout   ***'
603
 
604
  ! Copy input to output trajectory (apply scaling of output)
605
  do i=1,ntra
606
     do j=1,ntim
607
        do k=1,ncol+ntrace0
608
           if ( k.le.ncol ) then
609
              traout(i,j,k) = traint(i,j,k)
610
           elseif ( abs(traint(i,j,k)-mdv).gt.eps ) then
611
              traout(i,j,k) = fac(k-ncol) * traint(i,j,k)
612
           else
613
              traout(i,j,k) = mdv
614
           endif
615
        enddo
616
     enddo
617
  enddo
618
 
619
  ! Set the variable names for output trajectory
620
  do i=1,ncol+ntrace0
621
     varsout(i)      = varsint(i)
622
  enddo
623
 
624
 ! Transform coordinates into index space
625
   do i=1,ntra
626
     do j=1,ntim
627
         traout(i,j,2) = ( traout(i,j,2) - xmin ) / dx + 1.
628
         traout(i,j,3) = ( traout(i,j,3) - ymin ) / dy + 1.
629
     enddo
630
  enddo
631
 
632
  ! Write trajectories
633
  call wopen_tra(fod,outfile,ntra,ntim,ncol+ntrace0, &
634
       reftime,varsout,outmode)
635
  call write_tra(fod,traout ,ntra,ntim,ncol+ntrace0,outmode)
636
  call close_tra(fod,outmode)
637
 
638
  ! Write some status information, and end of program message
639
  print*
640
  print*,'---- STATUS INFORMATION --------------------------------'
641
  print*
642
  print*,' ok'
643
  print*
644
  print*,'              *** END OF PROGRAM TRACE ***'
645
  print*,'========================================================='
646
 
647
 
648
end program trace_label_circle
649