Subversion Repositories lagranto.wrf

Rev

Rev 15 | Rev 19 | Go to most recent revision | Only display areas with differences | Ignore whitespace | Details | Blame | Last modification | View Log | RSS feed

Rev 15 Rev 17
1
      PROGRAM caltra
1
      PROGRAM caltra
2
 
2
 
3
C     ********************************************************************
3
C     ********************************************************************
4
C     *                                                                  *
4
C     *                                                                  *
5
C     * Calculates trajectories                                          *
5
C     * Calculates trajectories                                          *
6
C     *                                                                  *
6
C     *                                                                  *
7
C     *	Heini Wernli	   first version:	April 1993                   *
7
C     *	Heini Wernli	   first version:	April 1993                   *
8
C     * Michael Sprenger   major upgrade:       2008-2009                *
8
C     * Michael Sprenger   major upgrade:       2008-2009                *
9
C     *                                                                  *
9
C     *                                                                  *
10
C     ********************************************************************
10
C     ********************************************************************
11
 
11
 
12
      implicit none      
12
      implicit none      
13
 
13
 
14
c     --------------------------------------------------------------------
14
c     --------------------------------------------------------------------
15
c     Declaration of parameters
15
c     Declaration of parameters
16
c     --------------------------------------------------------------------
16
c     --------------------------------------------------------------------
17
 
17
 
18
c     Maximum number of levels for input files
18
c     Maximum number of levels for input files
19
      integer   nlevmax
19
      integer   nlevmax
20
      parameter	(nlevmax=100)
20
      parameter	(nlevmax=100)
21
 
21
 
22
c     Maximum number of input files (dates, length of trajectories)
22
c     Maximum number of input files (dates, length of trajectories)
23
      integer   ndatmax
23
      integer   ndatmax
24
      parameter	(ndatmax=500)
24
      parameter	(ndatmax=500)
25
 
25
 
26
c     Numerical epsilon (for float comparison)
26
c     Numerical epsilon (for float comparison)
27
      real      eps
27
      real      eps
28
      parameter (eps=0.001)
28
      parameter (eps=0.001)
29
 
29
 
30
c     Distance in m between 2 lat circles 
30
c     Distance in m between 2 lat circles 
31
      real	deltay
31
      real	deltay
32
      parameter	(deltay=1.112E5)
32
      parameter	(deltay=1.112E5)
33
 
33
 
34
c     Number of iterations for iterative Euler scheme
34
c     Number of iterations for iterative Euler scheme
35
      integer   numit
35
      integer   numit
36
      parameter (numit=3)
36
      parameter (numit=3)
37
 
37
 
38
c     Input and output format for trajectories (see iotra.f)
38
c     Input and output format for trajectories (see iotra.f)
39
      integer   inpmode
39
      integer   inpmode
40
      integer   outmode
40
      integer   outmode
41
 
41
 
-
 
42
c     Missing data value
-
 
43
      real      mdv
-
 
44
      parameter (mdv = -999. )
-
 
45
 
42
c     Filename prefix (typically 'P')
46
c     Filename prefix (typically 'P')
43
      character*1 prefix
47
      character*1 prefix
44
      parameter   (prefix='P')
48
      parameter   (prefix='P')
45
 
49
 
46
c     Externals
50
c     Externals
47
      real       sdis
51
      real       sdis
48
      external   sdis
52
      external   sdis
49
 
53
 
50
c     --------------------------------------------------------------------
54
c     --------------------------------------------------------------------
51
c     Declaration of variables
55
c     Declaration of variables
52
c     --------------------------------------------------------------------
56
c     --------------------------------------------------------------------
53
 
57
 
54
c     Input parameters
58
c     Input parameters
55
      integer                                fbflag          ! Flag for forward/backward mode
59
      integer                                fbflag          ! Flag for forward/backward mode
56
      integer                                numdat          ! Number of input files
60
      integer                                numdat          ! Number of input files
57
      character*11                           dat(ndatmax)    ! Dates of input files
61
      character*11                           dat(ndatmax)    ! Dates of input files
58
      real                                   timeinc         ! Time increment between input files
62
      real                                   timeinc         ! Time increment between input files
59
      integer                                per             ! Periodicity flag (=0 if none)
63
      integer                                per             ! Periodicity flag (=0 if none)
60
      integer                                ntra            ! Number of trajectories
64
      integer                                ntra            ! Number of trajectories
61
      character*80                           cdfname         ! Name of output files
65
      character*80                           cdfname         ! Name of output files
62
      real                                   ts              ! Time step
66
      real                                   ts              ! Time step
63
      real                                   tst,ten         ! Shift of start and end time relative to first data file
67
      real                                   tst,ten         ! Shift of start and end time relative to first data file
64
      integer                                deltout         ! Output time interval (in minutes)
68
      integer                                deltout         ! Output time interval (in minutes)
65
      integer                                jflag           ! Jump flag (if =1 ground-touching trajectories reenter atmosphere)
69
      integer                                jflag           ! Jump flag (if =1 ground-touching trajectories reenter atmosphere)
66
      real                                   wfactor         ! Factor for vertical velocity field
70
      real                                   wfactor         ! Factor for vertical velocity field
67
      character*80                           strname         ! File with start positions
71
      character*80                           strname         ! File with start positions
68
      character*80                           timecheck       ! Either 'yes' or 'no'
72
      character*80                           timecheck       ! Either 'yes' or 'no'
69
      character*10		              	     gridflag        ! "ideal" or "real"
73
      character*10		              	     gridflag        ! "ideal" or "real"
70
 
74
 
71
c     Trajectories
75
c     Trajectories
72
      integer                                ncol            ! Number of columns for insput trajectories
76
      integer                                ncol            ! Number of columns for insput trajectories
73
      real,allocatable, dimension (:,:,:) :: trainp          ! Input start coordinates (ntra,1,ncol)
77
      real,allocatable, dimension (:,:,:) :: trainp          ! Input start coordinates (ntra,1,ncol)
74
      real,allocatable, dimension (:,:,:) :: traout          ! Output trajectories (ntra,ntim,4)
78
      real,allocatable, dimension (:,:,:) :: traout          ! Output trajectories (ntra,ntim,4)
75
      integer                                reftime(6)      ! Reference date
79
      integer                                reftime(6)      ! Reference date
76
      character*80                           vars(200)       ! Field names
80
      character*80                           vars(200)       ! Field names
77
      real,allocatable, dimension (:)     :: xx0,yy0,pp0     ! Position of air parcels
81
      real,allocatable, dimension (:)     :: xx0,yy0,pp0     ! Position of air parcels
78
      integer,allocatable, dimension (:)  :: leftflag        ! Flag for domain-leaving
82
      integer,allocatable, dimension (:)  :: leftflag        ! Flag for domain-leaving
79
      real                                   xx1,yy1,pp1     ! Updated position of air parcel
83
      real                                   xx1,yy1,pp1     ! Updated position of air parcel
80
      integer                                leftcount       ! Number of domain leaving trajectories
84
      integer                                leftcount       ! Number of domain leaving trajectories
81
      integer                                ntim            ! Number of output time steps
85
      integer                                ntim            ! Number of output time steps
82
 
86
 
83
c     Meteorological fields
87
c     Meteorological fields
84
      real,allocatable, dimension (:)     :: spt0,spt1       ! Surface pressure
88
      real,allocatable, dimension (:)     :: spt0,spt1       ! Surface pressure
85
      real,allocatable, dimension (:)     :: uut0,uut1       ! Zonal wind
89
      real,allocatable, dimension (:)     :: uut0,uut1       ! Zonal wind
86
      real,allocatable, dimension (:)     :: vvt0,vvt1       ! Meridional wind
90
      real,allocatable, dimension (:)     :: vvt0,vvt1       ! Meridional wind
87
      real,allocatable, dimension (:)     :: wwt0,wwt1       ! Vertical wind
91
      real,allocatable, dimension (:)     :: wwt0,wwt1       ! Vertical wind
88
      real,allocatable, dimension (:)     :: p3t0,p3t1       ! 3d-pressure 
92
      real,allocatable, dimension (:)     :: p3t0,p3t1       ! 3d-pressure 
89
 
93
 
90
c     Grid description
94
c     Grid description
91
      real                                   pollon,pollat   ! Longitude/latitude of pole
95
      real                                   pollon,pollat   ! Longitude/latitude of pole
92
      real                                   ak(nlevmax)     ! Vertical layers and levels
96
      real                                   ak(nlevmax)     ! Vertical layers and levels
93
      real                                   bk(nlevmax) 
97
      real                                   bk(nlevmax) 
94
      real                                   xmin,xmax       ! Zonal grid extension
98
      real                                   xmin,xmax       ! Zonal grid extension
95
      real                                   ymin,ymax       ! Meridional grid extension
99
      real                                   ymin,ymax       ! Meridional grid extension
96
      integer                                nx,ny,nz        ! Grid dimensions
100
      integer                                nx,ny,nz        ! Grid dimensions
97
      real                                   dx,dy           ! Horizontal grid resolution
101
      real                                   dx,dy           ! Horizontal grid resolution
98
      integer                                hem             ! Flag for hemispheric domain
102
      integer                                hem             ! Flag for hemispheric domain
99
      real                                   mdv             ! Missing data value
-
 
100
      real,allocatable, dimension (:,:)   :: mpx,mpy         ! Map scale factor in x,y direction
103
      real,allocatable, dimension (:,:)   :: mpx,mpy         ! Map scale factor in x,y direction
101
      real,allocatable, dimension (:,:)   :: lat,lon         ! Map scale factor in x,y direction
104
      real,allocatable, dimension (:,:)   :: lat,lon         ! Map scale factor in x,y direction
