2 |
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 Number of iterations for iterative Euler scheme
|
|
|
35 |
integer numit
|
|
|
36 |
parameter (numit=3)
|
|
|
37 |
|
|
|
38 |
c Input and output format for trajectories (see iotra.f)
|
|
|
39 |
integer inpmode
|
|
|
40 |
integer outmode
|
|
|
41 |
|
17 |
michaesp |
42 |
c Missing data value
|
|
|
43 |
real mdv
|
|
|
44 |
parameter (mdv = -999. )
|
|
|
45 |
|
2 |
michaesp |
46 |
c Filename prefix (typically 'P')
|
|
|
47 |
character*1 prefix
|
|
|
48 |
parameter (prefix='P')
|
|
|
49 |
|
|
|
50 |
c Externals
|
|
|
51 |
real sdis
|
|
|
52 |
external sdis
|
|
|
53 |
|
|
|
54 |
c --------------------------------------------------------------------
|
|
|
55 |
c Declaration of variables
|
|
|
56 |
c --------------------------------------------------------------------
|
|
|
57 |
|
|
|
58 |
c Input parameters
|
|
|
59 |
integer fbflag ! Flag for forward/backward mode
|
|
|
60 |
integer numdat ! Number of input files
|
23 |
michaesp |
61 |
character*13 dat(ndatmax) ! Dates of input files
|
2 |
michaesp |
62 |
real timeinc ! Time increment between input files
|
|
|
63 |
integer per ! Periodicity flag (=0 if none)
|
|
|
64 |
integer ntra ! Number of trajectories
|
|
|
65 |
character*80 cdfname ! Name of output files
|
|
|
66 |
real ts ! Time step
|
|
|
67 |
real tst,ten ! Shift of start and end time relative to first data file
|
|
|
68 |
integer deltout ! Output time interval (in minutes)
|
|
|
69 |
integer jflag ! Jump flag (if =1 ground-touching trajectories reenter atmosphere)
|
|
|
70 |
real wfactor ! Factor for vertical velocity field
|
|
|
71 |
character*80 strname ! File with start positions
|
|
|
72 |
character*80 timecheck ! Either 'yes' or 'no'
|
|
|
73 |
character*10 gridflag ! "ideal" or "real"
|
21 |
michaesp |
74 |
real aglconst ! constant level above ground
|
|
|
75 |
real isen ! isentropic level
|
|
|
76 |
character*80 vertmode ! Vertical mode (kinematic,aglconst)
|
2 |
michaesp |
77 |
|
|
|
78 |
c Trajectories
|
|
|
79 |
integer ncol ! Number of columns for insput trajectories
|
|
|
80 |
real,allocatable, dimension (:,:,:) :: trainp ! Input start coordinates (ntra,1,ncol)
|
|
|
81 |
real,allocatable, dimension (:,:,:) :: traout ! Output trajectories (ntra,ntim,4)
|
|
|
82 |
integer reftime(6) ! Reference date
|
|
|
83 |
character*80 vars(200) ! Field names
|
|
|
84 |
real,allocatable, dimension (:) :: xx0,yy0,pp0 ! Position of air parcels
|
|
|
85 |
integer,allocatable, dimension (:) :: leftflag ! Flag for domain-leaving
|
|
|
86 |
real xx1,yy1,pp1 ! Updated position of air parcel
|
|
|
87 |
integer leftcount ! Number of domain leaving trajectories
|
|
|
88 |
integer ntim ! Number of output time steps
|
|
|
89 |
|
|
|
90 |
c Meteorological fields
|
|
|
91 |
real,allocatable, dimension (:) :: spt0,spt1 ! Surface pressure
|
|
|
92 |
real,allocatable, dimension (:) :: uut0,uut1 ! Zonal wind
|
|
|
93 |
real,allocatable, dimension (:) :: vvt0,vvt1 ! Meridional wind
|
|
|
94 |
real,allocatable, dimension (:) :: wwt0,wwt1 ! Vertical wind
|
|
|
95 |
real,allocatable, dimension (:) :: p3t0,p3t1 ! 3d-pressure
|
21 |
michaesp |
96 |
real,allocatable, dimension (:) :: tht0,tht1 ! potential temperature
|
|
|
97 |
real,allocatable, dimension (:) :: sth0,sth1 ! surface potential temperature
|
2 |
michaesp |
98 |
|
|
|
99 |
c Grid description
|
|
|
100 |
real pollon,pollat ! Longitude/latitude of pole
|
|
|
101 |
real ak(nlevmax) ! Vertical layers and levels
|
|
|
102 |
real bk(nlevmax)
|
|
|
103 |
real xmin,xmax ! Zonal grid extension
|
|
|
104 |
real ymin,ymax ! Meridional grid extension
|
|
|
105 |
integer nx,ny,nz ! Grid dimensions
|
|
|
106 |
real dx,dy ! Horizontal grid resolution
|
|
|
107 |
integer hem ! Flag for hemispheric domain
|
|
|
108 |
real,allocatable, dimension (:,:) :: mpx,mpy ! Map scale factor in x,y direction
|
|
|
109 |
real,allocatable, dimension (:,:) :: lat,lon ! Map scale factor in x,y direction
|
|
|
110 |
|
|
|
111 |
c Auxiliary variables
|
|
|
112 |
real delta,rd
|
|
|
113 |
integer itm,iloop,i,j,k,filo,lalo
|
|
|
114 |
integer ierr,stat
|
|
|
115 |
integer cdfid,fid
|
|
|
116 |
real tstart,time0,time1,time
|
|
|
117 |
real reltpos0,reltpos1
|
|
|
118 |
real xind,yind,pind,pp,sp,stagz
|
|
|
119 |
character*80 filename,varname
|
|
|
120 |
integer reftmp(6)
|
|
|
121 |
character ch
|
|
|
122 |
real frac,tload
|
|
|
123 |
integer itim
|
|
|
124 |
integer wstep
|
|
|
125 |
character*80 radunit
|
|
|
126 |
real lon1,lat1,lon2,lat2
|
|
|
127 |
real scale_max,scale_min
|
|
|
128 |
|
21 |
michaesp |
129 |
c Externals
|
|
|
130 |
real int_index4
|
|
|
131 |
external int_index4
|
|
|
132 |
|
2 |
michaesp |
133 |
c --------------------------------------------------------------------
|
|
|
134 |
c Start of program, Read parameters
|
|
|
135 |
c --------------------------------------------------------------------
|
|
|
136 |
|
|
|
137 |
c Write start message
|
|
|
138 |
print*,'========================================================='
|
|
|
139 |
print*,' *** START OF PROGRAM CALTRA ***'
|
|
|
140 |
print*
|
|
|
141 |
|
|
|
142 |
c Open the parameter file
|
|
|
143 |
open(9,file='caltra.param')
|
|
|
144 |
|
|
|
145 |
c Read flag for forward/backward mode (fbflag)
|
|
|
146 |
read(9,*) fbflag
|
|
|
147 |
|
|
|
148 |
c Read number of input files (numdat)
|
|
|
149 |
read(9,*) numdat
|
|
|
150 |
if (numdat.gt.ndatmax) then
|
|
|
151 |
print*,' ERROR: too many input files ',numdat,ndatmax
|
|
|
152 |
goto 993
|
|
|
153 |
endif
|
|
|
154 |
|
|
|
155 |
c Read list of input dates (dat, sort depending on forward/backward mode)
|
|
|
156 |
if (fbflag.eq.1) then
|
|
|
157 |
do itm=1,numdat
|
23 |
michaesp |
158 |
read(9,'(a)') dat(itm)
|
2 |
michaesp |
159 |
enddo
|
|
|
160 |
else
|
|
|
161 |
do itm=numdat,1,-1
|
23 |
michaesp |
162 |
read(9,'(a)') dat(itm)
|
2 |
michaesp |
163 |
enddo
|
|
|
164 |
endif
|
|
|
165 |
|
|
|
166 |
c Read time increment between input files (timeinc)
|
|
|
167 |
read(9,*) timeinc
|
|
|
168 |
|
|
|
169 |
C Read if data domain is periodic and its periodicity
|
|
|
170 |
read(9,*) per
|
|
|
171 |
|
|
|
172 |
|
|
|
173 |
c Read the number of trajectories and name of position file
|
|
|
174 |
read(9,*) strname
|
|
|
175 |
read(9,*) ntra
|
|
|
176 |
read(9,*) ncol
|
|
|
177 |
if (ntra.eq.0) goto 991
|
|
|
178 |
|
|
|
179 |
C Read the name of the output trajectory file and set the constants file
|
|
|
180 |
read(9,*) cdfname
|
|
|
181 |
|
|
|
182 |
C Read the timestep for trajectory calculation (convert from minutes to hours)
|
|
|
183 |
read(9,*) ts
|
|
|
184 |
ts=ts/60.
|
|
|
185 |
|
|
|
186 |
C Read shift of start and end time relative to first data file
|
|
|
187 |
read(9,*) tst
|
|
|
188 |
read(9,*) ten
|
|
|
189 |
|
|
|
190 |
C Read output time interval (in minutes)
|
|
|
191 |
read(9,*) deltout
|
|
|
192 |
|
|
|
193 |
C Read jumpflag (if =1 ground-touching trajectories reenter the atmosphere)
|
|
|
194 |
read(9,*) jflag
|
|
|
195 |
|
|
|
196 |
C Read factor for vertical velocity field
|
|
|
197 |
read(9,*) wfactor
|
|
|
198 |
|
|
|
199 |
c Read the reference time and the time range
|
|
|
200 |
read(9,*) reftime(1) ! year
|
|
|
201 |
read(9,*) reftime(2) ! month
|
|
|
202 |
read(9,*) reftime(3) ! day
|
|
|
203 |
read(9,*) reftime(4) ! hour
|
|
|
204 |
read(9,*) reftime(5) ! min
|
|
|
205 |
read(9,*) reftime(6) ! time range (in min)
|
|
|
206 |
|
|
|
207 |
c Read flag for 'no time check'
|
|
|
208 |
read(9,*) timecheck
|
21 |
michaesp |
209 |
|
|
|
210 |
c Read flag for trajectories at constant level above ground
|
|
|
211 |
read(9,*) aglconst
|
|
|
212 |
|
|
|
213 |
c Read flag for trajectories at constant isentropic level
|
|
|
214 |
read(9,*) isen
|
2 |
michaesp |
215 |
|
|
|
216 |
c Close the input file
|
|
|
217 |
close(9)
|
|
|
218 |
|
21 |
michaesp |
219 |
c Set the vertical mode of the trajectory calculation
|
|
|
220 |
vertmode = 'kinematic'
|
|
|
221 |
if ( aglconst.gt.0. ) then
|
|
|
222 |
vertmode = 'aglconst'
|
|
|
223 |
endif
|
|
|
224 |
if ( isen.gt.0. ) then
|
|
|
225 |
vertmode = 'isen'
|
|
|
226 |
endif
|
|
|
227 |
|
2 |
michaesp |
228 |
c Calculate the number of output time steps
|
|
|
229 |
ntim = abs(reftime(6)/deltout) + 1
|
|
|
230 |
|
|
|
231 |
c Set the formats of the input and output files
|
|
|
232 |
call mode_tra(inpmode,strname)
|
|
|
233 |
call mode_tra(outmode,cdfname)
|
|
|
234 |
if (outmode.eq.-1) outmode=1
|
|
|
235 |
|
|
|
236 |
c Write some status information
|
|
|
237 |
print*,'---- INPUT PARAMETERS -----------------------------------'
|
|
|
238 |
print*
|
|
|
239 |
print*,' Forward/Backward : ',fbflag
|
|
|
240 |
print*,' #input files : ',numdat
|
|
|
241 |
print*,' First/last input file : ',trim(dat(1)),' ... ',
|
|
|
242 |
> trim(dat(numdat))
|
|
|
243 |
print*,' time increment : ',timeinc
|
|
|
244 |
print*,' Output file : ',trim(cdfname)
|
|
|
245 |
print*,' Time step (min) : ',60.*ts
|
|
|
246 |
write(*,'(a27,f7.2,f7.2)') ' Time shift (start,end) : ',tst,ten
|
|
|
247 |
print*,' Output time interval : ',deltout
|
|
|
248 |
print*,' Jump flag : ',jflag
|
|
|
249 |
print*,' Vertical wind (scale) : ',wfactor
|
|
|
250 |
print*,' Trajectory pos file : ',trim(strname)
|
|
|
251 |
print*,' # of trajectories : ',ntra
|
|
|
252 |
print*,' # of output timesteps : ',ntim
|
|
|
253 |
if ( inpmode.eq.-1) then
|
|
|
254 |
print*,' Input format : (x,y,z)-list'
|
|
|
255 |
else
|
|
|
256 |
print*,' Input format : ',inpmode
|
|
|
257 |
endif
|
|
|
258 |
print*,' Output format : ',outmode
|
|
|
259 |
print*,' Periodicity : ',per
|
|
|
260 |
print*,' Time check : ',trim(timecheck)
|
21 |
michaesp |
261 |
if ( vertmode.eq.'aglconst' ) then
|
|
|
262 |
print*,' Z = const AGL : ',aglconst
|
|
|
263 |
endif
|
|
|
264 |
if ( vertmode.eq.'isen' ) then
|
|
|
265 |
print*,' isentropic : ',isen
|
|
|
266 |
endif
|
|
|
267 |
|
2 |
michaesp |
268 |
print*
|
|
|
269 |
|
|
|
270 |
print*,'---- FIXED NUMERICAL PARAMETERS -------------------------'
|
|
|
271 |
print*
|
|
|
272 |
print*,' Number of iterations : ',numit
|
|
|
273 |
print*,' Filename prefix : ',prefix
|
|
|
274 |
print*,' Missing data value : ',mdv
|
|
|
275 |
print*
|
|
|
276 |
|
|
|
277 |
c --------------------------------------------------------------------
|
|
|
278 |
c Read grid parameters, checks and allocate memory
|
|
|
279 |
c --------------------------------------------------------------------
|
|
|
280 |
|
|
|
281 |
c Read the constant grid parameters (nx,ny,nz,xmin,xmax,ymin,ymax,pollon,pollat)
|
|
|
282 |
c The negative <-fid> of the file identifier is used as a flag for parameter retrieval
|
|
|
283 |
filename = prefix//dat(1)
|
|
|
284 |
varname = 'U'
|
|
|
285 |
nx = 1
|
|
|
286 |
ny = 1
|
|
|
287 |
nz = 1
|
|
|
288 |
tload = -tst
|
|
|
289 |
call input_open (fid,filename)
|
|
|
290 |
call input_grid_wrf(fid,xmin,xmax,ymin,ymax,dx,dy,nx,ny,nz)
|
|
|
291 |
call input_close(fid)
|
|
|
292 |
|
|
|
293 |
C Check if the number of levels is too large
|
|
|
294 |
if (nz.gt.nlevmax) goto 993
|
|
|
295 |
|
|
|
296 |
C Set logical flag for hemispheric mode (not yet supported)
|
|
|
297 |
hem = 0
|
|
|
298 |
|
|
|
299 |
C Allocate memory for some meteorological arrays
|
|
|
300 |
allocate(spt0(nx*ny),stat=stat)
|
|
|
301 |
if (stat.ne.0) print*,'*** error allocating array spt0 ***' ! Surface geopotential height
|
|
|
302 |
allocate(spt1(nx*ny),stat=stat)
|
|
|
303 |
if (stat.ne.0) print*,'*** error allocating array spt1 ***'
|
|
|
304 |
allocate(uut0(nx*ny*nz),stat=stat)
|
|
|
305 |
if (stat.ne.0) print*,'*** error allocating array uut0 ***' ! Zonal wind
|
|
|
306 |
allocate(uut1(nx*ny*nz),stat=stat)
|
|
|
307 |
if (stat.ne.0) print*,'*** error allocating array uut1 ***'
|
|
|
308 |
allocate(vvt0(nx*ny*nz),stat=stat)
|
|
|
309 |
if (stat.ne.0) print*,'*** error allocating array vvt0 ***' ! Meridional wind
|
|
|
310 |
allocate(vvt1(nx*ny*nz),stat=stat)
|
|
|
311 |
if (stat.ne.0) print*,'*** error allocating array vvt1 ***'
|
|
|
312 |
allocate(wwt0(nx*ny*nz),stat=stat)
|
|
|
313 |
if (stat.ne.0) print*,'*** error allocating array wwt0 ***' ! Vertical wind
|
|
|
314 |
allocate(wwt1(nx*ny*nz),stat=stat)
|
|
|
315 |
if (stat.ne.0) print*,'*** error allocating array wwt1 ***'
|
|
|
316 |
allocate(p3t0(nx*ny*nz),stat=stat)
|
|
|
317 |
if (stat.ne.0) print*,'*** error allocating array p3t0 ***' ! geopotential height
|
|
|
318 |
allocate(p3t1(nx*ny*nz),stat=stat)
|
|
|
319 |
if (stat.ne.0) print*,'*** error allocating array p3t1 ***'
|
|
|
320 |
allocate(mpx(nx,ny),stat=stat)
|
|
|
321 |
if (stat.ne.0) print*,'*** error allocating array mpx * **' ! Map scale factor in x
|
|
|
322 |
allocate(mpy(nx,ny),stat=stat)
|
|
|
323 |
if (stat.ne.0) print*,'*** error allocating array mpy ***' ! Map scale factor in y
|
|
|
324 |
allocate(lon(nx,ny),stat=stat)
|
|
|
325 |
if (stat.ne.0) print*,'*** error allocating array lon ***' ! Grid -> lon
|
|
|
326 |
allocate(lat(nx,ny),stat=stat)
|
|
|
327 |
if (stat.ne.0) print*,'*** error allocating array lat ***' ! Gridy -> lat
|
|
|
328 |
|
|
|
329 |
C Get memory for trajectory arrays
|
|
|
330 |
allocate(trainp(ntra,1,ncol),stat=stat)
|
|
|
331 |
if (stat.ne.0) print*,'*** error allocating array trainp ***' ! Input start coordinates
|
|
|
332 |
allocate(traout(ntra,ntim,4),stat=stat)
|
|
|
333 |
if (stat.ne.0) print*,'*** error allocating array traout ***' ! Output trajectories
|
|
|
334 |
allocate(xx0(ntra),stat=stat)
|
|
|
335 |
if (stat.ne.0) print*,'*** error allocating array xx0 ***' ! X position (longitude)
|
|
|
336 |
allocate(yy0(ntra),stat=stat)
|
|
|
337 |
if (stat.ne.0) print*,'*** error allocating array yy0 ***' ! Y position (latitude)
|
|
|
338 |
allocate(pp0(ntra),stat=stat)
|
|
|
339 |
if (stat.ne.0) print*,'*** error allocating array pp0 ***' ! Pressure
|
|
|
340 |
allocate(leftflag(ntra),stat=stat)
|
|
|
341 |
if (stat.ne.0) print*,'*** error allocating array leftflag ***' ! Leaving-domain flag
|
|
|
342 |
|
21 |
michaesp |
343 |
c Get memory for potential temperature (isentropic trajectories)
|
|
|
344 |
if ( vertmode.eq.'isen' ) then
|
|
|
345 |
allocate(tht0(nx*ny*nz),stat=stat)
|
|
|
346 |
if (stat.ne.0) print*,'*** error allocating array tht0 ***'
|
|
|
347 |
allocate(tht1(nx*ny*nz),stat=stat)
|
|
|
348 |
if (stat.ne.0) print*,'*** error allocating array tht1 ***'
|
|
|
349 |
allocate(sth0(nx*ny),stat=stat)
|
|
|
350 |
if (stat.ne.0) print*,'*** error allocating array sth0 ***'
|
|
|
351 |
allocate(sth1(nx*ny),stat=stat)
|
|
|
352 |
if (stat.ne.0) print*,'*** error allocating array sth1 ***'
|
|
|
353 |
endif
|
|
|
354 |
|
2 |
michaesp |
355 |
c Write some status information
|
|
|
356 |
print*,'---- CONSTANT GRID PARAMETERS ---------------------------'
|
|
|
357 |
print*
|
|
|
358 |
print*,' xmin,xmax : ',xmin,xmax
|
|
|
359 |
print*,' ymin,ymax : ',ymin,ymax
|
|
|
360 |
print*,' dx,dy : ',dx,dy
|
|
|
361 |
print*,' nx,ny,nz : ',nx,ny,nz
|
|
|
362 |
print*,' per, hem : ',per,hem
|
|
|
363 |
print*
|
|
|
364 |
|
|
|
365 |
c --------------------------------------------------------------------
|
|
|
366 |
c Calculate the map scale factors
|
|
|
367 |
c --------------------------------------------------------------------
|
|
|
368 |
|
|
|
369 |
c Read lon/lat on the model grid from first data file
|
|
|
370 |
filename = prefix//dat(1)
|
|
|
371 |
call input_open (fid,filename)
|
|
|
372 |
varname='XLONG'
|
|
|
373 |
call input_var_wrf(fid,varname,lon,nx,ny,1)
|
|
|
374 |
varname='XLAT'
|
|
|
375 |
call input_var_wrf(fid,varname,lat,nx,ny,1)
|
|
|
376 |
call input_close(fid)
|
|
|
377 |
|
|
|
378 |
c Get map scale for interior points
|
|
|
379 |
do i=2,nx-1
|
|
|
380 |
do j=2,ny-1
|
|
|
381 |
|
|
|
382 |
c Map scale in x direction
|
|
|
383 |
lon1 = lon(i-1,j)
|
|
|
384 |
lat1 = lat(i-1,j)
|
|
|
385 |
lon2 = lon(i+1,j)
|
|
|
386 |
lat2 = lat(i+1,j)
|
|
|
387 |
radunit = 'km.haversine'
|
|
|
388 |
mpx(i,j) = 0.5 * 1000. * sdis(lon1,lat1,lon2,lat2,radunit)
|
|
|
389 |
|
|
|
390 |
c Map scale in y direction
|
|
|
391 |
lon1 = lon(i,j-1)
|
|
|
392 |
lat1 = lat(i,j-1)
|
|
|
393 |
lon2 = lon(i,j+1)
|
|
|
394 |
lat2 = lat(i,j+1)
|
|
|
395 |
radunit = 'km.haversine'
|
|
|
396 |
mpy(i,j) = 0.5 * 1000. * sdis(lon1,lat1,lon2,lat2,radunit)
|
|
|
397 |
|
|
|
398 |
enddo
|
|
|
399 |
enddo
|
|
|
400 |
|
|
|
401 |
c Copy map scale for boundary points
|
|
|
402 |
do i=2,nx-1
|
|
|
403 |
mpx(i, 1) = mpx(i, 2)
|
|
|
404 |
mpx(i,ny) = mpx(i,ny-1)
|
|
|
405 |
mpy(i, 1) = mpy(i, 2)
|
|
|
406 |
mpy(i,ny) = mpy(i,ny-1)
|
|
|
407 |
enddo
|
|
|
408 |
do j=2,ny-1
|
|
|
409 |
mpx(1, j) = mpx(2, j)
|
|
|
410 |
mpx(nx,j) = mpx(nx-1,j)
|
|
|
411 |
mpy(1, j) = mpy(2, j)
|
|
|
412 |
mpy(nx,j) = mpy(nx-1,j)
|
|
|
413 |
enddo
|
|
|
414 |
mpx(1,1) = mpx(2,2)
|
|
|
415 |
mpy(1,1) = mpy(2,2)
|
|
|
416 |
mpx(nx,1) = mpx(nx-1,2)
|
|
|
417 |
mpy(nx,1) = mpy(nx-1,2)
|
|
|
418 |
mpx(nx,ny) = mpx(nx-1,ny-1)
|
|
|
419 |
mpy(nx,ny) = mpy(nx-1,ny-1)
|
|
|
420 |
mpx(1,ny) = mpx(2,ny-1)
|
|
|
421 |
mpy(1,ny) = mpy(2,ny-1)
|
|
|
422 |
|
|
|
423 |
c Write some status information
|
15 |
michaesp |
424 |
print*,'---- MAPPING SCALE FACTORS -----------------------------'
|
2 |
michaesp |
425 |
print*
|
|
|
426 |
|
|
|
427 |
scale_max = mpx(1,1)
|
|
|
428 |
scale_min = mpx(1,1)
|
|
|
429 |
do i=1,nx
|
|
|
430 |
do j=1,ny
|
|
|
431 |
if ( mpx(i,j).gt.scale_max ) scale_max = mpx(i,j)
|
|
|
432 |
if ( mpx(i,j).lt.scale_min ) scale_min = mpx(i,j)
|
|
|
433 |
enddo
|
|
|
434 |
enddo
|
|
|
435 |
print*,' map scale x (min,max) : ',scale_min,scale_max
|
|
|
436 |
scale_max = mpy(1,1)
|
|
|
437 |
scale_min = mpy(1,1)
|
|
|
438 |
do i=1,nx
|
|
|
439 |
do j=1,ny
|
|
|
440 |
if ( mpy(i,j).gt.scale_max ) scale_max = mpy(i,j)
|
|
|
441 |
if ( mpy(i,j).lt.scale_min ) scale_min = mpy(i,j)
|
|
|
442 |
enddo
|
|
|
443 |
enddo
|
|
|
444 |
print*,' map scale y (min,max) : ',scale_min,scale_max
|
|
|
445 |
print*
|
|
|
446 |
|
|
|
447 |
c --------------------------------------------------------------------
|
|
|
448 |
c Initialize the trajectory calculation
|
|
|
449 |
c --------------------------------------------------------------------
|
|
|
450 |
|
|
|
451 |
c Read start coordinates from file - Format (x,y,lev)
|
|
|
452 |
if (inpmode.eq.-1) then
|
|
|
453 |
open(fid,file=strname)
|
|
|
454 |
do i=1,ntra
|
|
|
455 |
read(fid,*) xx0(i),yy0(i),pp0(i)
|
|
|
456 |
enddo
|
|
|
457 |
close(fid)
|
|
|
458 |
|
|
|
459 |
c Read start coordinates from trajectory file - check consistency of ref time
|
|
|
460 |
else
|
|
|
461 |
call ropen_tra(cdfid,strname,ntra,1,ncol,reftmp,vars,inpmode)
|
|
|
462 |
call read_tra (cdfid,trainp,ntra,1,ncol,inpmode)
|
|
|
463 |
do i=1,ntra
|
|
|
464 |
time = trainp(i,1,1)
|
|
|
465 |
xx0(i) = trainp(i,1,2)
|
|
|
466 |
yy0(i) = trainp(i,1,3)
|
|
|
467 |
pp0(i) = trainp(i,1,4)
|
|
|
468 |
enddo
|
|
|
469 |
call close_tra(cdfid,inpmode)
|
|
|
470 |
|
|
|
471 |
if ( ( vars(2).ne.'x').and.(vars(3).ne.'y') ) then
|
|
|
472 |
print*,' ERROR: starting positions must be in x/y/z format'
|
|
|
473 |
stop
|
|
|
474 |
endif
|
|
|
475 |
|
|
|
476 |
if ( ( reftime(1).ne.reftmp(1) ).or.
|
|
|
477 |
> ( reftime(2).ne.reftmp(2) ).or.
|
|
|
478 |
> ( reftime(3).ne.reftmp(3) ).or.
|
|
|
479 |
> ( reftime(4).ne.reftmp(4) ).or.
|
|
|
480 |
> ( reftime(5).ne.reftmp(5) ) )
|
|
|
481 |
> then
|
|
|
482 |
print*,' WARNING: Inconsistent reference times'
|
|
|
483 |
write(*,'(5i8)') (reftime(i),i=1,5)
|
|
|
484 |
write(*,'(5i8)') (reftmp (i),i=1,5)
|
21 |
michaesp |
485 |
print*,' Enter a key to proceed...'
|
2 |
michaesp |
486 |
stop
|
|
|
487 |
endif
|
|
|
488 |
endif
|
|
|
489 |
|
|
|
490 |
c Transform start coordinates from index space to WRF grid
|
|
|
491 |
do i=1,ntra
|
|
|
492 |
xx0(i) = xmin + ( xx0(i) - 1. ) * dx
|
|
|
493 |
yy0(i) = ymin + ( yy0(i) - 1. ) * dy
|
|
|
494 |
enddo
|
|
|
495 |
|
|
|
496 |
c Set sign of time range
|
|
|
497 |
reftime(6) = fbflag * reftime(6)
|
|
|
498 |
|
|
|
499 |
c Write some status information
|
|
|
500 |
print*,'---- REFERENCE DATE---------- ---------------------------'
|
|
|
501 |
print*
|
|
|
502 |
print*,' Reference time (year) :',reftime(1)
|
|
|
503 |
print*,' (month) :',reftime(2)
|
|
|
504 |
print*,' (day) :',reftime(3)
|
|
|
505 |
print*,' (hour) :',reftime(4)
|
|
|
506 |
print*,' (min) :',reftime(5)
|
|
|
507 |
print*,' Time range :',reftime(6),' min'
|
|
|
508 |
print*
|
|
|
509 |
|
21 |
michaesp |
510 |
C Save starting positions - vertical position will be corrected below for <aglconst> and <isen>
|
2 |
michaesp |
511 |
itim = 1
|
|
|
512 |
do i=1,ntra
|
|
|
513 |
traout(i,itim,1) = 0.
|
|
|
514 |
traout(i,itim,2) = xx0(i)
|
|
|
515 |
traout(i,itim,3) = yy0(i)
|
|
|
516 |
traout(i,itim,4) = pp0(i)
|
|
|
517 |
enddo
|
|
|
518 |
|
|
|
519 |
c Init the flag and the counter for trajectories leaving the domain
|
|
|
520 |
leftcount=0
|
|
|
521 |
do i=1,ntra
|
|
|
522 |
leftflag(i)=0
|
|
|
523 |
enddo
|
|
|
524 |
|
|
|
525 |
C Convert time shifts <tst,ten> from <hh.mm> into fractional time
|
|
|
526 |
call hhmm2frac(tst,frac)
|
|
|
527 |
tst = frac
|
|
|
528 |
call hhmm2frac(ten,frac)
|
|
|
529 |
ten = frac
|
|
|
530 |
|
21 |
michaesp |
531 |
c Adjust the starting heights in mode <aglconst>
|
|
|
532 |
if ( vertmode.eq.'aglconst' ) then
|
|
|
533 |
filename = prefix//dat(1)
|
|
|
534 |
call input_open (fid,filename)
|
|
|
535 |
varname='geopot' ! geopot.height
|
|
|
536 |
call input_var_wrf(fid,varname,p3t0,nx,ny,nz)
|
|
|
537 |
varname='geopots' ! geopot.height at surface
|
|
|
538 |
call input_var_wrf(fid,varname,spt0,nx,ny,1)
|
|
|
539 |
call input_close(fid)
|
|
|
540 |
filename = prefix//dat(2)
|
|
|
541 |
call input_open (fid,filename)
|
|
|
542 |
varname='geopot' ! geopot.height
|
|
|
543 |
call input_var_wrf(fid,varname,p3t1,nx,ny,nz)
|
|
|
544 |
varname='geopots' ! geopot.height at surface
|
|
|
545 |
call input_var_wrf(fid,varname,spt1,nx,ny,1)
|
|
|
546 |
call input_close(fid)
|
|
|
547 |
do i=1,ntra
|
|
|
548 |
call get_index4 (xind,yind,pind,xx0(i),yy0(i),0.,reltpos0,
|
|
|
549 |
> p3t0,p3t1,spt0,spt1,3,
|
|
|
550 |
> nx,ny,nz,xmin,ymin,dx,dy,mdv)
|
|
|
551 |
sp = int_index4 (spt0,spt1,nx,ny,1,xind,yind,1.,reltpos0,mdv)
|
|
|
552 |
traout(i,1,4) = sp + aglconst
|
|
|
553 |
enddo
|
|
|
554 |
endif
|
|
|
555 |
|
|
|
556 |
c Adjust the starting heights in mode <isen>
|
|
|
557 |
if ( vertmode.eq.'isen' ) then
|
|
|
558 |
filename = prefix//dat(1)
|
|
|
559 |
call input_open (fid,filename)
|
|
|
560 |
varname='geopot' ! geopot.height
|
|
|
561 |
call input_var_wrf(fid,varname,p3t0,nx,ny,nz)
|
|
|
562 |
varname='geopots' ! geopot.height at surface
|
|
|
563 |
call input_var_wrf(fid,varname,spt0,nx,ny,1)
|
|
|
564 |
varname='TH' ! TH
|
|
|
565 |
call input_var_wrf(fid,varname,tht0,nx,ny,nz)
|
|
|
566 |
do i=1,nx*ny
|
|
|
567 |
sth0(i) = tht0(i)
|
|
|
568 |
enddo
|
|
|
569 |
call input_close(fid)
|
|
|
570 |
filename = prefix//dat(2)
|
|
|
571 |
call input_open (fid,filename)
|
|
|
572 |
varname='geopot' ! geopot.height
|
|
|
573 |
call input_var_wrf(fid,varname,p3t1,nx,ny,nz)
|
|
|
574 |
varname='geopots' ! geopot.height at surface
|
|
|
575 |
call input_var_wrf(fid,varname,spt1,nx,ny,1)
|
|
|
576 |
varname='TH' ! TH
|
|
|
577 |
call input_var_wrf(fid,varname,tht1,nx,ny,nz)
|
|
|
578 |
do i=1,nx*ny
|
|
|
579 |
sth1(i) = tht1(i)
|
|
|
580 |
enddo
|
|
|
581 |
call input_close(fid)
|
|
|
582 |
do i=1,ntra
|
|
|
583 |
|
|
|
584 |
call get_index4 (xind,yind,pind,xx0(i),yy0(i),isen,reltpos0,
|
|
|
585 |
> tht0,tht1,sth0,sth1,1,
|
|
|
586 |
> nx,ny,nz,xmin,ymin,dx,dy,mdv)
|
|
|
587 |
traout(i,1,4) = int_index4
|
|
|
588 |
> (p3t0,p3t1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
|
|
|
589 |
enddo
|
|
|
590 |
endif
|
|
|
591 |
|
2 |
michaesp |
592 |
c -----------------------------------------------------------------------
|
|
|
593 |
c Loop to calculate trajectories
|
|
|
594 |
c -----------------------------------------------------------------------
|
|
|
595 |
|
|
|
596 |
c Write some status information
|
|
|
597 |
print*,'---- TRAJECTORIES ----------- ---------------------------'
|
|
|
598 |
print*
|
|
|
599 |
|
|
|
600 |
C Set the time for the first data file (depending on forward/backward mode)
|
|
|
601 |
if (fbflag.eq.1) then
|
|
|
602 |
tstart = -tst
|
|
|
603 |
else
|
|
|
604 |
tstart = tst
|
|
|
605 |
endif
|
|
|
606 |
|
|
|
607 |
c Set the minute counter for output
|
|
|
608 |
wstep = 0
|
|
|
609 |
|
|
|
610 |
c Read wind fields and vertical grid from first file
|
|
|
611 |
filename = prefix//dat(1)
|
|
|
612 |
|
|
|
613 |
call frac2hhmm(tstart,tload)
|
|
|
614 |
|
25 |
michaesp |
615 |
write(*,'(a16,a20,2f7.2)') ' (file,time) : ',
|
2 |
michaesp |
616 |
> trim(filename),tload
|
|
|
617 |
|
|
|
618 |
call input_open (fid,filename)
|
|
|
619 |
|
|
|
620 |
varname='U' ! U
|
|
|
621 |
call input_var_wrf(fid,varname,uut1,nx,ny,nz)
|
|
|
622 |
|
|
|
623 |
varname='V' ! V
|
|
|
624 |
call input_var_wrf(fid,varname,vvt1,nx,ny,nz)
|
|
|
625 |
|
|
|
626 |
varname='W' ! W
|
|
|
627 |
call input_var_wrf(fid,varname,wwt1,nx,ny,nz)
|
|
|
628 |
|
|
|
629 |
varname='geopot' ! geopot.height
|
|
|
630 |
call input_var_wrf(fid,varname,p3t1,nx,ny,nz)
|
|
|
631 |
|
|
|
632 |
varname='geopots' ! geopot.height at surface
|
|
|
633 |
call input_var_wrf(fid,varname,spt1,nx,ny,1)
|
|
|
634 |
|
21 |
michaesp |
635 |
if (vertmode.eq.'isen') then
|
|
|
636 |
varname='TH' ! TH
|
|
|
637 |
call input_var_wrf(fid,varname,tht1,nx,ny,nz)
|
|
|
638 |
do i=1,nx*ny
|
|
|
639 |
sth1(i) = tht1(i)
|
|
|
640 |
enddo
|
|
|
641 |
endif
|
|
|
642 |
|
2 |
michaesp |
643 |
call input_close(fid)
|
|
|
644 |
|
|
|
645 |
c Loop over all input files (time step is <timeinc>)
|
|
|
646 |
do itm=1,numdat-1
|
|
|
647 |
|
|
|
648 |
c Calculate actual and next time
|
|
|
649 |
time0 = tstart+real(itm-1)*timeinc*fbflag
|
|
|
650 |
time1 = time0+timeinc*fbflag
|
|
|
651 |
|
|
|
652 |
c Copy old velocities and pressure fields to new ones
|
|
|
653 |
do i=1,nx*ny*nz
|
|
|
654 |
uut0(i)=uut1(i)
|
|
|
655 |
vvt0(i)=vvt1(i)
|
|
|
656 |
wwt0(i)=wwt1(i)
|
|
|
657 |
p3t0(i)=p3t1(i)
|
|
|
658 |
enddo
|
|
|
659 |
do i=1,nx*ny
|
|
|
660 |
spt0(i)=spt1(i)
|
|
|
661 |
enddo
|
|
|
662 |
|
21 |
michaesp |
663 |
c Copy potential temperature
|
|
|
664 |
if ( vertmode.eq.'isen' ) then
|
|
|
665 |
do i=1,nx*ny*nz
|
|
|
666 |
tht0(i)=tht1(i)
|
|
|
667 |
enddo
|
|
|
668 |
do i=1,nx*ny
|
|
|
669 |
sth0(i)=sth1(i)
|
|
|
670 |
enddo
|
|
|
671 |
endif
|
|
|
672 |
|
2 |
michaesp |
673 |
c Read wind fields and surface pressure at next time
|
|
|
674 |
filename = prefix//dat(itm+1)
|
|
|
675 |
|
|
|
676 |
call frac2hhmm(time1,tload)
|
25 |
michaesp |
677 |
write(*,'(a16,a20,2f7.2)') ' (file,time) : ',
|
2 |
michaesp |
678 |
> trim(filename),tload
|
|
|
679 |
|
|
|
680 |
call input_open (fid,filename)
|
|
|
681 |
|
|
|
682 |
varname='U' ! U
|
|
|
683 |
call input_var_wrf(fid,varname,uut1,nx,ny,nz)
|
|
|
684 |
|
|
|
685 |
varname='V' ! V
|
|
|
686 |
call input_var_wrf(fid,varname,vvt1,nx,ny,nz)
|
|
|
687 |
|
|
|
688 |
varname='W' ! W
|
|
|
689 |
call input_var_wrf(fid,varname,wwt1,nx,ny,nz)
|
|
|
690 |
|
|
|
691 |
varname='geopot' ! geopot.height
|
|
|
692 |
call input_var_wrf(fid,varname,p3t1,nx,ny,nz)
|
|
|
693 |
|
|
|
694 |
varname='geopots' ! geopot.height
|
|
|
695 |
call input_var_wrf(fid,varname,spt1,nx,ny,1)
|
21 |
michaesp |
696 |
|
|
|
697 |
if (vertmode.eq.'isen') then
|
|
|
698 |
varname='TH' ! TH
|
|
|
699 |
call input_var_wrf(fid,varname,tht1,nx,ny,nz)
|
|
|
700 |
do i=1,nx*ny
|
|
|
701 |
sth1(i) = tht1(i)
|
|
|
702 |
enddo
|
|
|
703 |
endif
|
|
|
704 |
|
2 |
michaesp |
705 |
call input_close(fid)
|
|
|
706 |
|
|
|
707 |
C Determine the first and last loop indices
|
|
|
708 |
if (numdat.eq.2) then
|
|
|
709 |
filo = nint(tst/ts)+1
|
|
|
710 |
lalo = nint((timeinc-ten)/ts)
|
|
|
711 |
elseif ( itm.eq.1 ) then
|
|
|
712 |
filo = nint(tst/ts)+1
|
|
|
713 |
lalo = nint(timeinc/ts)
|
|
|
714 |
else if (itm.eq.numdat-1) then
|
|
|
715 |
filo = 1
|
|
|
716 |
lalo = nint((timeinc-ten)/ts)
|
|
|
717 |
else
|
|
|
718 |
filo = 1
|
|
|
719 |
lalo = nint(timeinc/ts)
|
|
|
720 |
endif
|
|
|
721 |
|
|
|
722 |
c Split the interval <timeinc> into computational time steps <ts>
|
|
|
723 |
do iloop=filo,lalo
|
|
|
724 |
|
|
|
725 |
C Calculate relative time position in the interval timeinc (0=beginning, 1=end)
|
|
|
726 |
reltpos0 = ((real(iloop)-1.)*ts)/timeinc
|
|
|
727 |
reltpos1 = real(iloop)*ts/timeinc
|
|
|
728 |
|
|
|
729 |
C Initialize counter for domain leaving trajectories
|
|
|
730 |
leftcount=0
|
|
|
731 |
|
|
|
732 |
c Timestep for all trajectories
|
21 |
michaesp |
733 |
do i=1,ntra
|
2 |
michaesp |
734 |
|
|
|
735 |
C Check if trajectory has already left the data domain
|
|
|
736 |
if (leftflag(i).ne.1) then
|
|
|
737 |
|
|
|
738 |
c Iterative Euler timestep (x0,y0,p0 -> x1,y1,p1)
|
21 |
michaesp |
739 |
if ( vertmode.eq.'kinematic' ) then
|
|
|
740 |
call euler_kinematic(
|
2 |
michaesp |
741 |
> xx1,yy1,pp1,leftflag(i),
|
|
|
742 |
> xx0(i),yy0(i),pp0(i),reltpos0,reltpos1,
|
|
|
743 |
> ts*3600,numit,jflag,mdv,wfactor,fbflag,
|
|
|
744 |
> spt0,spt1,p3t0,p3t1,uut0,uut1,vvt0,vvt1,wwt0,wwt1,
|
|
|
745 |
> xmin,ymin,dx,dy,per,hem,nx,ny,nz,mpx,mpy)
|
21 |
michaesp |
746 |
elseif ( vertmode.eq.'aglconst' ) then
|
|
|
747 |
call euler_aglconst(
|
|
|
748 |
> xx1,yy1,pp1,leftflag(i),
|
|
|
749 |
> xx0(i),yy0(i),pp0(i),reltpos0,reltpos1,
|
|
|
750 |
> ts*3600,numit,jflag,mdv,wfactor,fbflag,
|
|
|
751 |
> spt0,spt1,p3t0,p3t1,uut0,uut1,vvt0,vvt1,
|
|
|
752 |
> xmin,ymin,dx,dy,per,hem,nx,ny,nz,mpx,mpy,aglconst)
|
|
|
753 |
elseif ( vertmode.eq.'isen' ) then
|
|
|
754 |
call euler_isen(
|
|
|
755 |
> xx1,yy1,pp1,leftflag(i),
|
|
|
756 |
> xx0(i),yy0(i),pp0(i),reltpos0,reltpos1,
|
|
|
757 |
> ts*3600,numit,jflag,mdv,wfactor,fbflag,
|
|
|
758 |
> spt0,spt1,p3t0,p3t1,uut0,uut1,vvt0,vvt1,
|
|
|
759 |
> tht0,tht1,sth0,sth1,
|
|
|
760 |
> xmin,ymin,dx,dy,per,hem,nx,ny,nz,mpx,mpy,isen)
|
|
|
761 |
|
|
|
762 |
endif
|
2 |
michaesp |
763 |
|
|
|
764 |
c Update trajectory position, or increase number of trajectories leaving domain
|
|
|
765 |
if (leftflag(i).eq.1) then
|
|
|
766 |
leftcount=leftcount+1
|
|
|
767 |
if ( leftcount.lt.10 ) then
|
|
|
768 |
print*,' -> Trajectory ',i,' leaves domain'
|
|
|
769 |
elseif ( leftcount.eq.10 ) then
|
|
|
770 |
print*,' -> N>=10 trajectories leave domain'
|
|
|
771 |
endif
|
|
|
772 |
else
|
|
|
773 |
xx0(i)=xx1
|
|
|
774 |
yy0(i)=yy1
|
|
|
775 |
pp0(i)=pp1
|
|
|
776 |
endif
|
|
|
777 |
|
|
|
778 |
c Trajectory has already left data domain (mark as <mdv>)
|
|
|
779 |
else
|
|
|
780 |
xx0(i)=mdv
|
|
|
781 |
yy0(i)=mdv
|
|
|
782 |
pp0(i)=mdv
|
|
|
783 |
|
|
|
784 |
endif
|
|
|
785 |
|
|
|
786 |
enddo
|
|
|
787 |
|
|
|
788 |
C Save positions only every deltout minutes
|
|
|
789 |
delta = aint(iloop*60*ts/deltout)-iloop*60*ts/deltout
|
|
|
790 |
if (abs(delta).lt.eps) then
|
|
|
791 |
time = time0+reltpos1*timeinc*fbflag
|
|
|
792 |
itim = itim + 1
|
|
|
793 |
do i=1,ntra
|
|
|
794 |
call frac2hhmm(time,tload)
|
|
|
795 |
traout(i,itim,1) = tload
|
|
|
796 |
traout(i,itim,2) = xx0(i)
|
|
|
797 |
traout(i,itim,3) = yy0(i)
|
|
|
798 |
traout(i,itim,4) = pp0(i)
|
|
|
799 |
enddo
|
|
|
800 |
endif
|
|
|
801 |
|
|
|
802 |
enddo
|
|
|
803 |
|
|
|
804 |
enddo
|
|
|
805 |
|
|
|
806 |
c Transform trajectory position from WRF grid do grid index
|
|
|
807 |
do i=1,ntra
|
|
|
808 |
do j=1,ntim
|
19 |
michaesp |
809 |
if ( abs(traout(i,j,2)-mdv).gt.eps ) then
|
|
|
810 |
traout(i,j,2) = ( traout(i,j,2) - xmin ) / dx + 1.
|
|
|
811 |
else
|
|
|
812 |
traout(i,j,2) = mdv
|
|
|
813 |
endif
|
|
|
814 |
if ( abs(traout(i,j,3)-mdv).gt.eps ) then
|
|
|
815 |
traout(i,j,3) = ( traout(i,j,3) - ymin ) / dy + 1.
|
|
|
816 |
else
|
|
|
817 |
traout(i,j,3) = mdv
|
|
|
818 |
endif
|
2 |
michaesp |
819 |
enddo
|
|
|
820 |
enddo
|
|
|
821 |
|
|
|
822 |
c Write trajectory file
|
|
|
823 |
vars(1) ='time'
|
|
|
824 |
vars(2) ='x'
|
|
|
825 |
vars(3) ='y'
|
|
|
826 |
vars(4) ='z'
|
|
|
827 |
call wopen_tra(cdfid,cdfname,ntra,ntim,4,reftime,vars,outmode)
|
|
|
828 |
call write_tra(cdfid,traout,ntra,ntim,4,outmode)
|
|
|
829 |
call close_tra(cdfid,outmode)
|
|
|
830 |
|
|
|
831 |
c Write some status information, and end of program message
|
|
|
832 |
print*
|
|
|
833 |
print*,'---- STATUS INFORMATION --------------------------------'
|
|
|
834 |
print*
|
|
|
835 |
print*,' #leaving domain ', leftcount
|
|
|
836 |
print*,' #staying in domain ', ntra-leftcount
|
|
|
837 |
print*
|
|
|
838 |
print*,' *** END OF PROGRAM CALTRA ***'
|
|
|
839 |
print*,'========================================================='
|
|
|
840 |
|
|
|
841 |
stop
|
|
|
842 |
|
|
|
843 |
c ------------------------------------------------------------------
|
|
|
844 |
c Exception handling
|
|
|
845 |
c ------------------------------------------------------------------
|
|
|
846 |
|
|
|
847 |
991 write(*,*) '*** ERROR: all start points outside the data domain'
|
|
|
848 |
call exit(1)
|
|
|
849 |
|
|
|
850 |
992 write(*,*) '*** ERROR: close arrays on files (prog. closear)'
|
|
|
851 |
call exit(1)
|
|
|
852 |
|
|
|
853 |
993 write(*,*) '*** ERROR: problems with array size'
|
|
|
854 |
call exit(1)
|
|
|
855 |
|
|
|
856 |
end
|
|
|
857 |
|
|
|
858 |
|
|
|
859 |
c *******************************************************************
|
|
|
860 |
c * Time step : either Euler or Runge-Kutta *
|
|
|
861 |
c *******************************************************************
|
|
|
862 |
|
|
|
863 |
C Time-step from (x0,y0,p0) to (x1,y1,p1)
|
|
|
864 |
C
|
|
|
865 |
C (x0,y0,p0) input coordinates (long,lat,p) for starting point
|
|
|
866 |
C (x1,y1,p1) output coordinates (long,lat,p) for end point
|
|
|
867 |
C deltat input timestep in seconds
|
|
|
868 |
C numit input number of iterations
|
|
|
869 |
C jump input flag (=1 trajectories don't enter the ground)
|
|
|
870 |
C left output flag (=1 if trajectory leaves data domain)
|
|
|
871 |
|
|
|
872 |
c -------------------------------------------------------------------
|
21 |
michaesp |
873 |
c Iterative Euler time step - kinematic trajectory
|
2 |
michaesp |
874 |
c -------------------------------------------------------------------
|
|
|
875 |
|
21 |
michaesp |
876 |
subroutine euler_kinematic
|
|
|
877 |
> (x1,y1,p1,left,x0,y0,p0,reltpos0,reltpos1,
|
|
|
878 |
> deltat,numit,jump,mdv,wfactor,fbflag,
|
|
|
879 |
> spt0,spt1,p3d0,p3d1,uut0,uut1,vvt0,vvt1,wwt0,wwt1,
|
|
|
880 |
> xmin,ymin,dx,dy,per,hem,nx,ny,nz,mpx,mpy)
|
2 |
michaesp |
881 |
|
|
|
882 |
implicit none
|
|
|
883 |
|
11 |
michaesp |
884 |
c Flag for test mode
|
|
|
885 |
integer test
|
|
|
886 |
parameter (test=0)
|
|
|
887 |
|
2 |
michaesp |
888 |
c Declaration of subroutine parameters
|
|
|
889 |
integer nx,ny,nz
|
|
|
890 |
real x1,y1,p1
|
|
|
891 |
integer left
|
|
|
892 |
real x0,y0,p0
|
|
|
893 |
real reltpos0,reltpos1
|
|
|
894 |
real deltat
|
|
|
895 |
integer numit
|
|
|
896 |
integer jump
|
|
|
897 |
real wfactor
|
|
|
898 |
integer fbflag
|
|
|
899 |
real spt0(nx*ny) ,spt1(nx*ny)
|
|
|
900 |
real uut0(nx*ny*nz),uut1(nx*ny*nz)
|
|
|
901 |
real vvt0(nx*ny*nz),vvt1(nx*ny*nz)
|
|
|
902 |
real wwt0(nx*ny*nz),wwt1(nx*ny*nz)
|
|
|
903 |
real p3d0(nx*ny*nz),p3d1(nx*ny*nz)
|
|
|
904 |
real xmin,ymin,dx,dy
|
|
|
905 |
integer per
|
|
|
906 |
integer hem
|
|
|
907 |
real mdv
|
|
|
908 |
real mpx(nx*ny),mpy(nx*ny)
|
|
|
909 |
|
|
|
910 |
c Auxiliary variables
|
|
|
911 |
real xmax,ymax
|
|
|
912 |
real xind,yind,pind
|
|
|
913 |
|
|
|
914 |
real u0,v0,w0,u1,v1,w1,u,v,w,sp
|
|
|
915 |
integer icount
|
|
|
916 |
character ch
|
|
|
917 |
real mpsc_x,mpsc_y
|
|
|
918 |
|
|
|
919 |
c Externals
|
|
|
920 |
real int_index4
|
|
|
921 |
external int_index4
|
|
|
922 |
|
|
|
923 |
c Reset the flag for domain-leaving
|
|
|
924 |
left=0
|
|
|
925 |
|
|
|
926 |
c Set the esat-north boundary of the domain
|
|
|
927 |
xmax = xmin+real(nx-1)*dx
|
|
|
928 |
ymax = ymin+real(ny-1)*dy
|
|
|
929 |
|
|
|
930 |
C Interpolate wind fields to starting position (x0,y0,p0)
|
|
|
931 |
call get_index4 (xind,yind,pind,x0,y0,p0,reltpos0,
|
|
|
932 |
> p3d0,p3d1,spt0,spt1,3,
|
|
|
933 |
> nx,ny,nz,xmin,ymin,dx,dy,mdv)
|
|
|
934 |
|
|
|
935 |
u0 = int_index4(uut0,uut1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
|
|
|
936 |
v0 = int_index4(vvt0,vvt1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
|
|
|
937 |
w0 = int_index4(wwt0,wwt1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
|
|
|
938 |
|
|
|
939 |
c Force the near-surface wind to zero
|
|
|
940 |
if (pind.lt.1.) w0=w0*pind
|
|
|
941 |
|
|
|
942 |
C For first iteration take ending position equal to starting position
|
|
|
943 |
x1=x0
|
|
|
944 |
y1=y0
|
|
|
945 |
p1=p0
|
|
|
946 |
|
|
|
947 |
C Iterative calculation of new position
|
|
|
948 |
do icount=1,numit
|
|
|
949 |
|
|
|
950 |
C Calculate new winds for advection
|
|
|
951 |
call get_index4 (xind,yind,pind,x1,y1,p1,reltpos1,
|
|
|
952 |
> p3d0,p3d1,spt0,spt1,3,
|
|
|
953 |
> nx,ny,nz,xmin,ymin,dx,dy,mdv)
|
|
|
954 |
u1 = int_index4(uut0,uut1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
|
|
|
955 |
v1 = int_index4(vvt0,vvt1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
|
|
|
956 |
w1 = int_index4(wwt0,wwt1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
|
|
|
957 |
|
|
|
958 |
c Force the near-surface wind to zero
|
|
|
959 |
if (pind.lt.1.) w1=w1*pind
|
|
|
960 |
|
|
|
961 |
c Get the new velocity in between
|
|
|
962 |
u=(u0+u1)/2.
|
|
|
963 |
v=(v0+v1)/2.
|
|
|
964 |
w=(w0+w1)/2.
|
|
|
965 |
|
|
|
966 |
c Get the mapping scale factors for this position
|
|
|
967 |
mpsc_x = int_index4 (mpx,mpx,nx,ny,1,xind,yind,1.,reltpos1,mdv)
|
|
|
968 |
mpsc_y = int_index4 (mpy,mpy,nx,ny,1,xind,yind,1.,reltpos1,mdv)
|
|
|
969 |
|
|
|
970 |
C Calculate new positions (adapted for cartesian grid)
|
|
|
971 |
x1 = x0 + fbflag * deltat * u * dx/mpsc_x
|
|
|
972 |
y1 = y0 + fbflag * deltat * v * dy/mpsc_y
|
|
|
973 |
p1 = p0 + fbflag * deltat * w * wfactor
|
|
|
974 |
|
11 |
michaesp |
975 |
c Write the update positions for tests
|
|
|
976 |
if ( (icount.eq.3).and.(test.eq.1) ) then
|
|
|
977 |
write(*,'(10f10.2)') x0,u,x1-x0,p0,w,p1-p0
|
|
|
978 |
endif
|
|
|
979 |
|
2 |
michaesp |
980 |
C Interpolate surface geop.height to actual position
|
|
|
981 |
call get_index4 (xind,yind,pind,x1,y1,0.,reltpos1,
|
|
|
982 |
> p3d0,p3d1,spt0,spt1,3,
|
|
|
983 |
> nx,ny,nz,xmin,ymin,dx,dy,mdv)
|
|
|
984 |
sp = int_index4 (spt0,spt1,nx,ny,1,xind,yind,1.,reltpos1,mdv)
|
|
|
985 |
|
|
|
986 |
c Handle trajectories which cross the lower boundary (jump flag)
|
|
|
987 |
if ((jump.eq.1).and.(p1.lt.sp)) p1=sp+10.
|
|
|
988 |
|
|
|
989 |
c Handle periodioc boundaries for 'ideal' mode
|
|
|
990 |
if ( per.eq.1 ) then
|
|
|
991 |
if ( x1.gt.xmax ) x1=x1-xmax
|
|
|
992 |
if ( x1.lt.xmin ) x1=x1+xmax
|
|
|
993 |
endif
|
|
|
994 |
|
|
|
995 |
C Check if trajectory leaves data domain
|
|
|
996 |
if ( (x1.lt.xmin).or.(x1.gt.xmax).or.
|
|
|
997 |
> (y1.lt.ymin).or.(y1.gt.ymax).or.(p1.lt.sp) ) ! geopot : .lt. ; pressure: .gt.
|
|
|
998 |
> then
|
|
|
999 |
left=1
|
|
|
1000 |
print*,x1,y1,p1
|
|
|
1001 |
print*,xmin,ymin
|
|
|
1002 |
print*,xmax,ymax
|
|
|
1003 |
print*,sp
|
|
|
1004 |
goto 100
|
|
|
1005 |
endif
|
|
|
1006 |
|
|
|
1007 |
enddo
|
|
|
1008 |
|
|
|
1009 |
c Exit point for subroutine
|
|
|
1010 |
100 continue
|
|
|
1011 |
|
|
|
1012 |
return
|
|
|
1013 |
|
|
|
1014 |
end
|
|
|
1015 |
|
21 |
michaesp |
1016 |
c -------------------------------------------------------------------
|
|
|
1017 |
c Iterative Euler time step - constant level above ground
|
|
|
1018 |
c -------------------------------------------------------------------
|
2 |
michaesp |
1019 |
|
21 |
michaesp |
1020 |
subroutine euler_aglconst
|
|
|
1021 |
> (x1,y1,p1,left,x0,y0,p0,reltpos0,reltpos1,
|
|
|
1022 |
> deltat,numit,jump,mdv,wfactor,fbflag,
|
|
|
1023 |
> spt0,spt1,p3d0,p3d1,uut0,uut1,vvt0,vvt1,
|
|
|
1024 |
> xmin,ymin,dx,dy,per,hem,nx,ny,nz,mpx,mpy,aglconst)
|
|
|
1025 |
|
|
|
1026 |
implicit none
|
|
|
1027 |
|
|
|
1028 |
c Flag for test mode
|
|
|
1029 |
integer test
|
|
|
1030 |
parameter (test=0)
|
|
|
1031 |
|
|
|
1032 |
c Declaration of subroutine parameters
|
|
|
1033 |
integer nx,ny,nz
|
|
|
1034 |
real x1,y1,p1
|
|
|
1035 |
integer left
|
|
|
1036 |
real x0,y0,p0
|
|
|
1037 |
real reltpos0,reltpos1
|
|
|
1038 |
real deltat
|
|
|
1039 |
integer numit
|
|
|
1040 |
integer jump
|
|
|
1041 |
real wfactor
|
|
|
1042 |
integer fbflag
|
|
|
1043 |
real spt0(nx*ny) ,spt1(nx*ny)
|
|
|
1044 |
real uut0(nx*ny*nz),uut1(nx*ny*nz)
|
|
|
1045 |
real vvt0(nx*ny*nz),vvt1(nx*ny*nz)
|
|
|
1046 |
real p3d0(nx*ny*nz),p3d1(nx*ny*nz)
|
|
|
1047 |
real xmin,ymin,dx,dy
|
|
|
1048 |
integer per
|
|
|
1049 |
integer hem
|
|
|
1050 |
real mdv
|
|
|
1051 |
real mpx(nx*ny),mpy(nx*ny)
|
|
|
1052 |
real aglconst
|
|
|
1053 |
|
|
|
1054 |
c Auxiliary variables
|
|
|
1055 |
real xmax,ymax
|
|
|
1056 |
real xind,yind,pind
|
|
|
1057 |
|
|
|
1058 |
real u0,v0,u1,v1,u,v,w,sp
|
|
|
1059 |
integer icount
|
|
|
1060 |
character ch
|
|
|
1061 |
real mpsc_x,mpsc_y
|
|
|
1062 |
|
|
|
1063 |
c Externals
|
|
|
1064 |
real int_index4
|
|
|
1065 |
external int_index4
|
|
|
1066 |
|
|
|
1067 |
c Reset the flag for domain-leaving
|
|
|
1068 |
left=0
|
|
|
1069 |
|
|
|
1070 |
c Set the esat-north boundary of the domain
|
|
|
1071 |
xmax = xmin+real(nx-1)*dx
|
|
|
1072 |
ymax = ymin+real(ny-1)*dy
|
|
|
1073 |
|
|
|
1074 |
C Set parcel height to aglconst above topography
|
|
|
1075 |
call get_index4 (xind,yind,pind,x1,y1,0.,reltpos0,
|
|
|
1076 |
> p3d0,p3d1,spt0,spt1,3,
|
|
|
1077 |
> nx,ny,nz,xmin,ymin,dx,dy,mdv)
|
|
|
1078 |
sp = int_index4 (spt0,spt1,nx,ny,1,xind,yind,1.,reltpos0,mdv)
|
|
|
1079 |
p0 = sp + aglconst
|
|
|
1080 |
|
|
|
1081 |
C Interpolate wind fields to starting position (x0,y0,p0)
|
|
|
1082 |
call get_index4 (xind,yind,pind,x0,y0,p0,reltpos0,
|
|
|
1083 |
> p3d0,p3d1,spt0,spt1,3,
|
|
|
1084 |
> nx,ny,nz,xmin,ymin,dx,dy,mdv)
|
|
|
1085 |
u0 = int_index4(uut0,uut1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
|
|
|
1086 |
v0 = int_index4(vvt0,vvt1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
|
|
|
1087 |
|
|
|
1088 |
C For first iteration take ending position equal to starting position
|
|
|
1089 |
x1 = x0
|
|
|
1090 |
y1 = y0
|
|
|
1091 |
p1 = p0
|
|
|
1092 |
|
|
|
1093 |
C Iterative calculation of new position
|
|
|
1094 |
do icount=1,numit
|
|
|
1095 |
|
|
|
1096 |
C Calculate new winds for advection
|
|
|
1097 |
call get_index4 (xind,yind,pind,x1,y1,p1,reltpos1,
|
|
|
1098 |
> p3d0,p3d1,spt0,spt1,3,
|
|
|
1099 |
> nx,ny,nz,xmin,ymin,dx,dy,mdv)
|
|
|
1100 |
u1 = int_index4(uut0,uut1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
|
|
|
1101 |
v1 = int_index4(vvt0,vvt1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
|
|
|
1102 |
|
|
|
1103 |
c Get the new velocity in between
|
|
|
1104 |
u=(u0+u1)/2.
|
|
|
1105 |
v=(v0+v1)/2.
|
|
|
1106 |
|
|
|
1107 |
c Get the mapping scale factors for this position
|
|
|
1108 |
mpsc_x = int_index4 (mpx,mpx,nx,ny,1,xind,yind,1.,reltpos1,mdv)
|
|
|
1109 |
mpsc_y = int_index4 (mpy,mpy,nx,ny,1,xind,yind,1.,reltpos1,mdv)
|
|
|
1110 |
|
|
|
1111 |
C Calculate new positions (adapted for cartesian grid)
|
|
|
1112 |
x1 = x0 + fbflag * deltat * u * dx/mpsc_x
|
|
|
1113 |
y1 = y0 + fbflag * deltat * v * dy/mpsc_y
|
|
|
1114 |
|
|
|
1115 |
c Write the update positions for tests
|
|
|
1116 |
if ( (icount.eq.3).and.(test.eq.1) ) then
|
|
|
1117 |
write(*,'(10f10.2)') x0,u,x1-x0,p0,w,p1-p0
|
|
|
1118 |
endif
|
|
|
1119 |
|
|
|
1120 |
C Set trajectory height to aglconst above topography
|
|
|
1121 |
call get_index4 (xind,yind,pind,x1,y1,0.,reltpos1,
|
|
|
1122 |
> p3d0,p3d1,spt0,spt1,3,
|
|
|
1123 |
> nx,ny,nz,xmin,ymin,dx,dy,mdv)
|
|
|
1124 |
sp = int_index4 (spt0,spt1,nx,ny,1,xind,yind,1.,reltpos1,mdv)
|
|
|
1125 |
p1 = sp + aglconst
|
|
|
1126 |
|
|
|
1127 |
c Handle trajectories which cross the lower boundary (jump flag)
|
|
|
1128 |
if ((jump.eq.1).and.(p1.lt.sp)) then
|
|
|
1129 |
print*,'jump'
|
|
|
1130 |
p1=sp+10.
|
|
|
1131 |
endif
|
|
|
1132 |
|
|
|
1133 |
c Handle periodioc boundaries for 'ideal' mode
|
|
|
1134 |
if ( per.eq.1 ) then
|
|
|
1135 |
if ( x1.gt.xmax ) x1=x1-xmax
|
|
|
1136 |
if ( x1.lt.xmin ) x1=x1+xmax
|
|
|
1137 |
endif
|
|
|
1138 |
|
|
|
1139 |
C Check if trajectory leaves data domain
|
|
|
1140 |
if ( (x1.lt.xmin).or.(x1.gt.xmax).or.
|
|
|
1141 |
> (y1.lt.ymin).or.(y1.gt.ymax).or.(p1.lt.sp) ) ! geopot : .lt. ; pressure: .gt.
|
|
|
1142 |
> then
|
|
|
1143 |
left=1
|
|
|
1144 |
goto 100
|
|
|
1145 |
endif
|
|
|
1146 |
|
|
|
1147 |
enddo
|
|
|
1148 |
|
|
|
1149 |
c Exit point for subroutine
|
|
|
1150 |
100 continue
|
|
|
1151 |
|
|
|
1152 |
return
|
|
|
1153 |
|
|
|
1154 |
end
|
|
|
1155 |
|
|
|
1156 |
c -------------------------------------------------------------------
|
|
|
1157 |
c Iterative Euler time step - isentropic
|
|
|
1158 |
c -------------------------------------------------------------------
|
|
|
1159 |
|
|
|
1160 |
subroutine euler_isen
|
|
|
1161 |
> (x1,y1,p1,left,x0,y0,p0,reltpos0,reltpos1,
|
|
|
1162 |
> deltat,numit,jump,mdv,wfactor,fbflag,
|
|
|
1163 |
> spt0,spt1,p3d0,p3d1,uut0,uut1,vvt0,vvt1,
|
|
|
1164 |
> tht0,tht1,sth0,sth1,
|
|
|
1165 |
> xmin,ymin,dx,dy,per,hem,nx,ny,nz,mpx,mpy,isen)
|
|
|
1166 |
|
|
|
1167 |
implicit none
|
|
|
1168 |
|
|
|
1169 |
c Flag for test mode
|
|
|
1170 |
integer test
|
|
|
1171 |
parameter (test=0)
|
|
|
1172 |
|
|
|
1173 |
c Declaration of subroutine parameters
|
|
|
1174 |
integer nx,ny,nz
|
|
|
1175 |
real x1,y1,p1
|
|
|
1176 |
integer left
|
|
|
1177 |
real x0,y0,p0
|
|
|
1178 |
real reltpos0,reltpos1
|
|
|
1179 |
real deltat
|
|
|
1180 |
integer numit
|
|
|
1181 |
integer jump
|
|
|
1182 |
real wfactor
|
|
|
1183 |
integer fbflag
|
|
|
1184 |
real spt0(nx*ny) ,spt1(nx*ny)
|
|
|
1185 |
real uut0(nx*ny*nz),uut1(nx*ny*nz)
|
|
|
1186 |
real vvt0(nx*ny*nz),vvt1(nx*ny*nz)
|
|
|
1187 |
real tht0(nx*ny*nz),tht1(nx*ny*nz),sth0(nx*ny),sth1(nx*ny)
|
|
|
1188 |
real p3d0(nx*ny*nz),p3d1(nx*ny*nz)
|
|
|
1189 |
real xmin,ymin,dx,dy
|
|
|
1190 |
integer per
|
|
|
1191 |
integer hem
|
|
|
1192 |
real mdv
|
|
|
1193 |
real mpx(nx*ny),mpy(nx*ny)
|
|
|
1194 |
real isen
|
|
|
1195 |
|
|
|
1196 |
c Auxiliary variables
|
|
|
1197 |
real xmax,ymax
|
|
|
1198 |
real xind,yind,pind
|
|
|
1199 |
|
|
|
1200 |
real u0,v0,u1,v1,u,v,w,sp
|
|
|
1201 |
integer icount
|
|
|
1202 |
character ch
|
|
|
1203 |
real mpsc_x,mpsc_y
|
|
|
1204 |
|
|
|
1205 |
c Externals
|
|
|
1206 |
real int_index4
|
|
|
1207 |
external int_index4
|
|
|
1208 |
|
|
|
1209 |
c Reset the flag for domain-leaving
|
|
|
1210 |
left=0
|
|
|
1211 |
|
|
|
1212 |
c Set the esat-north boundary of the domain
|
|
|
1213 |
xmax = xmin+real(nx-1)*dx
|
|
|
1214 |
ymax = ymin+real(ny-1)*dy
|
|
|
1215 |
|
|
|
1216 |
C Interpolate wind fields to starting position (x0,y0,p0)
|
|
|
1217 |
call get_index4 (xind,yind,pind,x0,y0,p0,reltpos0,
|
|
|
1218 |
> p3d0,p3d1,spt0,spt1,3,
|
|
|
1219 |
> nx,ny,nz,xmin,ymin,dx,dy,mdv)
|
|
|
1220 |
u0 = int_index4(uut0,uut1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
|
|
|
1221 |
v0 = int_index4(vvt0,vvt1,nx,ny,nz,xind,yind,pind,reltpos0,mdv)
|
|
|
1222 |
|
|
|
1223 |
C For first iteration take ending position equal to starting position
|
|
|
1224 |
x1 = x0
|
|
|
1225 |
y1 = y0
|
|
|
1226 |
p1 = p0
|
|
|
1227 |
|
|
|
1228 |
C Iterative calculation of new position
|
|
|
1229 |
do icount=1,numit
|
|
|
1230 |
|
|
|
1231 |
C Calculate new winds for advection
|
|
|
1232 |
call get_index4 (xind,yind,pind,x1,y1,p1,reltpos1,
|
|
|
1233 |
> p3d0,p3d1,spt0,spt1,3,
|
|
|
1234 |
> nx,ny,nz,xmin,ymin,dx,dy,mdv)
|
|
|
1235 |
u1 = int_index4(uut0,uut1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
|
|
|
1236 |
v1 = int_index4(vvt0,vvt1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
|
|
|
1237 |
|
|
|
1238 |
c Get the new velocity in between
|
|
|
1239 |
u=(u0+u1)/2.
|
|
|
1240 |
v=(v0+v1)/2.
|
|
|
1241 |
|
|
|
1242 |
c Get the mapping scale factors for this position
|
|
|
1243 |
mpsc_x = int_index4 (mpx,mpx,nx,ny,1,xind,yind,1.,reltpos1,mdv)
|
|
|
1244 |
mpsc_y = int_index4 (mpy,mpy,nx,ny,1,xind,yind,1.,reltpos1,mdv)
|
|
|
1245 |
|
|
|
1246 |
C Calculate new positions (adapted for cartesian grid)
|
|
|
1247 |
x1 = x0 + fbflag * deltat * u * dx/mpsc_x
|
|
|
1248 |
y1 = y0 + fbflag * deltat * v * dy/mpsc_y
|
|
|
1249 |
|
|
|
1250 |
c Write the update positions for tests
|
|
|
1251 |
if ( (icount.eq.3).and.(test.eq.1) ) then
|
|
|
1252 |
write(*,'(10f10.2)') x0,u,x1-x0,p0,w,p1-p0
|
|
|
1253 |
endif
|
|
|
1254 |
|
|
|
1255 |
c Set trajectory height to isentropic level
|
|
|
1256 |
call get_index4 (xind,yind,pind,x1,y1,isen,reltpos1,
|
|
|
1257 |
> tht0,tht1,sth0,sth1,1,
|
|
|
1258 |
> nx,ny,nz,xmin,ymin,dx,dy,mdv)
|
|
|
1259 |
p1 = int_index4
|
|
|
1260 |
> (p3d0,p3d1,nx,ny,nz,xind,yind,pind,reltpos1,mdv)
|
|
|
1261 |
|
|
|
1262 |
C Interpolate surface geop.height to actual position
|
|
|
1263 |
call get_index4 (xind,yind,pind,x1,y1,0.,reltpos1,
|
|
|
1264 |
> p3d0,p3d1,spt0,spt1,3,
|
|
|
1265 |
> nx,ny,nz,xmin,ymin,dx,dy,mdv)
|
|
|
1266 |
sp = int_index4 (spt0,spt1,nx,ny,1,xind,yind,1.,reltpos1,mdv)
|
|
|
1267 |
|
|
|
1268 |
c Handle trajectories which cross the lower boundary (jump flag)
|
|
|
1269 |
if ((jump.eq.1).and.(p1.lt.sp)) then
|
|
|
1270 |
print*,'jump'
|
|
|
1271 |
p1=sp+10.
|
|
|
1272 |
endif
|
|
|
1273 |
|
|
|
1274 |
c Handle periodioc boundaries for 'ideal' mode
|
|
|
1275 |
if ( per.eq.1 ) then
|
|
|
1276 |
if ( x1.gt.xmax ) x1=x1-xmax
|
|
|
1277 |
if ( x1.lt.xmin ) x1=x1+xmax
|
|
|
1278 |
endif
|
|
|
1279 |
|
|
|
1280 |
C Check if trajectory leaves data domain
|
|
|
1281 |
if ( (x1.lt.xmin).or.(x1.gt.xmax).or.
|
|
|
1282 |
> (y1.lt.ymin).or.(y1.gt.ymax).or.(p1.lt.sp) ) ! geopot : .lt. ; pressure: .gt.
|
|
|
1283 |
> then
|
|
|
1284 |
left=1
|
|
|
1285 |
goto 100
|
|
|
1286 |
endif
|
|
|
1287 |
|
|
|
1288 |
enddo
|
|
|
1289 |
|
|
|
1290 |
c Exit point for subroutine
|
|
|
1291 |
100 continue
|
|
|
1292 |
|
|
|
1293 |
return
|
|
|
1294 |
|
|
|
1295 |
end
|
|
|
1296 |
|
2 |
michaesp |
1297 |
c ----------------------------------------------------------------------
|
|
|
1298 |
c Get spherical distance between lat/lon points
|
|
|
1299 |
c ----------------------------------------------------------------------
|
|
|
1300 |
|
|
|
1301 |
real function sdis(xp,yp,xq,yq,unit)
|
|
|
1302 |
|
|
|
1303 |
c Calculates spherical distance (in km) between two points given
|
|
|
1304 |
c by their spherical coordinates (xp,yp) and (xq,yq), respectively.
|
|
|
1305 |
|
|
|
1306 |
real,parameter :: re=6371.229
|
|
|
1307 |
real,parameter :: rad2deg=180./3.14159265
|
|
|
1308 |
real,parameter :: deg2rad=3.141592654/180.
|
|
|
1309 |
|
|
|
1310 |
real xp,yp,xq,yq,arg
|
|
|
1311 |
character*80 unit
|
|
|
1312 |
real dlon,dlat,alpha,rlat1,ralt2
|
|
|
1313 |
|
|
|
1314 |
if ( unit.eq.'km' ) then
|
|
|
1315 |
|
15 |
michaesp |
1316 |
arg=sin(yp*deg2rad)*sin(yq*deg2rad)+
|
|
|
1317 |
> cos(yp*deg2rad)*cos(yq*deg2rad)*cos((xp-xq)*deg2rad)
|
2 |
michaesp |
1318 |
if (arg.lt.-1.) arg=-1.
|
|
|
1319 |
if (arg.gt.1.) arg=1.
|
|
|
1320 |
sdis=re*acos(arg)
|
|
|
1321 |
|
|
|
1322 |
elseif ( unit.eq.'deg' ) then
|
|
|
1323 |
|
|
|
1324 |
dlon = xp-xq
|
|
|
1325 |
if ( dlon.gt. 180. ) dlon = dlon - 360.
|
|
|
1326 |
if ( dlon.lt.-180. ) dlon = dlon + 360.
|
|
|
1327 |
sdis = sqrt( dlon**2 + (yp-yq)**2 )
|
|
|
1328 |
|
|
|
1329 |
elseif ( unit.eq.'km.haversine' ) then
|
|
|
1330 |
|
|
|
1331 |
dlat = (yp - yq)*deg2rad
|
|
|
1332 |
dlon = (xp - xq)*deg2rad
|
|
|
1333 |
rlat1 = yp*deg2rad
|
|
|
1334 |
rlat2 = yq*deg2rad
|
|
|
1335 |
|
|
|
1336 |
alpha = ( sin(0.5*dlat)**2 ) +
|
|
|
1337 |
> ( sin(0.5*dlon)**2 )*cos(rlat1)*cos(rlat2)
|
|
|
1338 |
|
|
|
1339 |
sdis = 4. * re * atan2( sqrt(alpha),1. + sqrt( 1. - alpha ) )
|
|
|
1340 |
|
|
|
1341 |
endif
|
|
|
1342 |
|
|
|
1343 |
end
|