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