102
 
105
 
103
c     Auxiliary variables                 
106
c     Auxiliary variables                 
104
      real                                   delta,rd
107
      real                                   delta,rd
105
      integer	                             itm,iloop,i,j,k,filo,lalo
108
      integer	                             itm,iloop,i,j,k,filo,lalo
106
      integer                                ierr,stat
109
      integer                                ierr,stat
107
      integer                                cdfid,fid
110
      integer                                cdfid,fid
108
      real	                                 tstart,time0,time1,time
111
      real	                                 tstart,time0,time1,time
109
      real                                   reltpos0,reltpos1
112
      real                                   reltpos0,reltpos1
110
      real                                   xind,yind,pind,pp,sp,stagz
113
      real                                   xind,yind,pind,pp,sp,stagz
111
      character*80                           filename,varname
114
      character*80                           filename,varname
112
      integer                                reftmp(6)
115
      integer                                reftmp(6)
113
      character                              ch
116
      character                              ch
114
      real                                   frac,tload
117
      real                                   frac,tload
115
      integer                                itim
118
      integer                                itim
116
      integer                                wstep
119
      integer                                wstep
117
      character*80                           radunit
120
      character*80                           radunit
118
      real                                   lon1,lat1,lon2,lat2
121
      real                                   lon1,lat1,lon2,lat2
119
      real                                   scale_max,scale_min
122
      real                                   scale_max,scale_min
120
 
123
 
121
c     --------------------------------------------------------------------
124
c     --------------------------------------------------------------------
122
c     Start of program, Read parameters
125
c     Start of program, Read parameters
123
c     --------------------------------------------------------------------
126
c     --------------------------------------------------------------------
124
 
127
 
125
c     Write start message
128
c     Write start message
126
      print*,'========================================================='
129
      print*,'========================================================='
127
      print*,'              *** START OF PROGRAM CALTRA ***'
130
      print*,'              *** START OF PROGRAM CALTRA ***'
128
      print*
131
      print*
129
 
132
 
130
c     Open the parameter file
133
c     Open the parameter file
131
      open(9,file='caltra.param')
134
      open(9,file='caltra.param')
132
 
135
 
133
c     Read flag for forward/backward mode (fbflag)
136
c     Read flag for forward/backward mode (fbflag)
134
      read(9,*) fbflag
137
      read(9,*) fbflag
135
 
138
 
136
c     Read number of input files (numdat)
139
c     Read number of input files (numdat)
137
      read(9,*) numdat
140
      read(9,*) numdat
138
      if (numdat.gt.ndatmax) then
141
      if (numdat.gt.ndatmax) then
139
        print*,' ERROR: too many input files ',numdat,ndatmax
142
        print*,' ERROR: too many input files ',numdat,ndatmax
140
        goto 993
143
        goto 993
141
      endif
144
      endif
142
 
145
 
143
c     Read list of input dates (dat, sort depending on forward/backward mode)
146
c     Read list of input dates (dat, sort depending on forward/backward mode)
144
      if (fbflag.eq.1) then
147
      if (fbflag.eq.1) then
145
        do itm=1,numdat
148
        do itm=1,numdat
146
          read(9,'(a11)') dat(itm)
149
          read(9,'(a11)') dat(itm)
147
        enddo
150
        enddo
148
      else
151
      else
149
        do itm=numdat,1,-1
152
        do itm=numdat,1,-1
150
          read(9,'(a11)') dat(itm)
153
          read(9,'(a11)') dat(itm)
151
        enddo
154
        enddo
152
      endif
155
      endif
153
 
156
 
154
c     Read time increment between input files (timeinc)
157
c     Read time increment between input files (timeinc)
155
      read(9,*) timeinc
158
      read(9,*) timeinc
156
 
159
 
157
C     Read if data domain is periodic and its periodicity
160
C     Read if data domain is periodic and its periodicity
158
      read(9,*) per
161
      read(9,*) per
159
 
162
 
160
 
163
 
161
c     Read the number of trajectories and name of position file
164
c     Read the number of trajectories and name of position file
162
      read(9,*) strname
165
      read(9,*) strname
163
      read(9,*) ntra
166
      read(9,*) ntra
164
      read(9,*) ncol 
167
      read(9,*) ncol 
165
      if (ntra.eq.0) goto 991
168
      if (ntra.eq.0) goto 991
166
 
169
 
167
C     Read the name of the output trajectory file and set the constants file
170
C     Read the name of the output trajectory file and set the constants file
168
      read(9,*) cdfname
171
      read(9,*) cdfname
169
 
172
 
170
C     Read the timestep for trajectory calculation (convert from minutes to hours)
173
C     Read the timestep for trajectory calculation (convert from minutes to hours)
171
      read(9,*) ts
174
      read(9,*) ts
172
      ts=ts/60.    
175
      ts=ts/60.    
173
 
176
 
174
C     Read shift of start and end time relative to first data file
177
C     Read shift of start and end time relative to first data file
175
      read(9,*) tst
178
      read(9,*) tst
176
      read(9,*) ten
179
      read(9,*) ten
177
 
180
 
178
C     Read output time interval (in minutes)
181
C     Read output time interval (in minutes)
179
      read(9,*) deltout
182
      read(9,*) deltout
180
 
183
 
181
C     Read jumpflag (if =1 ground-touching trajectories reenter the atmosphere)
184
C     Read jumpflag (if =1 ground-touching trajectories reenter the atmosphere)
182
      read(9,*) jflag
185
      read(9,*) jflag
183
 
186
 
184
C     Read factor for vertical velocity field
187
C     Read factor for vertical velocity field
185
      read(9,*) wfactor
188
      read(9,*) wfactor
186
 
189
 
187
c     Read the reference time and the time range
190
c     Read the reference time and the time range
188
      read(9,*) reftime(1)                  ! year
191
      read(9,*) reftime(1)                  ! year
189
      read(9,*) reftime(2)                  ! month
192
      read(9,*) reftime(2)                  ! month
190
      read(9,*) reftime(3)                  ! day
193
      read(9,*) reftime(3)                  ! day
191
      read(9,*) reftime(4)                  ! hour
194
      read(9,*) reftime(4)                  ! hour
192
      read(9,*) reftime(5)                  ! min
195
      read(9,*) reftime(5)                  ! min
193
      read(9,*) reftime(6)                  ! time range (in min)
196
      read(9,*) reftime(6)                  ! time range (in min)
194
 
197
 
195
c     Read flag for 'no time check'
198
c     Read flag for 'no time check'
196
      read(9,*) timecheck
199
      read(9,*) timecheck
197
 
200
 
198
c     Close the input file
201
c     Close the input file
199
      close(9)
202
      close(9)
200
 
203
 
201
c     Calculate the number of output time steps
204
c     Calculate the number of output time steps
202
      ntim = abs(reftime(6)/deltout) + 1
205
      ntim = abs(reftime(6)/deltout) + 1
203
 
206
 
204
c     Set the formats of the input and output files
207
c     Set the formats of the input and output files
205
      call mode_tra(inpmode,strname)
208
      call mode_tra(inpmode,strname)
206
      call mode_tra(outmode,cdfname)
209
      call mode_tra(outmode,cdfname)
207
      if (outmode.eq.-1) outmode=1
210
      if (outmode.eq.-1) outmode=1
208
 
211
 
209
c     Write some status information
212
c     Write some status information
210
      print*,'---- INPUT PARAMETERS -----------------------------------'
213
      print*,'---- INPUT PARAMETERS -----------------------------------'
211
      print* 
214
      print* 
212
      print*,'  Forward/Backward       : ',fbflag
215
      print*,'  Forward/Backward       : ',fbflag
213
      print*,'  #input files           : ',numdat
216
      print*,'  #input files           : ',numdat
214
      print*,'  First/last input file  : ',trim(dat(1)),' ... ',
217
      print*,'  First/last input file  : ',trim(dat(1)),' ... ',
215
     >                                     trim(dat(numdat))
218
     >                                     trim(dat(numdat))
216
      print*,'  time increment         : ',timeinc
219
      print*,'  time increment         : ',timeinc
217
      print*,'  Output file            : ',trim(cdfname)
220
      print*,'  Output file            : ',trim(cdfname)
218
      print*,'  Time step (min)        : ',60.*ts
221
      print*,'  Time step (min)        : ',60.*ts
219
      write(*,'(a27,f7.2,f7.2)') '   Time shift (start,end) : ',tst,ten
222
      write(*,'(a27,f7.2,f7.2)') '   Time shift (start,end) : ',tst,ten
220
      print*,'  Output time interval   : ',deltout
223
      print*,'  Output time interval   : ',deltout
221
      print*,'  Jump flag              : ',jflag
224
      print*,'  Jump flag              : ',jflag
222
      print*,'  Vertical wind (scale)  : ',wfactor
225
      print*,'  Vertical wind (scale)  : ',wfactor
223
      print*,'  Trajectory pos file    : ',trim(strname)
226
      print*,'  Trajectory pos file    : ',trim(strname)
224
      print*,'  # of trajectories      : ',ntra
227
      print*,'  # of trajectories      : ',ntra
225
      print*,'  # of output timesteps  : ',ntim
228
      print*,'  # of output timesteps  : ',ntim
226
      if ( inpmode.eq.-1) then
229
      if ( inpmode.eq.-1) then
227
         print*,'  Input format           : (x,y,z)-list'
230
         print*,'  Input format           : (x,y,z)-list'
228
      else
231
      else
229
         print*,'  Input format           : ',inpmode
232
         print*,'  Input format           : ',inpmode
230
      endif
233
      endif
231
      print*,'  Output format          : ',outmode
