11 subroutine advect_pv_tl (qnew,q,q_traj,q_north,q_south, &
13 & nx,ny,deltax,deltay,dt)
20 real(kind=kind_real),
intent(out) :: qnew(nx,ny,2)
21 real(kind=kind_real),
intent(in) :: q(nx,ny,2)
22 real(kind=kind_real),
intent(in) :: q_traj(nx,ny,2)
23 real(kind=kind_real),
intent(in) :: q_north(nx,2)
24 real(kind=kind_real),
intent(in) :: q_south(nx,2)
25 real(kind=kind_real),
intent(in) :: u(nx,ny,2)
26 real(kind=kind_real),
intent(in) :: u_traj(nx,ny,2)
27 real(kind=kind_real),
intent(in) :: v(nx,ny,2)
28 real(kind=kind_real),
intent(in) :: v_traj(nx,ny,2)
29 integer,
intent(in) :: nx
30 integer,
intent(in) :: ny
31 real(kind=kind_real),
intent(in) :: deltax
32 real(kind=kind_real),
intent(in) :: deltay
33 real(kind=kind_real),
intent(in) :: dt
35 integer :: ii,jj,kk,ixm1,ix,ixp1,ixp2,jym1,jy,jyp1,jyp2
36 real(kind=kind_real) :: ax_traj,ay_traj,qjm1_traj,qj_traj,qjp1_traj,qjp2_traj
37 real(kind=kind_real) :: ax,ay,qjm1,qj,qjp1,qjp2
38 real(kind_real),
parameter :: one=1.0_kind_real
39 real(kind_real),
parameter :: two=2.0_kind_real
40 real(kind_real),
parameter :: half=0.5_kind_real
41 real(kind_real),
parameter :: sixth=1.0_kind_real/6.0_kind_real
51 ax_traj =
real(ii,kind_real) - u_traj(ii,jj,kk)*dt/deltax
52 ax = -u(ii,jj,kk)*dt/deltax
55 ax_traj = ax_traj-
real(ix,kind_real)
56 ixm1 = 1 + modulo(ix-2,nx)
57 ixp1 = 1 + modulo(ix ,nx)
58 ixp2 = 1 + modulo(ix+1,nx)
59 ix = 1 + modulo(ix-1,nx)
61 ay_traj =
real(jj,kind_real) - v_traj(ii,jj,kk)*dt/deltay
62 ay = -v(ii,jj,kk)*dt/deltay
65 ay_traj = ay_traj-
real(jy,kind_real)
74 & ax_traj *(ax_traj-one)*(ax_traj-two) &
75 & *q_south(ixm1,kk)*(-sixth) + &
76 & (ax_traj+one)*(ax_traj-one)*(ax_traj-two) &
77 & *q_south(ix, kk)*half + &
78 & (ax_traj+one)* ax_traj *(ax_traj-two) &
79 & *q_south(ixp1,kk)*(-half) + &
80 & (ax_traj+one) * ax_traj *(ax_traj-one) &
81 & *q_south(ixp2,kk)*(sixth)
84 & ax *(ax_traj-one)*(ax_traj-two) &
85 & *q_south(ixm1,kk)*(-sixth) + &
86 & ax_traj * ax *(ax_traj-two) &
87 & *q_south(ixm1,kk)*(-sixth) + &
88 & ax_traj *(ax_traj-one)* ax &
89 & *q_south(ixm1,kk)*(-sixth) + &
90 & ax *(ax_traj-one)*(ax_traj-two) &
91 & *q_south(ix, kk)*half + &
92 & (ax_traj+one)* ax *(ax_traj-two) &
93 & *q_south(ix, kk)*half + &
94 & (ax_traj+one)*(ax_traj-one)* ax &
95 & *q_south(ix, kk)*half + &
96 & ax * ax_traj *(ax_traj-two) &
97 & *q_south(ixp1,kk)*(-half) + &
98 & (ax_traj+one)* ax *(ax_traj-two) &
99 & *q_south(ixp1,kk)*(-half) + &
100 & (ax_traj+one)* ax_traj * ax &
101 & *q_south(ixp1,kk)*(-half) + &
102 & ax * ax_traj *(ax_traj-one) &
103 & *q_south(ixp2,kk)*(sixth) + &
104 & (ax_traj+one) * ax *(ax_traj-one) &
105 & *q_south(ixp2,kk)*(sixth) + &
106 & (ax_traj+one) * ax_traj * ax &
107 & *q_south(ixp2,kk)*(sixth)
108 else if (jym1 > ny)
then 110 & ax_traj *(ax_traj-one)*(ax_traj-two) &
111 & *q_north(ixm1,kk)*(-sixth) + &
112 & (ax_traj+one)*(ax_traj-one)*(ax_traj-two) &
113 & *q_north(ix, kk)*half + &
114 & (ax_traj+one)* ax_traj *(ax_traj-two) &
115 & *q_north(ixp1,kk)*(-half) + &
116 & (ax_traj+one)* ax_traj *(ax_traj-one) &
117 & *q_north(ixp2,kk)*(sixth)
120 & ax *(ax_traj-one)*(ax_traj-two) &
121 & *q_north(ixm1,kk)*(-sixth) + &
122 & ax_traj * ax *(ax_traj-two) &
123 & *q_north(ixm1,kk)*(-sixth) + &
124 & ax_traj *(ax_traj-one)* ax &
125 & *q_north(ixm1,kk)*(-sixth) + &
126 & ax *(ax_traj-one)*(ax_traj-two) &
127 & *q_north(ix, kk)*half + &
128 & (ax_traj+one)* ax *(ax_traj-two) &
129 & *q_north(ix, kk)*half + &
130 & (ax_traj+one)*(ax_traj-one)* ax &
131 & *q_north(ix, kk)*half + &
132 & ax * ax_traj *(ax_traj-two) &
133 & *q_north(ixp1,kk)*(-half) + &
134 & (ax_traj+one)* ax *(ax_traj-two) &
135 & *q_north(ixp1,kk)*(-half) + &
136 & (ax_traj+one)* ax_traj * ax &
137 & *q_north(ixp1,kk)*(-half) + &
138 & ax * ax_traj *(ax_traj-one) &
139 & *q_north(ixp2,kk)*(sixth) + &
140 & (ax_traj+one)* ax *(ax_traj-one) &
141 & *q_north(ixp2,kk)*(sixth) + &
142 & (ax_traj+one)* ax_traj * ax &
143 & *q_north(ixp2,kk)*(sixth)
146 & ax_traj *(ax_traj-one)*(ax_traj-two) &
147 & *q_traj(ixm1,jym1,kk)*(-sixth) + &
148 & (ax_traj+one)*(ax_traj-one)*(ax_traj-two) &
149 & *q_traj(ix, jym1,kk)*half + &
150 & (ax_traj+one)* ax_traj *(ax_traj-two) &
151 & *q_traj(ixp1,jym1,kk)*(-half) + &
152 & (ax_traj+one)* ax_traj *(ax_traj-one) &
153 & *q_traj(ixp2,jym1,kk)*(sixth)
156 & ax *(ax_traj-one)*(ax_traj-two) &
157 & *q_traj(ixm1,jym1,kk)*(-sixth) + &
158 & ax_traj * ax *(ax_traj-two) &
159 & *q_traj(ixm1,jym1,kk)*(-sixth) + &
160 & ax_traj *(ax_traj-one)* ax &
161 & *q_traj(ixm1,jym1,kk)*(-sixth) + &
162 & ax_traj *(ax_traj-one)*(ax_traj-two) &
163 & *q(ixm1,jym1,kk)*(-sixth) + &
164 & ax *(ax_traj-one)*(ax_traj-two) &
165 & *q_traj(ix, jym1,kk)*half + &
166 & (ax_traj+one)* ax *(ax_traj-two) &
167 & *q_traj(ix, jym1,kk)*half + &
168 & (ax_traj+one)*(ax_traj-one)* ax &
169 & *q_traj(ix, jym1,kk)*half + &
170 & (ax_traj+one)*(ax_traj-one)*(ax_traj-two) &
171 & *q(ix, jym1,kk)*half + &
172 & ax * ax_traj *(ax_traj-two) &
173 & *q_traj(ixp1,jym1,kk)*(-half) + &
174 & (ax_traj+one)* ax *(ax_traj-two) &
175 & *q_traj(ixp1,jym1,kk)*(-half) + &
176 & (ax_traj+one)* ax_traj * ax &
177 & *q_traj(ixp1,jym1,kk)*(-half) + &
178 & (ax_traj+one)* ax_traj *(ax_traj-two) &
179 & *q(ixp1,jym1,kk)*(-half) + &
180 & ax * ax_traj *(ax_traj-one) &
181 & *q_traj(ixp2,jym1,kk)*(sixth) + &
182 & (ax_traj+one)* ax *(ax_traj-one) &
183 & *q_traj(ixp2,jym1,kk)*(sixth) + &
184 & (ax_traj+one)* ax_traj * ax &
185 & *q_traj(ixp2,jym1,kk)*(sixth) + &
186 & (ax_traj+one)* ax_traj *(ax_traj-one) &
187 & *q(ixp2,jym1,kk)*(sixth)
192 & ax_traj *(ax_traj-one)*(ax_traj-two) &
193 & *q_south(ixm1,kk)*(-sixth) + &
194 & (ax_traj+one)*(ax_traj-one)*(ax_traj-two) &
195 & *q_south(ix ,kk)*half + &
196 & (ax_traj+one)* ax_traj *(ax_traj-two) &
197 & *q_south(ixp1,kk)*(-half) + &
198 & (ax_traj+one)* ax_traj *(ax_traj-one) &
199 & *q_south(ixp2,kk)*(sixth)
202 & ax *(ax_traj-one)*(ax_traj-two) &
203 & *q_south(ixm1,kk)*(-sixth) + &
204 & ax_traj * ax *(ax_traj-two) &
205 & *q_south(ixm1,kk)*(-sixth) + &
206 & ax_traj *(ax_traj-one)* ax &
207 & *q_south(ixm1,kk)*(-sixth) + &
208 & ax *(ax_traj-one)*(ax_traj-two) &
209 & *q_south(ix ,kk)*half + &
210 & (ax_traj+one)* ax *(ax_traj-two) &
211 & *q_south(ix ,kk)*half + &
212 & (ax_traj+one)*(ax_traj-one)* ax &
213 & *q_south(ix ,kk)*half + &
214 & ax * ax_traj *(ax_traj-two) &
215 & *q_south(ixp1,kk)*(-half) + &
216 & (ax_traj+one)* ax *(ax_traj-two) &
217 & *q_south(ixp1,kk)*(-half) + &
218 & (ax_traj+one)* ax_traj * ax &
219 & *q_south(ixp1,kk)*(-half) + &
220 & ax * ax_traj *(ax_traj-one) &
221 & *q_south(ixp2,kk)*(sixth) + &
222 & (ax_traj+one)* ax *(ax_traj-one) &
223 & *q_south(ixp2,kk)*(sixth) + &
224 & (ax_traj+one)* ax_traj * ax &
225 & *q_south(ixp2,kk)*(sixth)
226 else if (jy > ny)
then 228 & ax_traj *(ax_traj-one)*(ax_traj-two) &
229 & *q_north(ixm1,kk)*(-sixth) + &
230 & (ax_traj+one)*(ax_traj-one)*(ax_traj-two) &
231 & *q_north(ix ,kk)*half + &
232 & (ax_traj+one)* ax_traj *(ax_traj-two) &
233 & *q_north(ixp1,kk)*(-half) + &
234 & (ax_traj+one)* ax_traj *(ax_traj-one) &
235 & *q_north(ixp2,kk)*(sixth)
238 & ax *(ax_traj-one)*(ax_traj-two) &
239 & *q_north(ixm1,kk)*(-sixth) + &
240 & ax_traj * ax *(ax_traj-two) &
241 & *q_north(ixm1,kk)*(-sixth) + &
242 & ax_traj *(ax_traj-one)* ax &
243 & *q_north(ixm1,kk)*(-sixth) + &
244 & ax *(ax_traj-one)*(ax_traj-two) &
245 & *q_north(ix ,kk)*half + &
246 & (ax_traj+one)* ax *(ax_traj-two) &
247 & *q_north(ix ,kk)*half + &
248 & (ax_traj+one)*(ax_traj-one)* ax &
249 & *q_north(ix ,kk)*half + &
250 & ax * ax_traj *(ax_traj-two) &
251 & *q_north(ixp1,kk)*(-half) + &
252 & (ax_traj+one)* ax *(ax_traj-two) &
253 & *q_north(ixp1,kk)*(-half) + &
254 & (ax_traj+one)* ax_traj * ax &
255 & *q_north(ixp1,kk)*(-half) + &
256 & ax * ax_traj *(ax_traj-one) &
257 & *q_north(ixp2,kk)*(sixth) + &
258 & (ax_traj+one)* ax *(ax_traj-one) &
259 & *q_north(ixp2,kk)*(sixth) + &
260 & (ax_traj+one)* ax_traj * ax &
261 & *q_north(ixp2,kk)*(sixth)
264 & ax_traj *(ax_traj-one)*(ax_traj-two) &
265 & *q_traj(ixm1,jy,kk)*(-sixth) + &
266 & (ax_traj+one)*(ax_traj-one)*(ax_traj-two) &
267 & *q_traj(ix ,jy,kk)*half + &
268 & (ax_traj+one)* ax_traj *(ax_traj-two) &
269 & *q_traj(ixp1,jy,kk)*(-half) + &
270 & (ax_traj+one)* ax_traj *(ax_traj-one) &
271 & *q_traj(ixp2,jy,kk)*(sixth)
274 & ax *(ax_traj-one)*(ax_traj-two) &
275 & *q_traj(ixm1,jy,kk)*(-sixth) + &
276 & ax_traj * ax *(ax_traj-two) &
277 & *q_traj(ixm1,jy,kk)*(-sixth) + &
278 & ax_traj *(ax_traj-one)* ax &
279 & *q_traj(ixm1,jy,kk)*(-sixth) + &
280 & ax_traj *(ax_traj-one)*(ax_traj-two) &
281 & *q(ixm1,jy,kk)*(-sixth) + &
282 & ax *(ax_traj-one)*(ax_traj-two) &
283 & *q_traj(ix ,jy,kk)*half + &
284 & (ax_traj+one)* ax *(ax_traj-two) &
285 & *q_traj(ix ,jy,kk)*half + &
286 & (ax_traj+one)*(ax_traj-one)* ax &
287 & *q_traj(ix ,jy,kk)*half + &
288 & (ax_traj+one)*(ax_traj-one)*(ax_traj-two) &
289 & *q(ix ,jy,kk)*half + &
290 & ax * ax_traj *(ax_traj-two) &
291 & *q_traj(ixp1,jy,kk)*(-half) + &
292 & (ax_traj+one)* ax *(ax_traj-two) &
293 & *q_traj(ixp1,jy,kk)*(-half) + &
294 & (ax_traj+one)* ax_traj * ax &
295 & *q_traj(ixp1,jy,kk)*(-half) + &
296 & (ax_traj+one)* ax_traj *(ax_traj-two) &
297 & *q(ixp1,jy,kk)*(-half) + &
298 & ax * ax_traj *(ax_traj-one) &
299 & *q_traj(ixp2,jy,kk)*(sixth) + &
300 & (ax_traj+one)* ax *(ax_traj-one) &
301 & *q_traj(ixp2,jy,kk)*(sixth) + &
302 & (ax_traj+one)* ax_traj * ax &
303 & *q_traj(ixp2,jy,kk)*(sixth) + &
304 & (ax_traj+one)* ax_traj *(ax_traj-one) &
305 & *q(ixp2,jy,kk)*(sixth)
310 & ax_traj *(ax_traj-one)*(ax_traj-two) &
311 & *q_south(ixm1,kk)*(-sixth) + &
312 & (ax_traj+one)*(ax_traj-one)*(ax_traj-two) &
313 & *q_south(ix ,kk)*half + &
314 & (ax_traj+one)* ax_traj *(ax_traj-two) &
315 & *q_south(ixp1,kk)*(-half) + &
316 & (ax_traj+one)* ax_traj *(ax_traj-one) &
317 & *q_south(ixp2,kk)*(sixth)
320 & ax *(ax_traj-one)*(ax_traj-two) &
321 & *q_south(ixm1,kk)*(-sixth) + &
322 & ax_traj * ax *(ax_traj-two) &
323 & *q_south(ixm1,kk)*(-sixth) + &
324 & ax_traj *(ax_traj-one)* ax &
325 & *q_south(ixm1,kk)*(-sixth) + &
326 & ax *(ax_traj-one)*(ax_traj-two) &
327 & *q_south(ix ,kk)*half + &
328 & (ax_traj+one)* ax *(ax_traj-two) &
329 & *q_south(ix ,kk)*half + &
330 & (ax_traj+one)*(ax_traj-one)* ax &
331 & *q_south(ix ,kk)*half + &
332 & ax * ax_traj *(ax_traj-two) &
333 & *q_south(ixp1,kk)*(-half) + &
334 & (ax_traj+one)* ax *(ax_traj-two) &
335 & *q_south(ixp1,kk)*(-half) + &
336 & (ax_traj+one)* ax_traj * ax &
337 & *q_south(ixp1,kk)*(-half) + &
338 & ax * ax_traj *(ax_traj-one) &
339 & *q_south(ixp2,kk)*(sixth) + &
340 & (ax_traj+one)* ax *(ax_traj-one) &
341 & *q_south(ixp2,kk)*(sixth) + &
342 & (ax_traj+one)* ax_traj * ax &
343 & *q_south(ixp2,kk)*(sixth)
344 else if (jyp1 > ny)
then 346 & ax_traj *(ax_traj-one)*(ax_traj-two) &
347 & *q_north(ixm1,kk)*(-sixth) + &
348 & (ax_traj+one)*(ax_traj-one)*(ax_traj-two) &
349 & *q_north(ix ,kk)*half + &
350 & (ax_traj+one)* ax_traj *(ax_traj-two) &
351 & *q_north(ixp1,kk)*(-half) + &
352 & (ax_traj+one)* ax_traj *(ax_traj-one) &
353 & *q_north(ixp2,kk)*(sixth)
356 & ax *(ax_traj-one)*(ax_traj-two) &
357 & *q_north(ixm1,kk)*(-sixth) + &
358 & ax_traj * ax *(ax_traj-two) &
359 & *q_north(ixm1,kk)*(-sixth) + &
360 & ax_traj *(ax_traj-one)* ax &
361 & *q_north(ixm1,kk)*(-sixth) + &
362 & ax *(ax_traj-one)*(ax_traj-two) &
363 & *q_north(ix ,kk)*half + &
364 & (ax_traj+one)* ax *(ax_traj-two) &
365 & *q_north(ix ,kk)*half + &
366 & (ax_traj+one)*(ax_traj-one)* ax &
367 & *q_north(ix ,kk)*half + &
368 & ax * ax_traj *(ax_traj-two) &
369 & *q_north(ixp1,kk)*(-half) + &
370 & (ax_traj+one)* ax *(ax_traj-two) &
371 & *q_north(ixp1,kk)*(-half) + &
372 & (ax_traj+one)* ax_traj * ax &
373 & *q_north(ixp1,kk)*(-half) + &
374 & ax * ax_traj *(ax_traj-one) &
375 & *q_north(ixp2,kk)*(sixth) + &
376 & (ax_traj+one)* ax *(ax_traj-one) &
377 & *q_north(ixp2,kk)*(sixth) + &
378 & (ax_traj+one)* ax_traj * ax &
379 & *q_north(ixp2,kk)*(sixth)
382 & ax_traj *(ax_traj-one)*(ax_traj-two) &
383 & *q_traj(ixm1,jyp1,kk)*(-sixth) + &
384 & (ax_traj+one)*(ax_traj-one)*(ax_traj-two) &
385 & *q_traj(ix ,jyp1,kk)*half + &
386 & (ax_traj+one)* ax_traj *(ax_traj-two) &
387 & *q_traj(ixp1,jyp1,kk)*(-half) + &
388 & (ax_traj+one)* ax_traj *(ax_traj-one) &
389 & *q_traj(ixp2,jyp1,kk)*(sixth)
392 & ax *(ax_traj-one)*(ax_traj-two) &
393 & *q_traj(ixm1,jyp1,kk)*(-sixth) + &
394 & ax_traj * ax *(ax_traj-two) &
395 & *q_traj(ixm1,jyp1,kk)*(-sixth) + &
396 & ax_traj *(ax_traj-one)* ax &
397 & *q_traj(ixm1,jyp1,kk)*(-sixth) + &
398 & ax_traj *(ax_traj-one)*(ax_traj-two) &
399 & *q(ixm1,jyp1,kk)*(-sixth) + &
400 & ax *(ax_traj-one)*(ax_traj-two) &
401 & *q_traj(ix ,jyp1,kk)*half + &
402 & (ax_traj+one)* ax *(ax_traj-two) &
403 & *q_traj(ix ,jyp1,kk)*half + &
404 & (ax_traj+one)*(ax_traj-one)* ax &
405 & *q_traj(ix ,jyp1,kk)*half + &
406 & (ax_traj+one)*(ax_traj-one)*(ax_traj-two) &
407 & *q(ix ,jyp1,kk)*half + &
408 & ax * ax_traj *(ax_traj-two) &
409 & *q_traj(ixp1,jyp1,kk)*(-half) + &
410 & (ax_traj+one)* ax *(ax_traj-two) &
411 & *q_traj(ixp1,jyp1,kk)*(-half) + &
412 & (ax_traj+one)* ax_traj * ax &
413 & *q_traj(ixp1,jyp1,kk)*(-half) + &
414 & (ax_traj+one)* ax_traj *(ax_traj-two) &
415 & *q(ixp1,jyp1,kk)*(-half) + &
416 & ax * ax_traj *(ax_traj-one) &
417 & *q_traj(ixp2,jyp1,kk)*(sixth) + &
418 & (ax_traj+one)* ax *(ax_traj-one) &
419 & *q_traj(ixp2,jyp1,kk)*(sixth) + &
420 & (ax_traj+one)* ax_traj * ax &
421 & *q_traj(ixp2,jyp1,kk)*(sixth) + &
422 & (ax_traj+one)* ax_traj *(ax_traj-one) &
423 & *q(ixp2,jyp1,kk)*(sixth)
428 & ax_traj *(ax_traj-one)*(ax_traj-two) &
429 & *q_south(ixm1,kk)*(-sixth) + &
430 & (ax_traj+one)*(ax_traj-one)*(ax_traj-two) &
431 & *q_south(ix ,kk)*half + &
432 & (ax_traj+one)* ax_traj *(ax_traj-two) &
433 & *q_south(ixp1,kk)*(-half) + &
434 & (ax_traj+one)* ax_traj *(ax_traj-one) &
435 & *q_south(ixp2,kk)*(sixth)
438 & ax *(ax_traj-one)*(ax_traj-two) &
439 & *q_south(ixm1,kk)*(-sixth) + &
440 & ax_traj * ax *(ax_traj-two) &
441 & *q_south(ixm1,kk)*(-sixth) + &
442 & ax_traj *(ax_traj-one)* ax &
443 & *q_south(ixm1,kk)*(-sixth) + &
444 & ax *(ax_traj-one)*(ax_traj-two) &
445 & *q_south(ix ,kk)*half + &
446 & (ax_traj+one)* ax *(ax_traj-two) &
447 & *q_south(ix ,kk)*half + &
448 & (ax_traj+one)*(ax_traj-one)* ax &
449 & *q_south(ix ,kk)*half + &
450 & ax * ax_traj *(ax_traj-two) &
451 & *q_south(ixp1,kk)*(-half) + &
452 & (ax_traj+one)* ax *(ax_traj-two) &
453 & *q_south(ixp1,kk)*(-half) + &
454 & (ax_traj+one)* ax_traj * ax &
455 & *q_south(ixp1,kk)*(-half) + &
456 & ax * ax_traj *(ax_traj-one) &
457 & *q_south(ixp2,kk)*(sixth) + &
458 & (ax_traj+one)* ax *(ax_traj-one) &
459 & *q_south(ixp2,kk)*(sixth) + &
460 & (ax_traj+one)* ax_traj * ax &
461 & *q_south(ixp2,kk)*(sixth)
462 else if (jyp2 > ny)
then 464 & ax_traj *(ax_traj-one)*(ax_traj-two) &
465 & *q_north(ixm1,kk)*(-sixth) + &
466 & (ax_traj+one)*(ax_traj-one)*(ax_traj-two) &
467 & *q_north(ix ,kk)*half + &
468 & (ax_traj+one)* ax_traj *(ax_traj-two) &
469 & *q_north(ixp1,kk)*(-half) + &
470 & (ax_traj+one)* ax_traj *(ax_traj-one) &
471 & *q_north(ixp2,kk)*(sixth)
474 & ax *(ax_traj-one)*(ax_traj-two) &
475 & *q_north(ixm1,kk)*(-sixth) + &
476 & ax_traj * ax *(ax_traj-two) &
477 & *q_north(ixm1,kk)*(-sixth) + &
478 & ax_traj *(ax_traj-one)* ax &
479 & *q_north(ixm1,kk)*(-sixth) + &
480 & ax *(ax_traj-one)*(ax_traj-two) &
481 & *q_north(ix ,kk)*half + &
482 & (ax_traj+one)* ax *(ax_traj-two) &
483 & *q_north(ix ,kk)*half + &
484 & (ax_traj+one)*(ax_traj-one)* ax &
485 & *q_north(ix ,kk)*half + &
486 & ax * ax_traj *(ax_traj-two) &
487 & *q_north(ixp1,kk)*(-half) + &
488 & (ax_traj+one)* ax *(ax_traj-two) &
489 & *q_north(ixp1,kk)*(-half) + &
490 & (ax_traj+one)* ax_traj * ax &
491 & *q_north(ixp1,kk)*(-half) + &
492 & ax * ax_traj *(ax_traj-one) &
493 & *q_north(ixp2,kk)*(sixth) + &
494 & (ax_traj+one)* ax *(ax_traj-one) &
495 & *q_north(ixp2,kk)*(sixth) + &
496 & (ax_traj+one)* ax_traj * ax &
497 & *q_north(ixp2,kk)*(sixth)
500 & ax_traj *(ax_traj-one)*(ax_traj-two) &
501 & *q_traj(ixm1,jyp2,kk)*(-sixth) + &
502 & (ax_traj+one)*(ax_traj-one)*(ax_traj-two) &
503 & *q_traj(ix ,jyp2,kk)*half + &
504 & (ax_traj+one)* ax_traj *(ax_traj-two) &
505 & *q_traj(ixp1,jyp2,kk)*(-half) + &
506 & (ax_traj+one)* ax_traj *(ax_traj-one) &
507 &*q_traj(ixp2,jyp2,kk)*(sixth)
510 & ax *(ax_traj-one)*(ax_traj-two) &
511 & *q_traj(ixm1,jyp2,kk)*(-sixth) + &
512 & ax_traj * ax *(ax_traj-two) &
513 & *q_traj(ixm1,jyp2,kk)*(-sixth) + &
514 & ax_traj *(ax_traj-one)* ax &
515 & *q_traj(ixm1,jyp2,kk)*(-sixth) + &
516 & ax_traj *(ax_traj-one)*(ax_traj-two) &
517 & *q(ixm1,jyp2,kk)*(-sixth) + &
518 & ax *(ax_traj-one)*(ax_traj-two) &
519 & *q_traj(ix ,jyp2,kk)*half + &
520 & (ax_traj+one)* ax *(ax_traj-two) &
521 & *q_traj(ix ,jyp2,kk)*half + &
522 & (ax_traj+one)*(ax_traj-one)* ax &
523 & *q_traj(ix ,jyp2,kk)*half + &
524 & (ax_traj+one)*(ax_traj-one)*(ax_traj-two) &
525 & *q(ix ,jyp2,kk)*half + &
526 & ax * ax_traj *(ax_traj-two) &
527 & *q_traj(ixp1,jyp2,kk)*(-half) + &
528 & (ax_traj+one)* ax *(ax_traj-two) &
529 & *q_traj(ixp1,jyp2,kk)*(-half) + &
530 & (ax_traj+one)* ax_traj * ax &
531 & *q_traj(ixp1,jyp2,kk)*(-half) + &
532 & (ax_traj+one)* ax_traj *(ax_traj-two) &
533 & *q(ixp1,jyp2,kk)*(-half) + &
534 & ax * ax_traj *(ax_traj-one) &
535 &*q_traj(ixp2,jyp2,kk)*(sixth) + &
536 & (ax_traj+one)* ax *(ax_traj-one) &
537 &*q_traj(ixp2,jyp2,kk)*(sixth) + &
538 & (ax_traj+one)* ax_traj * ax &
539 &*q_traj(ixp2,jyp2,kk)*(sixth) + &
540 & (ax_traj+one)* ax_traj *(ax_traj-one) &
541 &*q(ixp2,jyp2,kk)*(sixth)
553 & ay *(ay_traj-one)*(ay_traj-two)*(-sixth)*qjm1_traj + &
554 & ay_traj * ay *(ay_traj-two)*(-sixth)*qjm1_traj + &
555 & ay_traj *(ay_traj-one)* ay *(-sixth)*qjm1_traj + &
556 & ay_traj *(ay_traj-one)*(ay_traj-two)*(-sixth)*qjm1 + &
557 & ay *(ay_traj-one)*(ay_traj-two)*half *qj_traj + &
558 & (ay_traj+one)* ay *(ay_traj-two)*half *qj_traj + &
559 & (ay_traj+one)*(ay_traj-one)* ay *half *qj_traj + &
560 & (ay_traj+one)*(ay_traj-one)*(ay_traj-two)*half *qj + &
561 & ay * ay_traj *(ay_traj-two)*(-half) *qjp1_traj + &
562 & (ay_traj+one)* ay *(ay_traj-two)*(-half) *qjp1_traj + &
563 & (ay_traj+one)* ay_traj * ay *(-half) *qjp1_traj + &
564 & (ay_traj+one)* ay_traj *(ay_traj-two)*(-half) *qjp1 + &
565 & ay * ay_traj *(ay_traj-one)*(sixth) *qjp2_traj + &
566 & (ay_traj+one)* ay *(ay_traj-one)*(sixth) *qjp2_traj + &
567 & (ay_traj+one)* ay_traj * ay *(sixth) *qjp2_traj + &
568 & (ay_traj+one)* ay_traj *(ay_traj-one)*(sixth) *qjp2
subroutine advect_pv_tl(qnew, q, q_traj, q_north, q_south, u, u_traj, v, v_traj, nx, ny, deltax, deltay, dt)
Advect potential vorticity - Tangent Linear.