234
      print*,'  Output format          : ',outmode
232
      print*,'  Periodicity            : ',per
235
      print*,'  Periodicity            : ',per
233
      print*,'  Time check             : ',trim(timecheck)
236
      print*,'  Time check             : ',trim(timecheck)
234
      print*
237
      print*
235
 
238
 
236
      print*,'---- FIXED NUMERICAL PARAMETERS -------------------------'
239
      print*,'---- FIXED NUMERICAL PARAMETERS -------------------------'
237
      print*
240
      print*
238
      print*,'  Number of iterations   : ',numit
241
      print*,'  Number of iterations   : ',numit
239
      print*,'  Filename prefix        : ',prefix
242
      print*,'  Filename prefix        : ',prefix
240
      print*,'  Missing data value     : ',mdv
243
      print*,'  Missing data value     : ',mdv
241
      print*
244
      print*
242
 
245
 
243
c     --------------------------------------------------------------------
246
c     --------------------------------------------------------------------
244
c     Read grid parameters, checks and allocate memory
247
c     Read grid parameters, checks and allocate memory
245
c     --------------------------------------------------------------------
248
c     --------------------------------------------------------------------
246
 
249
 
247
c     Read the constant grid parameters (nx,ny,nz,xmin,xmax,ymin,ymax,pollon,pollat)
250
c     Read the constant grid parameters (nx,ny,nz,xmin,xmax,ymin,ymax,pollon,pollat)
248
c     The negative <-fid> of the file identifier is used as a flag for parameter retrieval  
251
c     The negative <-fid> of the file identifier is used as a flag for parameter retrieval  
249
      filename = prefix//dat(1)
252
      filename = prefix//dat(1)
250
      varname  = 'U'
253
      varname  = 'U'
251
      nx       = 1
254
      nx       = 1
252
      ny       = 1
255
      ny       = 1
253
      nz       = 1
256
      nz       = 1
254
      tload    = -tst
257
      tload    = -tst
255
      call input_open (fid,filename)
258
      call input_open (fid,filename)
256
	  call input_grid_wrf(fid,xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz)
259
	  call input_grid_wrf(fid,xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz)
257
      call input_close(fid)
260
      call input_close(fid)
258
 
261
 
259
C     Check if the number of levels is too large
262
C     Check if the number of levels is too large
260
      if (nz.gt.nlevmax) goto 993
263
      if (nz.gt.nlevmax) goto 993
261
 
264
 
262
C     Set logical flag for hemispheric mode (not yet supported)
265
C     Set logical flag for hemispheric mode (not yet supported)
263
      hem = 0
266
      hem = 0
264
 
267
 
265
C     Allocate memory for some meteorological arrays
268
C     Allocate memory for some meteorological arrays
266
      allocate(spt0(nx*ny),stat=stat)
269
      allocate(spt0(nx*ny),stat=stat)
267
      if (stat.ne.0) print*,'*** error allocating array spt0 ***'   ! Surface geopotential height
270
      if (stat.ne.0) print*,'*** error allocating array spt0 ***'   ! Surface geopotential height
268
      allocate(spt1(nx*ny),stat=stat)
271
      allocate(spt1(nx*ny),stat=stat)
269
      if (stat.ne.0) print*,'*** error allocating array spt1 ***'
272
      if (stat.ne.0) print*,'*** error allocating array spt1 ***'
270
      allocate(uut0(nx*ny*nz),stat=stat)
273
      allocate(uut0(nx*ny*nz),stat=stat)
271
      if (stat.ne.0) print*,'*** error allocating array uut0 ***'   ! Zonal wind
274
      if (stat.ne.0) print*,'*** error allocating array uut0 ***'   ! Zonal wind
272
      allocate(uut1(nx*ny*nz),stat=stat)
275
      allocate(uut1(nx*ny*nz),stat=stat)
273
      if (stat.ne.0) print*,'*** error allocating array uut1 ***'
276
      if (stat.ne.0) print*,'*** error allocating array uut1 ***'
274
      allocate(vvt0(nx*ny*nz),stat=stat)
277
      allocate(vvt0(nx*ny*nz),stat=stat)
275
      if (stat.ne.0) print*,'*** error allocating array vvt0 ***'   ! Meridional wind
278
      if (stat.ne.0) print*,'*** error allocating array vvt0 ***'   ! Meridional wind
276
      allocate(vvt1(nx*ny*nz),stat=stat)
279
      allocate(vvt1(nx*ny*nz),stat=stat)
277
      if (stat.ne.0) print*,'*** error allocating array vvt1 ***'
280
      if (stat.ne.0) print*,'*** error allocating array vvt1 ***'
278
      allocate(wwt0(nx*ny*nz),stat=stat)
281
      allocate(wwt0(nx*ny*nz),stat=stat)
279
      if (stat.ne.0) print*,'*** error allocating array wwt0 ***'   ! Vertical wind
282
      if (stat.ne.0) print*,'*** error allocating array wwt0 ***'   ! Vertical wind
280
      allocate(wwt1(nx*ny*nz),stat=stat)
283
      allocate(wwt1(nx*ny*nz),stat=stat)
281
      if (stat.ne.0) print*,'*** error allocating array wwt1 ***'
284
      if (stat.ne.0) print*,'*** error allocating array wwt1 ***'
282
      allocate(p3t0(nx*ny*nz),stat=stat)
285
      allocate(p3t0(nx*ny*nz),stat=stat)
283
      if (stat.ne.0) print*,'*** error allocating array p3t0 ***'   ! geopotential height
286
      if (stat.ne.0) print*,'*** error allocating array p3t0 ***'   ! geopotential height
284
      allocate(p3t1(nx*ny*nz),stat=stat)
287
      allocate(p3t1(nx*ny*nz),stat=stat)
285
      if (stat.ne.0) print*,'*** error allocating array p3t1 ***'
288
      if (stat.ne.0) print*,'*** error allocating array p3t1 ***'
286
      allocate(mpx(nx,ny),stat=stat)
289
      allocate(mpx(nx,ny),stat=stat)
287
      if (stat.ne.0) print*,'*** error allocating array mpx * **'   ! Map scale factor in x
290
      if (stat.ne.0) print*,'*** error allocating array mpx * **'   ! Map scale factor in x
288
      allocate(mpy(nx,ny),stat=stat)
291
      allocate(mpy(nx,ny),stat=stat)
289
      if (stat.ne.0) print*,'*** error allocating array mpy ***'    ! Map scale factor in y
292
      if (stat.ne.0) print*,'*** error allocating array mpy ***'    ! Map scale factor in y
290
      allocate(lon(nx,ny),stat=stat)
293
      allocate(lon(nx,ny),stat=stat)
291
      if (stat.ne.0) print*,'*** error allocating array lon ***'    ! Grid -> lon
294
      if (stat.ne.0) print*,'*** error allocating array lon ***'    ! Grid -> lon
292
      allocate(lat(nx,ny),stat=stat)
295
      allocate(lat(nx,ny),stat=stat)
293
      if (stat.ne.0) print*,'*** error allocating array lat ***'    ! Gridy -> lat
296
      if (stat.ne.0) print*,'*** error allocating array lat ***'    ! Gridy -> lat
294
 
297
 
295
C     Get memory for trajectory arrays
298
C     Get memory for trajectory arrays
296
      allocate(trainp(ntra,1,ncol),stat=stat)
299
      allocate(trainp(ntra,1,ncol),stat=stat)
297
      if (stat.ne.0) print*,'*** error allocating array trainp   ***' ! Input start coordinates
300
      if (stat.ne.0) print*,'*** error allocating array trainp   ***' ! Input start coordinates
298
      allocate(traout(ntra,ntim,4),stat=stat)
301
      allocate(traout(ntra,ntim,4),stat=stat)
299
      if (stat.ne.0) print*,'*** error allocating array traout   ***' ! Output trajectories
302
      if (stat.ne.0) print*,'*** error allocating array traout   ***' ! Output trajectories
300
      allocate(xx0(ntra),stat=stat)
303
      allocate(xx0(ntra),stat=stat)
301
      if (stat.ne.0) print*,'*** error allocating array xx0      ***' ! X position (longitude)
304
      if (stat.ne.0) print*,'*** error allocating array xx0      ***' ! X position (longitude)
302
      allocate(yy0(ntra),stat=stat)
305
      allocate(yy0(ntra),stat=stat)
303
      if (stat.ne.0) print*,'*** error allocating array yy0      ***' ! Y position (latitude)
306
      if (stat.ne.0) print*,'*** error allocating array yy0      ***' ! Y position (latitude)
304
      allocate(pp0(ntra),stat=stat)
307
      allocate(pp0(ntra),stat=stat)
305
      if (stat.ne.0) print*,'*** error allocating array pp0      ***' ! Pressure
308
      if (stat.ne.0) print*,'*** error allocating array pp0      ***' ! Pressure
306
      allocate(leftflag(ntra),stat=stat)
309
      allocate(leftflag(ntra),stat=stat)
307
      if (stat.ne.0) print*,'*** error allocating array leftflag ***' ! Leaving-domain flag
310
      if (stat.ne.0) print*,'*** error allocating array leftflag ***' ! Leaving-domain flag
308
 
311
 
309
c     Write some status information
312
c     Write some status information
310
      print*,'---- CONSTANT GRID PARAMETERS ---------------------------'
313
      print*,'---- CONSTANT GRID PARAMETERS ---------------------------'
311
      print*
314
      print*
312
      print*,'  xmin,xmax     : ',xmin,xmax
315
      print*,'  xmin,xmax     : ',xmin,xmax
313
      print*,'  ymin,ymax     : ',ymin,ymax
316
      print*,'  ymin,ymax     : ',ymin,ymax
314
      print*,'  dx,dy         : ',dx,dy
317
      print*,'  dx,dy         : ',dx,dy
315
      print*,'  nx,ny,nz      : ',nx,ny,nz
318
      print*,'  nx,ny,nz      : ',nx,ny,nz
316
      print*,'  per, hem      : ',per,hem
319
      print*,'  per, hem      : ',per,hem
317
      print*
320
      print*
318
 
321
 
319
c     --------------------------------------------------------------------
322
c     --------------------------------------------------------------------
320
c     Calculate the map scale factors
323
c     Calculate the map scale factors
321
c     --------------------------------------------------------------------
324
c     --------------------------------------------------------------------
322
 
325
 
323
c     Read lon/lat on the model grid from first data file
326
c     Read lon/lat on the model grid from first data file
324
      filename = prefix//dat(1)
327
      filename = prefix//dat(1)
325
      call input_open (fid,filename)
328
      call input_open (fid,filename)
326
      varname='XLONG'
329
      varname='XLONG'
327
      call input_var_wrf(fid,varname,lon,nx,ny,1)
330
      call input_var_wrf(fid,varname,lon,nx,ny,1)
328
      varname='XLAT'
331
      varname='XLAT'
329
      call input_var_wrf(fid,varname,lat,nx,ny,1)
332
      call input_var_wrf(fid,varname,lat,nx,ny,1)
330
      call input_close(fid)
333
      call input_close(fid)
331
 
334
 
332
c     Get map scale for interior points
335
c     Get map scale for interior points
333
      do i=2,nx-1
336
      do i=2,nx-1
334
         do j=2,ny-1
337
         do j=2,ny-1
335
 
338
 
336
c           Map scale in x direction
339
c           Map scale in x direction
337
            lon1     = lon(i-1,j)
340
            lon1     = lon(i-1,j)
338
            lat1     = lat(i-1,j)
341
            lat1     = lat(i-1,j)
339
            lon2     = lon(i+1,j)
342
            lon2     = lon(i+1,j)
340
            lat2     = lat(i+1,j)
343
            lat2     = lat(i+1,j)
341
            radunit  = 'km.haversine'
344
            radunit  = 'km.haversine'
342
            mpx(i,j) = 0.5 * 1000. * sdis(lon1,lat1,lon2,lat2,radunit)
345
            mpx(i,j) = 0.5 * 1000. * sdis(lon1,lat1,lon2,lat2,radunit)
343
 
346
 
344
c           Map scale in y direction
347
c           Map scale in y direction
345
            lon1     = lon(i,j-1)
348
            lon1     = lon(i,j-1)
346
            lat1     = lat(i,j-1)
349
            lat1     = lat(i,j-1)
347
            lon2     = lon(i,j+1)
350
            lon2     = lon(i,j+1)
348
            lat2     = lat(i,j+1)
351
            lat2     = lat(i,j+1)
349
            radunit  = 'km.haversine'
352
            radunit  = 'km.haversine'
350
            mpy(i,j) = 0.5 * 1000. * sdis(lon1,lat1,lon2,lat2,radunit)
353
            mpy(i,j) = 0.5 * 1000. * sdis(lon1,lat1,lon2,lat2,radunit)
351
 
354
 
352
         enddo
355
         enddo
353
      enddo
356
      enddo
354
 
357
 
355
c     Copy map scale for boundary points
358
c     Copy map scale for boundary points
356
      do i=2,nx-1
359
      do i=2,nx-1
357
        mpx(i, 1) = mpx(i,   2)
360
        mpx(i, 1) = mpx(i,   2)
358
        mpx(i,ny) = mpx(i,ny-1)
361
        mpx(i,ny) = mpx(i,ny-1)
359
        mpy(i, 1) = mpy(i,   2)
362
        mpy(i, 1) = mpy(i,   2)
360
        mpy(i,ny) = mpy(i,ny-1)
363
        mpy(i,ny) = mpy(i,ny-1)
361
      enddo
364
      enddo
362
      do j=2,ny-1
365
      do j=2,ny-1
363
        mpx(1, j) = mpx(2,   j)
366
        mpx(1, j) = mpx(2,   j)
364
        mpx(nx,j) = mpx(nx-1,j)
367
        mpx(nx,j) = mpx(nx-1,j)
365
        mpy(1, j) = mpy(2,   j)
368
        mpy(1, j) = mpy(2,   j)
366
        mpy(nx,j) = mpy(nx-1,j)
369
        mpy(nx,j) = mpy(nx-1,j)
367
      enddo
370
      enddo
368
      mpx(1,1)   = mpx(2,2)
371
      mpx(1,1)   = mpx(2,2)
369
      mpy(1,1)   = mpy(2,2)
372
      mpy(1,1)   = mpy(2,2)
370
      mpx(nx,1)  = mpx(nx-1,2)
373
      mpx(nx,1)  = mpx(nx-1,2)
371
      mpy(nx,1)  = mpy(nx-1,2)
374
      mpy(nx,1)  = mpy(nx-1,2)
372
      mpx(nx,ny) = mpx(nx-1,ny-1)
375
      mpx(nx,ny) = mpx(nx-1,ny-1)
373
      mpy(nx,ny) = mpy(nx-1,ny-1)
376
      mpy(nx,ny) = mpy(nx-1,ny-1)
374
      mpx(1,ny)  = mpx(2,ny-1)
377
      mpx(1,ny)  = mpx(2,ny-1)
375
      mpy(1,ny)  = mpy(2,ny-1)
378
      mpy(1,ny)  = mpy(2,ny-1)
376
 
379
 
377
c	  Write some status information
380
c	  Write some status information
378
      print*,'---- MAPPING SCALE FACTORS -----------------------------'
381
      print*,'---- MAPPING SCALE FACTORS -----------------------------'
379
      print*
382
      print*
380
 
383
 
381
      scale_max = mpx(1,1)
384
      scale_max = mpx(1,1)
382
      scale_min = mpx(1,1)
385
      scale_min = mpx(1,1)
383
      do i=1,nx
386
      do i=1,nx
384
      	do j=1,ny
387
      	do j=1,ny
385
      	   if ( mpx(i,j).gt.scale_max ) scale_max = mpx(i,j)
388
      	   if ( mpx(i,j).gt.scale_max ) scale_max = mpx(i,j)
386
      	   if ( mpx(i,j).lt.scale_min ) scale_min = mpx(i,j)
389
      	   if ( mpx(i,j).lt.scale_min ) scale_min = mpx(i,j)
387
      	 enddo
390
      	 enddo
388
      enddo
391
      enddo
389
      print*,'  map scale x (min,max) : ',scale_min,scale_max
392
      print*,'  map scale x (min,max) : ',scale_min,scale_max
390
      scale_max = mpy(1,1)
393
      scale_max = mpy(1,1)
391
      scale_min = mpy(1,1)
394
      scale_min = mpy(1,1)
392
      do i=1,nx
395
      do i=1,nx
393
      	do j=1,ny
396
      	do j=1,ny
394
      	   if ( mpy(i,j).gt.scale_max ) scale_max = mpy(i,j)
397
      	   if ( mpy(i,j).gt.scale_max ) scale_max = mpy(i,j)
395
      	   if ( mpy(i,j).lt.scale_min ) scale_min = mpy(i,j)
398
      	   if ( mpy(i,j).lt.scale_min ) scale_min = mpy(i,j)
396
      	 enddo
399
      	 enddo
397
      enddo
400
      enddo
398
      print*,'  map scale y (min,max) : ',scale_min,scale_max
401
      print*,'  map scale y (min,max) : ',scale_min,scale_max
399
      print*
402
      print*
400
 
403
 
401
c     --------------------------------------------------------------------
404
c     --------------------------------------------------------------------
402
c     Initialize the trajectory calculation
405
c     Initialize the trajectory calculation
403
c     --------------------------------------------------------------------
406
c     --------------------------------------------------------------------
404
 
407
 
405
c     Read start coordinates from file - Format (x,y,lev)
408
c     Read start coordinates from file - Format (x,y,lev)
406
      if (inpmode.eq.-1) then
409
      if (inpmode.eq.-1) then
407
         open(fid,file=strname)
410
         open(fid,file=strname)
408
          do i=1,ntra
411
          do i=1,ntra
409
             read(fid,*) xx0(i),yy0(i),pp0(i)
412
             read(fid,*) xx0(i),yy0(i),pp0(i)
410
          enddo
413
          enddo
411
         close(fid)
414
         close(fid)
412
 
415
 
413
c     Read start coordinates from trajectory file - check consistency of ref time
416
c     Read start coordinates from trajectory file - check consistency of ref time
414
      else
417
      else
415
         call ropen_tra(cdfid,strname,ntra,1,ncol,reftmp,vars,inpmode)
418
         call ropen_tra(cdfid,strname,ntra,1,ncol,reftmp,vars,inpmode)
416
         call read_tra (cdfid,trainp,ntra,1,ncol,inpmode)
419
         call read_tra (cdfid,trainp,ntra,1,ncol,inpmode)
417
         do i=1,ntra
420
         do i=1,ntra
418
            time   = trainp(i,1,1)
421
            time   = trainp(i,1,1)
419
            xx0(i) = trainp(i,1,2) 
422
            xx0(i) = trainp(i,1,2) 
420
            yy0(i) = trainp(i,1,3) 
423
            yy0(i) = trainp(i,1,3) 
421
            pp0(i) = trainp(i,1,4) 
424
            pp0(i) = trainp(i,1,4) 
422
         enddo
425
         enddo
423
         call close_tra(cdfid,inpmode)
426
         call close_tra(cdfid,inpmode)
424
 
427
 
425
         if ( ( vars(2).ne.'x').and.(vars(3).ne.'y') ) then
428
         if ( ( vars(2).ne.'x').and.(vars(3).ne.'y') ) then
426
            print*,' ERROR: starting positions must be in x/y/z format'
429
            print*,' ERROR: starting positions must be in x/y/z format'
427
            stop
430
            stop
428
         endif
431
         endif
429
 
432
 
430
         if ( ( reftime(1).ne.reftmp(1) ).or.
433
         if ( ( reftime(1).ne.reftmp(1) ).or.
431
     >        ( reftime(2).ne.reftmp(2) ).or.
434
     >        ( reftime(2).ne.reftmp(2) ).or.
432
     >        ( reftime(3).ne.reftmp(3) ).or.
435
     >        ( reftime(3).ne.reftmp(3) ).or.
433
     >        ( reftime(4).ne.reftmp(4) ).or.
436
     >        ( reftime(4).ne.reftmp(4) ).or.
434
     >        ( reftime(5).ne.reftmp(5) ) )
437
     >        ( reftime(5).ne.reftmp(5) ) )
435
     >   then
438
     >   then
436
            print*,' WARNING: Inconsistent reference times'
439
            print*,' WARNING: Inconsistent reference times'
437
            write(*,'(5i8)') (reftime(i),i=1,5)
440
            write(*,'(5i8)') (reftime(i),i=1,5)
438
            write(*,'(5i8)') (reftmp (i),i=1,5)
441
            write(*,'(5i8)') (reftmp (i),i=1,5)
439
            print*,'Enter a key to proceed...'
442
            print*,'Enter a key to proceed...'
440
            stop
443
            stop
441
         endif
444
         endif
442
      endif
445
      endif
443
 
446
 
444
c     Transform start coordinates from index space to WRF grid
447
c     Transform start coordinates from index space to WRF grid
445
      do i=1,ntra
448
      do i=1,ntra
446
         xx0(i) = xmin + ( xx0(i) - 1. ) * dx
449
         xx0(i) = xmin + ( xx0(i) - 1. ) * dx
447
         yy0(i) = ymin + ( yy0(i) - 1. ) * dy
450
         yy0(i) = ymin + ( yy0(i) - 1. ) * dy
448
      enddo
451
      enddo
449
 
452
 
450
c     Set sign of time range
453
c     Set sign of time range
451
      reftime(6) = fbflag * reftime(6)
454
      reftime(6) = fbflag * reftime(6)
452
         
455
         
453
c     Write some status information
456
c     Write some status information
454
      print*,'---- REFERENCE DATE---------- ---------------------------'
457
      print*,'---- REFERENCE DATE---------- ---------------------------'
455
      print*
458
      print*
456
      print*,' Reference time (year)  :',reftime(1)
459
      print*,' Reference time (year)  :',reftime(1)
457
      print*,'                (month) :',reftime(2)
460
      print*,'                (month) :',reftime(2)
458
      print*,'                (day)   :',reftime(3)
461
      print*,'                (day)   :',reftime(3)
459
      print*,'                (hour)  :',reftime(4)
462
      print*,'                (hour)  :',reftime(4)
460
      print*,'                (min)   :',reftime(5)
463
      print*,'                (min)   :',reftime(5)
461
      print*,' Time range             :',reftime(6),' min'
464
      print*,' Time range             :',reftime(6),' min'
462
      print*
465
      print*
463
 
466
 
464
C     Save starting positions 
467
C     Save starting positions 
465
      itim = 1
468
      itim = 1
466
      do i=1,ntra
469
      do i=1,ntra
467
         traout(i,itim,1) = 0.
470
         traout(i,itim,1) = 0.
468
         traout(i,itim,2) = xx0(i)
471
         traout(i,itim,2) = xx0(i)
469
         traout(i,itim,3) = yy0(i)
472
         traout(i,itim,3) = yy0(i)
470
         traout(i,itim,4) = pp0(i)
473
         traout(i,itim,4) = pp0(i)
471
      enddo
474
      enddo
472
      
475
      
473
c     Init the flag and the counter for trajectories leaving the domain
476
c     Init the flag and the counter for trajectories leaving the domain
474
      leftcount=0
477
      leftcount=0
475
      do i=1,ntra
478
      do i=1,ntra
476
         leftflag(i)=0
479
         leftflag(i)=0
477
      enddo
480
      enddo
478
 
481
 
479
C     Convert time shifts <tst,ten> from <hh.mm> into fractional time
482
C     Convert time shifts <tst,ten> from <hh.mm> into fractional time
480
      call hhmm2frac(tst,frac)
483
      call hhmm2frac(tst,frac)
481
      tst = frac
484
      tst = frac
482
      call hhmm2frac(ten,frac)
485
      call hhmm2frac(ten,frac)
483
      ten = frac
486
      ten = frac
484
 
487
 
485
c     -----------------------------------------------------------------------
488
c     -----------------------------------------------------------------------
486
c     Loop to calculate trajectories
489
c     Loop to calculate trajectories
487
c     -----------------------------------------------------------------------
490
c     -----------------------------------------------------------------------
488
 
491
 
489
c     Write some status information
492
c     Write some status information
490
      print*,'---- TRAJECTORIES ----------- ---------------------------'
493
      print*,'---- TRAJECTORIES ----------- ---------------------------'
491
      print*    
494
      print*    
492
 
495
 
493
C     Set the time for the first data file (depending on forward/backward mode)
496
C     Set the time for the first data file (depending on forward/backward mode)
494
      if (fbflag.eq.1) then
497
      if (fbflag.eq.1) then
495
        tstart = -tst
498
        tstart = -tst
496
      else
499
      else
497
        tstart = tst
500
        tstart = tst
498
      endif
501
      endif
499
 
502
 
500
c     Set the minute counter for output
503
c     Set the minute counter for output
501
      wstep = 0
504
      wstep = 0
502
 
505
 
503
c     Read wind fields and vertical grid from first file
506
c     Read wind fields and vertical grid from first file
504
      filename = prefix//dat(1)
507
      filename = prefix//dat(1)
505
 
508
 
506
      call frac2hhmm(tstart,tload)
509
      call frac2hhmm(tstart,tload)
507
 
510
 
508
      write(*,'(a16,a20,f7.2)') '  (file,time) : ',
511
      write(*,'(a16,a20,f7.2)') '  (file,time) : ',
509
     >                       trim(filename),tload
512
     >                       trim(filename),tload
510
 
513
 
511
      call input_open (fid,filename)
514
      call input_open (fid,filename)
512
 
515
 
513
      varname='U'					! U
516
      varname='U'					! U
514
 	  call input_var_wrf(fid,varname,uut1,nx,ny,nz)
517
 	  call input_var_wrf(fid,varname,uut1,nx,ny,nz)
515
 
518
 
516
	  varname='V'					! V
519
	  varname='V'					! V
517
	  call input_var_wrf(fid,varname,vvt1,nx,ny,nz)
520
	  call input_var_wrf(fid,varname,vvt1,nx,ny,nz)
518
	
521
	
519
	  varname='W'					! W
522
	  varname='W'					! W
520
	  call input_var_wrf(fid,varname,wwt1,nx,ny,nz)
523
	  call input_var_wrf(fid,varname,wwt1,nx,ny,nz)
521
 
524
 
522
	  varname='geopot'				! geopot.height
525
	  varname='geopot'				! geopot.height
523
	  call input_var_wrf(fid,varname,p3t1,nx,ny,nz)
526
	  call input_var_wrf(fid,varname,p3t1,nx,ny,nz)
524
 
527
 
525
	  varname='geopots'				! geopot.height at surface
528
	  varname='geopots'				! geopot.height at surface
526
	  call input_var_wrf(fid,varname,spt1,nx,ny,1)
529
	  call input_var_wrf(fid,varname,spt1,nx,ny,1)
527
 
530
 
528
      call input_close(fid)
531
      call input_close(fid)
529
 
532
 
530
c     Loop over all input files (time step is <timeinc>)
533
c     Loop over all input files (time step is <timeinc>)
531
      do itm=1,numdat-1
534
      do itm=1,numdat-1
532
 
535
 
533
c       Calculate actual and next time
536
c       Calculate actual and next time
534
        time0 = tstart+real(itm-1)*timeinc*fbflag
537
        time0 = tstart+real(itm-1)*timeinc*fbflag
535
        time1 = time0+timeinc*fbflag
538
        time1 = time0+timeinc*fbflag
536
 
539
 
537
 
540
 
538
c       Copy old velocities and pressure fields to new ones       
541
c       Copy old velocities and pressure fields to new ones       
539
        do i=1,nx*ny*nz
542
        do i=1,nx*ny*nz
540
           uut0(i)=uut1(i)
543
           uut0(i)=uut1(i)
541
           vvt0(i)=vvt1(i)
544
           vvt0(i)=vvt1(i)
542
           wwt0(i)=wwt1(i)
545
           wwt0(i)=wwt1(i)
543
           p3t0(i)=p3t1(i)
546
           p3t0(i)=p3t1(i)
544
        enddo
547
        enddo
545
        do i=1,nx*ny
548
        do i=1,nx*ny
546
           spt0(i)=spt1(i)
549
           spt0(i)=spt1(i)
547
        enddo
550
        enddo
548
 
551
 
549
c       Read wind fields and surface pressure at next time
552
c       Read wind fields and surface pressure at next time
550
        filename = prefix//dat(itm+1)
553
        filename = prefix//dat(itm+1)
551
 
554
 
552
        call frac2hhmm(time1,tload)
555
        call frac2hhmm(time1,tload)
553
        write(*,'(a16,a20,f7.2)') '  (file,time) : ',
556
        write(*,'(a16,a20,f7.2)') '  (file,time) : ',
554
     >                          trim(filename),tload
557
     >                          trim(filename),tload
555
 
558
 
556
        call input_open (fid,filename)
559
        call input_open (fid,filename)
557
 
560
 
558
  	    varname='U'					! U
561
  	    varname='U'					! U
559
	    call input_var_wrf(fid,varname,uut1,nx,ny,nz)
562
	    call input_var_wrf(fid,varname,uut1,nx,ny,nz)
560
 
563
 
561
	    varname='V'					! V
564
	    varname='V'					! V
562
	    call input_var_wrf(fid,varname,vvt1,nx,ny,nz)
565
	    call input_var_wrf(fid,varname,vvt1,nx,ny,nz)
563
	
566
	
564
	    varname='W'					! W
567
	    varname='W'					! W
565
	    call input_var_wrf(fid,varname,wwt1,nx,ny,nz)
568
	    call input_var_wrf(fid,varname,wwt1,nx,ny,nz)
566
 
569
 
567
	    varname='geopot'				! geopot.height
570
	    varname='geopot'				! geopot.height
568
	    call input_var_wrf(fid,varname,p3t1,nx,ny,nz)
571
	    call input_var_wrf(fid,varname,p3t1,nx,ny,nz)
569
 
572
 
570
	    varname='geopots'				! geopot.height
573
	    varname='geopots'				! geopot.height
571
	    call input_var_wrf(fid,varname,spt1,nx,ny,1)
574
	    call input_var_wrf(fid,varname,spt1,nx,ny,1)
572
 
575
 
573
        call input_close(fid)
576
        call input_close(fid)
574
        
577
        
575
C       Determine the first and last loop indices
578
C       Determine the first and last loop indices
576
        if (numdat.eq.2) then
579
        if (numdat.eq.2) then
577
          filo = nint(tst/ts)+1
580
          filo = nint(tst/ts)+1
578
          lalo = nint((timeinc-ten)/ts) 
581
          lalo = nint((timeinc-ten)/ts) 
579
        elseif ( itm.eq.1 ) then
582
        elseif ( itm.eq.1 ) then
580
          filo = nint(tst/ts)+1
583
          filo = nint(tst/ts)+1
581
          lalo = nint(timeinc/ts)
584
          lalo = nint(timeinc/ts)
582
        else if (itm.eq.numdat-1) then
585
        else if (itm.eq.numdat-1) then
583
          filo = 1
586
          filo = 1
584
          lalo = nint((timeinc-ten)/ts)
587
          lalo = nint((timeinc-ten)/ts)
585
        else
588
        else
586
          filo = 1
589
          filo = 1
587
          lalo = nint(timeinc/ts)
590
          lalo = nint(timeinc/ts)
588
        endif
591
        endif
589
 
592
 
590
c       Split the interval <timeinc> into computational time steps <ts>
593
c       Split the interval <timeinc> into computational time steps <ts>
591
        do iloop=filo,lalo
594
        do iloop=filo,lalo
592
 
595
 
593
C         Calculate relative time position in the interval timeinc (0=beginning, 1=end)
596
C         Calculate relative time position in the interval timeinc (0=beginning, 1=end)
594
          reltpos0 = ((real(iloop)-1.)*ts)/timeinc
597
          reltpos0 = ((real(iloop)-1.)*ts)/timeinc
595
          reltpos1 = real(iloop)*ts/timeinc
598
          reltpos1 = real(iloop)*ts/timeinc
596
 
599
 
597
C         Initialize counter for domain leaving trajectories
600
C         Initialize counter for domain leaving trajectories
598
          leftcount=0
601
          leftcount=0
599
 
602
 
600
c         Timestep for all trajectories
603
c         Timestep for all trajectories
601
          do i=1,ntra
604
          do i=1,ntra
602
 
605
 
603
C           Check if trajectory has already left the data domain
606
C           Check if trajectory has already left the data domain
604
            if (leftflag(i).ne.1) then	
607
            if (leftflag(i).ne.1) then	
605
 
608
 
606
c             Iterative Euler timestep (x0,y0,p0 -> x1,y1,p1)
609
c             Iterative Euler timestep (x0,y0,p0 -> x1,y1,p1)
607
              call euler(
610
              call euler(
608
     >               xx1,yy1,pp1,leftflag(i),
611
     >               xx1,yy1,pp1,leftflag(i),
609
     >               xx0(i),yy0(i),pp0(i),reltpos0,reltpos1,
612
     >               xx0(i),yy0(i),pp0(i),reltpos0,reltpos1,
610
     >               ts*3600,numit,jflag,mdv,wfactor,fbflag,
613
     >               ts*3600,numit,jflag,mdv,wfactor,fbflag,
611
     >               spt0,spt1,p3t0,p3t1,uut0,uut1,vvt0,vvt1,wwt0,wwt1,
614
     >               spt0,spt1,p3t0,p3t1,uut0,uut1,vvt0,vvt1,wwt0,wwt1,
612
     >               xmin,ymin,dx,dy,per,hem,nx,ny,nz,mpx,mpy)
615
     >               xmin,ymin,dx,dy,per,hem,nx,ny,nz,mpx,mpy)
613
 
616
 
614
c             Update trajectory position, or increase number of trajectories leaving domain
617
c             Update trajectory position, or increase number of trajectories leaving domain
615
              if (leftflag(i).eq.1) then
618
              if (leftflag(i).eq.1) then
616
                leftcount=leftcount+1
619
                leftcount=leftcount+1
617
                if ( leftcount.lt.10 ) then
620
                if ( leftcount.lt.10 ) then
618
                   print*,'     -> Trajectory ',i,' leaves domain'
621
                   print*,'     -> Trajectory ',i,' leaves domain'
619
                elseif ( leftcount.eq.10 ) then
622
                elseif ( leftcount.eq.10 ) then
620
                   print*,'     -> N>=10 trajectories leave domain'
623
                   print*,'     -> N>=10 trajectories leave domain'
621
                endif
624
                endif
622
              else
625
              else
623
                xx0(i)=xx1      
626
                xx0(i)=xx1      
624
                yy0(i)=yy1
627
                yy0(i)=yy1
625
                pp0(i)=pp1
628
                pp0(i)=pp1
626
              endif
629
              endif
627
 
630
 
628
c          Trajectory has already left data domain (mark as <mdv>)
631
c          Trajectory has already left data domain (mark as <mdv>)
629
           else
632
           else
630
              xx0(i)=mdv
633
              xx0(i)=mdv
631
              yy0(i)=mdv
634
              yy0(i)=mdv
632
              pp0(i)=mdv
635
              pp0(i)=mdv
633
              
636
              
634
           endif
637
           endif
635
 
638
 
636
          enddo
639
          enddo
637
 
640
 
638
C         Save positions only every deltout minutes
641
C         Save positions only every deltout minutes
639
          delta = aint(iloop*60*ts/deltout)-iloop*60*ts/deltout
642
          delta = aint(iloop*60*ts/deltout)-iloop*60*ts/deltout
640
          if (abs(delta).lt.eps) then
643
          if (abs(delta).lt.eps) then
641
            time = time0+reltpos1*timeinc*fbflag
644
            time = time0+reltpos1*timeinc*fbflag
642
            itim = itim + 1
645
            itim = itim + 1
643
            do i=1,ntra
646
            do i=1,ntra
644
               call frac2hhmm(time,tload)
647
               call frac2hhmm(time,tload)
645
               traout(i,itim,1) = tload
648
               traout(i,itim,1) = tload
646
               traout(i,itim,2) = xx0(i)
649
               traout(i,itim,2) = xx0(i)
647
               traout(i,itim,3) = yy0(i) 
650
               traout(i,itim,3) = yy0(i) 
648
               traout(i,itim,4) = pp0(i)
651
               traout(i,itim,4) = pp0(i)
649
            enddo
652
            enddo
650
          endif
653
          endif
651
 
654
 
652
        enddo
655
        enddo
653
 
656
 
654
      enddo
657
      enddo
655
 
658
 
656
c     Transform trajectory position from WRF grid do grid index
659
c     Transform trajectory position from WRF grid do grid index
657
      do i=1,ntra
660
      do i=1,ntra
658
        do j=1,ntim
661
        do j=1,ntim
659
           traout(i,j,2) = ( traout(i,j,2) - xmin ) / dx + 1.
662
           traout(i,j,2) = ( traout(i,j,2) - xmin ) / dx + 1.
660
           traout(i,j,3) = ( traout(i,j,3) - ymin ) / dy + 1.
663
           traout(i,j,3) = ( traout(i,j,3) - ymin ) / dy + 1.
661
        enddo
664
        enddo
662
      enddo
665
      enddo
663
 
666
 
664
c     Write trajectory file
667
c     Write trajectory file
665
      vars(1)  ='time'
668
      vars(1)  ='time'
666
      vars(2)  ='x'
669
      vars(2)  ='x'
667
      vars(3)  ='y'
670
      vars(3)  ='y'
668
      vars(4)  ='z'
671
      vars(4)  ='z'
669
      call wopen_tra(cdfid,cdfname,ntra,ntim,4,reftime,vars,outmode)
672
      call wopen_tra(cdfid,cdfname,ntra,ntim,4,reftime,vars,outmode)
670
      call write_tra(cdfid,traout,ntra,ntim,4,outmode)
673
      call write_tra(cdfid,traout,ntra,ntim,4,outmode)
671
      call close_tra(cdfid,outmode)   
674
      call close_tra(cdfid,outmode)   
672
 
675
 
673
c     Write some status information, and end of program message
676
c     Write some status information, and end of program message
674
      print*  
677
      print*  
675
      print*,'---- STATUS INFORMATION --------------------------------'
678
      print*,'---- STATUS INFORMATION --------------------------------'
676
      print*
679
      print*
677
      print*,'  #leaving domain    ', leftcount
680
      print*,'  #leaving domain    ', leftcount
678
      print*,'  #staying in domain ', ntra-leftcount
681
      print*,'  #staying in domain ', ntra-leftcount
679
      print*
682
      print*
680
      print*,'              *** END OF PROGRAM CALTRA ***'
683
      print*,'              *** END OF PROGRAM CALTRA ***'
681
      print*,'========================================================='
684
      print*,'========================================================='
682
 
685
 
683
      stop
686
      stop
684
 
687
 
685
c     ------------------------------------------------------------------
688
c     ------------------------------------------------------------------
686
c     Exception handling
689
c     Exception handling
687
c     ------------------------------------------------------------------
690
c     ------------------------------------------------------------------
688
 
691
 
689
 991  write(*,*) '*** ERROR: all start points outside the data domain'
692
 991  write(*,*) '*** ERROR: all start points outside the data domain'
690
      call exit(1)
693
      call exit(1)
691
      
694
      
692
 992  write(*,*) '*** ERROR: close arrays on files (prog. closear)'
695
 992  write(*,*) '*** ERROR: close arrays on files (prog. closear)'
693
      call exit(1)
696
      call exit(1)
694
 
697
 
695
 993  write(*,*) '*** ERROR: problems with array size'
698
 993  write(*,*) '*** ERROR: problems with array size'
696
      call exit(1)
699
      call exit(1)
697
 
700
 
698
      end 
701
      end 
699
 
702
 
700
 
703
 
701
c     *******************************************************************
704
c     *******************************************************************
702
c     * Time step : either Euler or Runge-Kutta                         *
705
c     * Time step : either Euler or Runge-Kutta                         *
703
c     *******************************************************************
706
c     *******************************************************************
704
 
707
 
705
C     Time-step from (x0,y0,p0) to (x1,y1,p1)
708
C     Time-step from (x0,y0,p0) to (x1,y1,p1)
706
C
709
C
707
C     (x0,y0,p0) input	coordinates (long,lat,p) for starting point
710
C     (x0,y0,p0) input	coordinates (long,lat,p) for starting point
708
C     (x1,y1,p1) output	coordinates (long,lat,p) for end point
711
C     (x1,y1,p1) output	coordinates (long,lat,p) for end point
709
C     deltat	 input	timestep in seconds
712
C     deltat	 input	timestep in seconds
710
C     numit	 input	number of iterations
713
C     numit	 input	number of iterations
711
C     jump	 input  flag (=1 trajectories don't enter the ground)
714
C     jump	 input  flag (=1 trajectories don't enter the ground)
712
C     left	 output	flag (=1 if trajectory leaves data domain)
715
C     left	 output	flag (=1 if trajectory leaves data domain)
713
 
716
 
714
c     -------------------------------------------------------------------
717
c     -------------------------------------------------------------------
715
c     Iterative Euler time step
718
c     Iterative Euler time step
716
c     -------------------------------------------------------------------
719
c     -------------------------------------------------------------------
717
 
720
 
718
      subroutine euler(x1,y1,p1,left,x0,y0,p0,reltpos0,reltpos1,
721
      subroutine euler(x1,y1,p1,left,x0,y0,p0,reltpos0,reltpos1,
719
     >                 deltat,numit,jump,mdv,wfactor,fbflag,
722
     >                 deltat,numit,jump,mdv,wfactor,fbflag,
720
     >		           spt0,spt1,p3d0,p3d1,uut0,uut1,vvt0,vvt1,wwt0,wwt1,
723
     >		           spt0,spt1,p3d0,p3d1,uut0,uut1,vvt0,vvt1,wwt0,wwt1,
721
     >                 xmin,ymin,dx,dy,per,hem,nx,ny,nz,mpx,mpy)
724
     >                 xmin,ymin,dx,dy,per,hem,nx,ny,nz,mpx,mpy)
722
 
725
 
723
      implicit none
726
      implicit none
724
 
727
 
725
c     Flag for test mode
728
c     Flag for test mode
726
      integer      test
729
      integer      test
727
      parameter    (test=0)
730
      parameter    (test=0)
728
 
731
 
729
c     Declaration of subroutine parameters
732
c     Declaration of subroutine parameters
730
      integer      nx,ny,nz
733
      integer      nx,ny,nz
731
      real         x1,y1,p1
734
      real         x1,y1,p1
732
      integer      left
735
      integer      left
733
      real	       x0,y0,p0
736
      real	       x0,y0,p0
734
      real         reltpos0,reltpos1
737
      real         reltpos0,reltpos1
735
      real   	   deltat
738
      real   	   deltat
736
      integer      numit
739
      integer      numit
737
      integer      jump
740
      integer      jump
738
      real         wfactor
741
      real         wfactor
739
      integer      fbflag
742
      integer      fbflag
740
      real     	   spt0(nx*ny)   ,spt1(nx*ny)
743
      real     	   spt0(nx*ny)   ,spt1(nx*ny)
741
      real         uut0(nx*ny*nz),uut1(nx*ny*nz)
744
      real         uut0(nx*ny*nz),uut1(nx*ny*nz)
742
      real 	       vvt0(nx*ny*nz),vvt1(nx*ny*nz)
745
      real 	       vvt0(nx*ny*nz),vvt1(nx*ny*nz)
743
      real         wwt0(nx*ny*nz),wwt1(nx*ny*nz)
746
      real         wwt0(nx*ny*nz),wwt1(nx*ny*nz)
744
      real         p3d0(nx*ny*nz),p3d1(nx*ny*nz)
747
      real         p3d0(nx*ny*nz),p3d1(nx*ny*nz)
745
      real         xmin,ymin,dx,dy
748
      real         xmin,ymin,dx,dy
746
      integer      per
749
      integer      per
747
      integer      hem
750
      integer      hem
748
      real         mdv
751
      real         mdv
749
      real         mpx(nx*ny),mpy(nx*ny)
752
      real         mpx(nx*ny),mpy(nx*ny)
750
 
753
 
751
c     Auxiliary variables
754
c     Auxiliary variables
752
      real         xmax,ymax
755
      real         xmax,ymax
753
      real	       xind,yind,pind
756
      real	       xind,yind,pind
754
 
757
 
755
      real	       u0,v0,w0,u1,v1,w1,u,v,w,sp
758
      real	       u0,v0,w0,u1,v1,w1,u,v,w,sp
756
      integer	   icount
759
      integer	   icount
757
      character    ch
760
      character    ch
758
      real         mpsc_x,mpsc_y
761
      real         mpsc_x,mpsc_y
759
 
762
 
760
c     Externals    
763
c     Externals    
761
      real         int_index4
764
      real         int_index4
762
      external     int_index4
765
      external     int_index4
763
 
766
 
764
c     Reset the flag for domain-leaving
767
c     Reset the flag for domain-leaving
765
      left=0
768
      left=0
766
 
769
 
767
c     Set the esat-north boundary of the domain
770
c     Set the esat-north boundary of the domain
768
      xmax = xmin+real(nx-1)*dx
771
      xmax = xmin+real(nx-1)*dx
769
      ymax = ymin+real(ny-1)*dy
772
      ymax = ymin+real(ny-1)*dy
770
 
773
 
771
C     Interpolate wind fields to starting position (x0,y0,p0)
774
C     Interpolate wind fields to starting position (x0,y0,p0)
772
      call get_index4 (xind,yind,pind,x0,y0,p0,reltpos0,
775
      call get_index4 (xind,yind,pind,x0,y0,p0,reltpos0,
773
     >                 p3d0,p3d1,spt0,spt1,3,
776
     >                 p3d0,p3d1,spt0,spt1,3,
774
     >                 nx,ny,nz,xmin,ymin,dx,dy,mdv)
777
     >                 nx,ny,nz,xmin,ymin,dx,dy,mdv)
775
 
778
 
776
      u0 = int_index4(uut0,uut1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
779
      u0 = int_index4(uut0,uut1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
777
      v0 = int_index4(vvt0,vvt1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
780
      v0 = int_index4(vvt0,vvt1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
778
      w0 = int_index4(wwt0,wwt1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
781
      w0 = int_index4(wwt0,wwt1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
779
 
782
 
780
c     Force the near-surface wind to zero
783
c     Force the near-surface wind to zero
781
      if (pind.lt.1.) w0=w0*pind
784
      if (pind.lt.1.) w0=w0*pind
782
 
785
 
783
C     For first iteration take ending position equal to starting position
786
C     For first iteration take ending position equal to starting position
784
      x1=x0
787
      x1=x0
785
      y1=y0
788
      y1=y0
786
      p1=p0
789
      p1=p0
787
 
790
 
788
C     Iterative calculation of new position
791
C     Iterative calculation of new position
789
      do icount=1,numit
792
      do icount=1,numit
790
 
793
 
791
C        Calculate new winds for advection
794
C        Calculate new winds for advection
792
         call get_index4 (xind,yind,pind,x1,y1,p1,reltpos1,
795
         call get_index4 (xind,yind,pind,x1,y1,p1,reltpos1,
793
     >                    p3d0,p3d1,spt0,spt1,3,
796
     >                    p3d0,p3d1,spt0,spt1,3,
794
     >                    nx,ny,nz,xmin,ymin,dx,dy,mdv)
797
     >                    nx,ny,nz,xmin,ymin,dx,dy,mdv)
795
         u1 = int_index4(uut0,uut1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
798
         u1 = int_index4(uut0,uut1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
796
         v1 = int_index4(vvt0,vvt1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
799
         v1 = int_index4(vvt0,vvt1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
797
         w1 = int_index4(wwt0,wwt1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
800
         w1 = int_index4(wwt0,wwt1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
798
 
801
 
799
c        Force the near-surface wind to zero
802
c        Force the near-surface wind to zero
800
         if (pind.lt.1.) w1=w1*pind
803
         if (pind.lt.1.) w1=w1*pind
801
 
804
 
802
c        Get the new velocity in between
805
c        Get the new velocity in between
803
         u=(u0+u1)/2.
806
         u=(u0+u1)/2.
804
         v=(v0+v1)/2.
807
         v=(v0+v1)/2.
805
         w=(w0+w1)/2.
808
         w=(w0+w1)/2.
806
         
809
         
807
c        Get the mapping scale factors for this position
810
c        Get the mapping scale factors for this position
808
         mpsc_x = int_index4 (mpx,mpx,nx,ny,1,xind,yind,1.,reltpos1,mdv)
811
         mpsc_x = int_index4 (mpx,mpx,nx,ny,1,xind,yind,1.,reltpos1,mdv)
809
         mpsc_y = int_index4 (mpy,mpy,nx,ny,1,xind,yind,1.,reltpos1,mdv)
812
         mpsc_y = int_index4 (mpy,mpy,nx,ny,1,xind,yind,1.,reltpos1,mdv)
810
 
813
 
811
C        Calculate new positions (adapted for cartesian grid)
814
C        Calculate new positions (adapted for cartesian grid)
812
         x1 = x0 + fbflag * deltat * u * dx/mpsc_x
815
         x1 = x0 + fbflag * deltat * u * dx/mpsc_x
813
         y1 = y0 + fbflag * deltat * v * dy/mpsc_y
816
         y1 = y0 + fbflag * deltat * v * dy/mpsc_y
814
         p1 = p0 + fbflag * deltat * w * wfactor
817
         p1 = p0 + fbflag * deltat * w * wfactor
815
 
818
 
816
c        Write the update positions for tests
819
c        Write the update positions for tests
817
         if ( (icount.eq.3).and.(test.eq.1) ) then
820
         if ( (icount.eq.3).and.(test.eq.1) ) then
818
            write(*,'(10f10.2)') x0,u,x1-x0,p0,w,p1-p0
821
            write(*,'(10f10.2)') x0,u,x1-x0,p0,w,p1-p0
819
         endif
822
         endif
820
 
823
 
821
C       Interpolate surface geop.height to actual position
824
C       Interpolate surface geop.height to actual position
822
        call get_index4 (xind,yind,pind,x1,y1,0.,reltpos1,
825
        call get_index4 (xind,yind,pind,x1,y1,0.,reltpos1,
823
     >                   p3d0,p3d1,spt0,spt1,3,
826
     >                   p3d0,p3d1,spt0,spt1,3,
824
     >                   nx,ny,nz,xmin,ymin,dx,dy,mdv)
827
     >                   nx,ny,nz,xmin,ymin,dx,dy,mdv)
825
        sp = int_index4 (spt0,spt1,nx,ny,1,xind,yind,1.,reltpos1,mdv)
828
        sp = int_index4 (spt0,spt1,nx,ny,1,xind,yind,1.,reltpos1,mdv)
826
 
829
 
827
c       Handle trajectories which cross the lower boundary (jump flag)
830
c       Handle trajectories which cross the lower boundary (jump flag)
828
        if ((jump.eq.1).and.(p1.lt.sp)) p1=sp+10.
831
        if ((jump.eq.1).and.(p1.lt.sp)) p1=sp+10.
829
 
832
 
830
c       Handle periodioc boundaries for 'ideal' mode
833
c       Handle periodioc boundaries for 'ideal' mode
831
        if ( per.eq.1 ) then
834
        if ( per.eq.1 ) then
832
	        if ( x1.gt.xmax ) x1=x1-xmax
835
	        if ( x1.gt.xmax ) x1=x1-xmax
833
	        if ( x1.lt.xmin ) x1=x1+xmax
836
	        if ( x1.lt.xmin ) x1=x1+xmax
834
        endif
837
        endif
835
 
838
 
836
C       Check if trajectory leaves data domain
839
C       Check if trajectory leaves data domain
837
        if ( (x1.lt.xmin).or.(x1.gt.xmax).or.
840
        if ( (x1.lt.xmin).or.(x1.gt.xmax).or.
838
     >       (y1.lt.ymin).or.(y1.gt.ymax).or.(p1.lt.sp) )   ! geopot : .lt. ; pressure: .gt.
841
     >       (y1.lt.ymin).or.(y1.gt.ymax).or.(p1.lt.sp) )   ! geopot : .lt. ; pressure: .gt.
839
     >  then
842
     >  then
840
          left=1
843
          left=1
841
          print*,x1,y1,p1
844
          print*,x1,y1,p1
842
          print*,xmin,ymin
845
          print*,xmin,ymin
843
	      print*,xmax,ymax
846
	      print*,xmax,ymax
844
	      print*,sp
847
	      print*,sp
845
          goto 100
848
          goto 100
846
        endif
849
        endif
847
 
850
 
848
      enddo
851
      enddo
849
 
852
 
850
c     Exit point for subroutine
853
c     Exit point for subroutine
851
 100  continue
854
 100  continue
852
 
855
 
853
      return
856
      return
854
 
857
 
855
      end
858
      end
856
 
859
 
857
 
860
 
858
c     ----------------------------------------------------------------------
861
c     ----------------------------------------------------------------------
859
c     Get spherical distance between lat/lon points
862
c     Get spherical distance between lat/lon points
860
c     ----------------------------------------------------------------------
863
c     ----------------------------------------------------------------------
861
 
864
 
862
      real function sdis(xp,yp,xq,yq,unit)
865
      real function sdis(xp,yp,xq,yq,unit)
863
 
866
 
864
c     Calculates spherical distance (in km) between two points given
867
c     Calculates spherical distance (in km) between two points given
865
c     by their spherical coordinates (xp,yp) and (xq,yq), respectively.
868
c     by their spherical coordinates (xp,yp) and (xq,yq), respectively.
866
 
869
 
867
      real,parameter ::       re=6371.229
870
      real,parameter ::       re=6371.229
868
      real,parameter ::       rad2deg=180./3.14159265
871
      real,parameter ::       rad2deg=180./3.14159265
869
      real,parameter ::       deg2rad=3.141592654/180.
872
      real,parameter ::       deg2rad=3.141592654/180.
870
 
873
 
871
      real         xp,yp,xq,yq,arg
874
      real         xp,yp,xq,yq,arg
872
      character*80 unit
875
      character*80 unit
873
      real         dlon,dlat,alpha,rlat1,ralt2
876
      real         dlon,dlat,alpha,rlat1,ralt2
874
 
877
 
875
      if ( unit.eq.'km' ) then
878
      if ( unit.eq.'km' ) then
876
 
879
 
877
         arg=sin(yp*deg2rad)*sin(yq*deg2rad)+
880
         arg=sin(yp*deg2rad)*sin(yq*deg2rad)+
878
     >              cos(yp*deg2rad)*cos(yq*deg2rad)*cos((xp-xq)*deg2rad)
881
     >              cos(yp*deg2rad)*cos(yq*deg2rad)*cos((xp-xq)*deg2rad)
879
         if (arg.lt.-1.) arg=-1.
882
         if (arg.lt.-1.) arg=-1.
880
         if (arg.gt.1.) arg=1.
883
         if (arg.gt.1.) arg=1.
881
         sdis=re*acos(arg)
884
         sdis=re*acos(arg)
882
 
885
 
883
      elseif ( unit.eq.'deg' ) then
886
      elseif ( unit.eq.'deg' ) then
884
 
887
 
885
         dlon = xp-xq
888
         dlon = xp-xq
886
         if ( dlon.gt. 180. ) dlon = dlon - 360.
889
         if ( dlon.gt. 180. ) dlon = dlon - 360.
887
         if ( dlon.lt.-180. ) dlon = dlon + 360.
890
         if ( dlon.lt.-180. ) dlon = dlon + 360.
888
         sdis = sqrt( dlon**2 + (yp-yq)**2 )
891
         sdis = sqrt( dlon**2 + (yp-yq)**2 )
889
 
892
 
890
      elseif ( unit.eq.'km.haversine' ) then
893
      elseif ( unit.eq.'km.haversine' ) then
891
 
894
 
892
        dlat   =  (yp - yq)*deg2rad
895
        dlat   =  (yp - yq)*deg2rad
893
        dlon   =  (xp - xq)*deg2rad
896
        dlon   =  (xp - xq)*deg2rad
894
        rlat1   =  yp*deg2rad
897
        rlat1   =  yp*deg2rad
895
        rlat2   =  yq*deg2rad
898
        rlat2   =  yq*deg2rad
896
 
899
 
897
        alpha  = ( sin(0.5*dlat)**2 ) +
900
        alpha  = ( sin(0.5*dlat)**2 ) +
898
     >           ( sin(0.5*dlon)**2 )*cos(rlat1)*cos(rlat2)
901
     >           ( sin(0.5*dlon)**2 )*cos(rlat1)*cos(rlat2)
899
 
902
 
900
        sdis = 4. * re * atan2( sqrt(alpha),1. + sqrt( 1. - alpha ) )
903
        sdis = 4. * re * atan2( sqrt(alpha),1. + sqrt( 1. - alpha ) )
901
 
904
 
902
      endif
905
      endif
903
 
906
 
904
      end
907
      end