1 |
! |
2 |
! $Header: /home/cvsroot/LMDZ4/libf/dyn3d/ppm3d.F,v 1.1.1.1 2004/05/19 12:53:07 lmdzadmin Exp $ |
3 |
! |
4 |
|
5 |
cFrom lin@explorer.gsfc.nasa.gov Wed Apr 15 17:44:44 1998 |
6 |
cDate: Wed, 15 Apr 1998 11:37:03 -0400 |
7 |
cFrom: lin@explorer.gsfc.nasa.gov |
8 |
cTo: Frederic.Hourdin@lmd.jussieu.fr |
9 |
cSubject: 3D transport module of the GSFC CTM and GEOS GCM |
10 |
|
11 |
|
12 |
cThis code is sent to you by S-J Lin, DAO, NASA-GSFC |
13 |
|
14 |
cNote: this version is intended for machines like CRAY |
15 |
C-90. No multitasking directives implemented. |
16 |
|
17 |
|
18 |
C ******************************************************************** |
19 |
C |
20 |
C TransPort Core for Goddard Chemistry Transport Model (G-CTM), Goddard |
21 |
C Earth Observing System General Circulation Model (GEOS-GCM), and Data |
22 |
C Assimilation System (GEOS-DAS). |
23 |
C |
24 |
C ******************************************************************** |
25 |
C |
26 |
C Purpose: given horizontal winds on a hybrid sigma-p surfaces, |
27 |
C one call to tpcore updates the 3-D mixing ratio |
28 |
C fields one time step (NDT). [vertical mass flux is computed |
29 |
C internally consistent with the discretized hydrostatic mass |
30 |
C continuity equation of the C-Grid GEOS-GCM (for IGD=1)]. |
31 |
C |
32 |
C Schemes: Multi-dimensional Flux Form Semi-Lagrangian (FFSL) scheme based |
33 |
C on the van Leer or PPM. |
34 |
C (see Lin and Rood 1996). |
35 |
C Version 4.5 |
36 |
C Last modified: Dec. 5, 1996 |
37 |
C Major changes from version 4.0: a more general vertical hybrid sigma- |
38 |
C pressure coordinate. |
39 |
C Subroutines modified: xtp, ytp, fzppm, qckxyz |
40 |
C Subroutines deleted: vanz |
41 |
C |
42 |
C Author: Shian-Jiann Lin |
43 |
C mail address: |
44 |
C Shian-Jiann Lin* |
45 |
C Code 910.3, NASA/GSFC, Greenbelt, MD 20771 |
46 |
C Phone: 301-286-9540 |
47 |
C E-mail: lin@dao.gsfc.nasa.gov |
48 |
C |
49 |
C *affiliation: |
50 |
C Joint Center for Earth Systems Technology |
51 |
C The University of Maryland Baltimore County |
52 |
C NASA - Goddard Space Flight Center |
53 |
C References: |
54 |
C |
55 |
C 1. Lin, S.-J., and R. B. Rood, 1996: Multidimensional flux form semi- |
56 |
C Lagrangian transport schemes. Mon. Wea. Rev., 124, 2046-2070. |
57 |
C |
58 |
C 2. Lin, S.-J., W. C. Chao, Y. C. Sud, and G. K. Walker, 1994: A class of |
59 |
C the van Leer-type transport schemes and its applications to the moist- |
60 |
C ure transport in a General Circulation Model. Mon. Wea. Rev., 122, |
61 |
C 1575-1593. |
62 |
C |
63 |
C ****6***0*********0*********0*********0*********0*********0**********72 |
64 |
C |
65 |
subroutine ppm3d(IGD,Q,PS1,PS2,U,V,W,NDT,IORD,JORD,KORD,NC,IMR, |
66 |
& JNP,j1,NLAY,AP,BP,PT,AE,fill,dum,Umax) |
67 |
|
68 |
c implicit none |
69 |
|
70 |
c rajout de déclarations |
71 |
c integer Jmax,kmax,ndt0,nstep,k,j,i,ic,l,js,jn,imh,iad,jad,krd |
72 |
c integer iu,iiu,j2,jmr,js0,jt |
73 |
c real dtdy,dtdy5,rcap,iml,jn0,imjm,pi,dl,dp |
74 |
c real dt,cr1,maxdt,ztc,d5,sum1,sum2,ru |
75 |
C |
76 |
C ******************************************************************** |
77 |
C |
78 |
C ============= |
79 |
C INPUT: |
80 |
C ============= |
81 |
C |
82 |
C Q(IMR,JNP,NLAY,NC): mixing ratios at current time (t) |
83 |
C NC: total number of constituents |
84 |
C IMR: first dimension (E-W); number of Grid intervals in E-W is IMR |
85 |
C JNP: 2nd dimension (N-S); number of Grid intervals in N-S is JNP-1 |
86 |
C NLAY: 3rd dimension (number of layers); vertical index increases from 1 at |
87 |
C the model top to NLAY near the surface (see fig. below). |
88 |
C It is assumed that 6 <= NLAY <= JNP (for dynamic memory allocation) |
89 |
C |
90 |
C PS1(IMR,JNP): surface pressure at current time (t) |
91 |
C PS2(IMR,JNP): surface pressure at mid-time-level (t+NDT/2) |
92 |
C PS2 is replaced by the predicted PS (at t+NDT) on output. |
93 |
C Note: surface pressure can have any unit or can be multiplied by any |
94 |
C const. |
95 |
C |
96 |
C The pressure at layer edges are defined as follows: |
97 |
C |
98 |
C p(i,j,k) = AP(k)*PT + BP(k)*PS(i,j) (1) |
99 |
C |
100 |
C Where PT is a constant having the same unit as PS. |
101 |
C AP and BP are unitless constants given at layer edges |
102 |
C defining the vertical coordinate. |
103 |
C BP(1) = 0., BP(NLAY+1) = 1. |
104 |
C The pressure at the model top is PTOP = AP(1)*PT |
105 |
C |
106 |
C For pure sigma system set AP(k) = 1 for all k, PT = PTOP, |
107 |
C BP(k) = sige(k) (sigma at edges), PS = Psfc - PTOP. |
108 |
C |
109 |
C Note: the sigma-P coordinate is a subset of Eq. 1, which in turn |
110 |
C is a subset of the following even more general sigma-P-thelta coord. |
111 |
C currently under development. |
112 |
C p(i,j,k) = (AP(k)*PT + BP(k)*PS(i,j))/(D(k)-C(k)*TE**(-1/kapa)) |
113 |
C |
114 |
C ///////////////////////////////// |
115 |
C / \ ------------- PTOP -------------- AP(1), BP(1) |
116 |
C | |
117 |
C delp(1) | ........... Q(i,j,1) ............ |
118 |
C | |
119 |
C W(1) \ / --------------------------------- AP(2), BP(2) |
120 |
C |
121 |
C |
122 |
C |
123 |
C W(k-1) / \ --------------------------------- AP(k), BP(k) |
124 |
C | |
125 |
C delp(K) | ........... Q(i,j,k) ............ |
126 |
C | |
127 |
C W(k) \ / --------------------------------- AP(k+1), BP(k+1) |
128 |
C |
129 |
C |
130 |
C |
131 |
C / \ --------------------------------- AP(NLAY), BP(NLAY) |
132 |
C | |
133 |
C delp(NLAY) | ........... Q(i,j,NLAY) ......... |
134 |
C | |
135 |
C W(NLAY)=0 \ / ------------- surface ----------- AP(NLAY+1), BP(NLAY+1) |
136 |
C ////////////////////////////////// |
137 |
C |
138 |
C U(IMR,JNP,NLAY) & V(IMR,JNP,NLAY):winds (m/s) at mid-time-level (t+NDT/2) |
139 |
C U and V may need to be polar filtered in advance in some cases. |
140 |
C |
141 |
C IGD: grid type on which winds are defined. |
142 |
C IGD = 0: A-Grid [all variables defined at the same point from south |
143 |
C pole (j=1) to north pole (j=JNP) ] |
144 |
C |
145 |
C IGD = 1 GEOS-GCM C-Grid |
146 |
C [North] |
147 |
C |
148 |
C V(i,j) |
149 |
C | |
150 |
C | |
151 |
C | |
152 |
C U(i-1,j)---Q(i,j)---U(i,j) [EAST] |
153 |
C | |
154 |
C | |
155 |
C | |
156 |
C V(i,j-1) |
157 |
C |
158 |
C U(i, 1) is defined at South Pole. |
159 |
C V(i, 1) is half grid north of the South Pole. |
160 |
C V(i,JMR) is half grid south of the North Pole. |
161 |
C |
162 |
C V must be defined at j=1 and j=JMR if IGD=1 |
163 |
C V at JNP need not be given. |
164 |
C |
165 |
C NDT: time step in seconds (need not be constant during the course of |
166 |
C the integration). Suggested value: 30 min. for 4x5, 15 min. for 2x2.5 |
167 |
C (Lat-Lon) resolution. Smaller values are recommanded if the model |
168 |
C has a well-resolved stratosphere. |
169 |
C |
170 |
C J1 defines the size of the polar cap: |
171 |
C South polar cap edge is located at -90 + (j1-1.5)*180/(JNP-1) deg. |
172 |
C North polar cap edge is located at 90 - (j1-1.5)*180/(JNP-1) deg. |
173 |
C There are currently only two choices (j1=2 or 3). |
174 |
C IMR must be an even integer if j1 = 2. Recommended value: J1=3. |
175 |
C |
176 |
C IORD, JORD, and KORD are integers controlling various options in E-W, N-S, |
177 |
C and vertical transport, respectively. Recommended values for positive |
178 |
C definite scalars: IORD=JORD=3, KORD=5. Use KORD=3 for non- |
179 |
C positive definite scalars or when linear correlation between constituents |
180 |
C is to be maintained. |
181 |
C |
182 |
C _ORD= |
183 |
C 1: 1st order upstream scheme (too diffusive, not a useful option; it |
184 |
C can be used for debugging purposes; this is THE only known "linear" |
185 |
C monotonic advection scheme.). |
186 |
C 2: 2nd order van Leer (full monotonicity constraint; |
187 |
C see Lin et al 1994, MWR) |
188 |
C 3: monotonic PPM* (slightly improved PPM of Collela & Woodward 1984) |
189 |
C 4: semi-monotonic PPM (same as 3, but overshoots are allowed) |
190 |
C 5: positive-definite PPM (constraint on the subgrid distribution is |
191 |
C only strong enough to prevent generation of negative values; |
192 |
C both overshoots & undershoots are possible). |
193 |
C 6: un-constrained PPM (nearly diffusion free; slightly faster but |
194 |
C positivity not quaranteed. Use this option only when the fields |
195 |
C and winds are very smooth). |
196 |
C |
197 |
C *PPM: Piece-wise Parabolic Method |
198 |
C |
199 |
C Note that KORD <=2 options are no longer supported. DO not use option 4 or 5. |
200 |
C for non-positive definite scalars (such as Ertel Potential Vorticity). |
201 |
C |
202 |
C The implicit numerical diffusion decreases as _ORD increases. |
203 |
C The last two options (ORDER=5, 6) should only be used when there is |
204 |
C significant explicit diffusion (such as a turbulence parameterization). You |
205 |
C might get dispersive results otherwise. |
206 |
C No filter of any kind is applied to the constituent fields here. |
207 |
C |
208 |
C AE: Radius of the sphere (meters). |
209 |
C Recommended value for the planet earth: 6.371E6 |
210 |
C |
211 |
C fill(logical): flag to do filling for negatives (see note below). |
212 |
C |
213 |
C Umax: Estimate (upper limit) of the maximum U-wind speed (m/s). |
214 |
C (220 m/s is a good value for troposphere model; 280 m/s otherwise) |
215 |
C |
216 |
C ============= |
217 |
C Output |
218 |
C ============= |
219 |
C |
220 |
C Q: mixing ratios at future time (t+NDT) (original values are over-written) |
221 |
C W(NLAY): large-scale vertical mass flux as diagnosed from the hydrostatic |
222 |
C relationship. W will have the same unit as PS1 and PS2 (eg, mb). |
223 |
C W must be divided by NDT to get the correct mass-flux unit. |
224 |
C The vertical Courant number C = W/delp_UPWIND, where delp_UPWIND |
225 |
C is the pressure thickness in the "upwind" direction. For example, |
226 |
C C(k) = W(k)/delp(k) if W(k) > 0; |
227 |
C C(k) = W(k)/delp(k+1) if W(k) < 0. |
228 |
C ( W > 0 is downward, ie, toward surface) |
229 |
C PS2: predicted PS at t+NDT (original values are over-written) |
230 |
C |
231 |
C ******************************************************************** |
232 |
C NOTES: |
233 |
C This forward-in-time upstream-biased transport scheme reduces to |
234 |
C the 2nd order center-in-time center-in-space mass continuity eqn. |
235 |
C if Q = 1 (constant fields will remain constant). This also ensures |
236 |
C that the computed vertical velocity to be identical to GEOS-1 GCM |
237 |
C for on-line transport. |
238 |
C |
239 |
C A larger polar cap is used if j1=3 (recommended for C-Grid winds or when |
240 |
C winds are noisy near poles). |
241 |
C |
242 |
C Flux-Form Semi-Lagrangian transport in the East-West direction is used |
243 |
C when and where Courant number is greater than one. |
244 |
C |
245 |
C The user needs to change the parameter Jmax or Kmax if the resolution |
246 |
C is greater than 0.5 deg in N-S or 150 layers in the vertical direction. |
247 |
C (this TransPort Core is otherwise resolution independent and can be used |
248 |
C as a library routine). |
249 |
C |
250 |
C PPM is 4th order accurate when grid spacing is uniform (x & y); 3rd |
251 |
C order accurate for non-uniform grid (vertical sigma coord.). |
252 |
C |
253 |
C Time step is limitted only by transport in the meridional direction. |
254 |
C (the FFSL scheme is not implemented in the meridional direction). |
255 |
C |
256 |
C Since only 1-D limiters are applied, negative values could |
257 |
C potentially be generated when large time step is used and when the |
258 |
C initial fields contain discontinuities. |
259 |
C This does not necessarily imply the integration is unstable. |
260 |
C These negatives are typically very small. A filling algorithm is |
261 |
C activated if the user set "fill" to be true. |
262 |
C |
263 |
C The van Leer scheme used here is nearly as accurate as the original PPM |
264 |
C due to the use of a 4th order accurate reference slope. The PPM imple- |
265 |
C mented here is an improvement over the original and is also based on |
266 |
C the 4th order reference slope. |
267 |
C |
268 |
C ****6***0*********0*********0*********0*********0*********0**********72 |
269 |
C |
270 |
C User modifiable parameters |
271 |
C |
272 |
parameter (Jmax = 361, kmax = 150) |
273 |
C |
274 |
C ****6***0*********0*********0*********0*********0*********0**********72 |
275 |
C |
276 |
C Input-Output arrays |
277 |
C |
278 |
|
279 |
real Q(IMR,JNP,NLAY,NC),PS1(IMR,JNP),PS2(IMR,JNP), |
280 |
& U(IMR,JNP,NLAY),V(IMR,JNP,NLAY),AP(NLAY+1), |
281 |
& BP(NLAY+1),W(IMR,JNP,NLAY),NDT,val(NLAY),Umax |
282 |
integer IGD,IORD,JORD,KORD,NC,IMR,JNP,j1,NLAY,AE |
283 |
integer IMRD2 |
284 |
real PT |
285 |
logical cross, fill, dum |
286 |
C |
287 |
C Local dynamic arrays |
288 |
C |
289 |
real CRX(IMR,JNP),CRY(IMR,JNP),xmass(IMR,JNP),ymass(IMR,JNP), |
290 |
& fx1(IMR+1),DPI(IMR,JNP,NLAY),delp1(IMR,JNP,NLAY), |
291 |
& WK1(IMR,JNP,NLAY),PU(IMR,JNP),PV(IMR,JNP),DC2(IMR,JNP), |
292 |
& delp2(IMR,JNP,NLAY),DQ(IMR,JNP,NLAY,NC),VA(IMR,JNP), |
293 |
& UA(IMR,JNP),qtmp(-IMR:2*IMR) |
294 |
C |
295 |
C Local static arrays |
296 |
C |
297 |
real DTDX(Jmax), DTDX5(Jmax), acosp(Jmax), |
298 |
& cosp(Jmax), cose(Jmax), DAP(kmax),DBK(Kmax) |
299 |
data NDT0, NSTEP /0, 0/ |
300 |
data cross /.true./ |
301 |
SAVE DTDY, DTDY5, RCAP, JS0, JN0, IML, |
302 |
& DTDX, DTDX5, ACOSP, COSP, COSE, DAP,DBK |
303 |
C |
304 |
|
305 |
JMR = JNP -1 |
306 |
IMJM = IMR*JNP |
307 |
j2 = JNP - j1 + 1 |
308 |
NSTEP = NSTEP + 1 |
309 |
C |
310 |
C *********** Initialization ********************** |
311 |
if(NSTEP.eq.1) then |
312 |
c |
313 |
write(6,*) '------------------------------------ ' |
314 |
write(6,*) 'NASA/GSFC Transport Core Version 4.5' |
315 |
write(6,*) '------------------------------------ ' |
316 |
c |
317 |
WRITE(6,*) 'IMR=',IMR,' JNP=',JNP,' NLAY=',NLAY,' j1=',j1 |
318 |
WRITE(6,*) 'NC=',NC,IORD,JORD,KORD,NDT |
319 |
C |
320 |
C controles sur les parametres |
321 |
if(NLAY.LT.6) then |
322 |
write(6,*) 'NLAY must be >= 6' |
323 |
stop |
324 |
endif |
325 |
if (JNP.LT.NLAY) then |
326 |
write(6,*) 'JNP must be >= NLAY' |
327 |
stop |
328 |
endif |
329 |
IMRD2=mod(IMR,2) |
330 |
if (j1.eq.2.and.IMRD2.NE.0) then |
331 |
write(6,*) 'if j1=2 IMR must be an even integer' |
332 |
stop |
333 |
endif |
334 |
|
335 |
C |
336 |
if(Jmax.lt.JNP .or. Kmax.lt.NLAY) then |
337 |
write(6,*) 'Jmax or Kmax is too small' |
338 |
stop |
339 |
endif |
340 |
C |
341 |
DO k=1,NLAY |
342 |
DAP(k) = (AP(k+1) - AP(k))*PT |
343 |
DBK(k) = BP(k+1) - BP(k) |
344 |
ENDDO |
345 |
C |
346 |
PI = 4. * ATAN(1.) |
347 |
DL = 2.*PI / float(IMR) |
348 |
DP = PI / float(JMR) |
349 |
C |
350 |
if(IGD.eq.0) then |
351 |
C Compute analytic cosine at cell edges |
352 |
call cosa(cosp,cose,JNP,PI,DP) |
353 |
else |
354 |
C Define cosine consistent with GEOS-GCM (using dycore2.0 or later) |
355 |
call cosc(cosp,cose,JNP,PI,DP) |
356 |
endif |
357 |
C |
358 |
do 15 J=2,JMR |
359 |
15 acosp(j) = 1. / cosp(j) |
360 |
C |
361 |
C Inverse of the Scaled polar cap area. |
362 |
C |
363 |
RCAP = DP / (IMR*(1.-COS((j1-1.5)*DP))) |
364 |
acosp(1) = RCAP |
365 |
acosp(JNP) = RCAP |
366 |
endif |
367 |
C |
368 |
if(NDT0 .ne. NDT) then |
369 |
DT = NDT |
370 |
NDT0 = NDT |
371 |
|
372 |
if(Umax .lt. 180.) then |
373 |
write(6,*) 'Umax may be too small!' |
374 |
endif |
375 |
CR1 = abs(Umax*DT)/(DL*AE) |
376 |
MaxDT = DP*AE / abs(Umax) + 0.5 |
377 |
write(6,*)'Largest time step for max(V)=',Umax,' is ',MaxDT |
378 |
if(MaxDT .lt. abs(NDT)) then |
379 |
write(6,*) 'Warning!!! NDT maybe too large!' |
380 |
endif |
381 |
C |
382 |
if(CR1.ge.0.95) then |
383 |
JS0 = 0 |
384 |
JN0 = 0 |
385 |
IML = IMR-2 |
386 |
ZTC = 0. |
387 |
else |
388 |
ZTC = acos(CR1) * (180./PI) |
389 |
C |
390 |
JS0 = float(JMR)*(90.-ZTC)/180. + 2 |
391 |
JS0 = max(JS0, J1+1) |
392 |
IML = min(6*JS0/(J1-1)+2, 4*IMR/5) |
393 |
JN0 = JNP-JS0+1 |
394 |
endif |
395 |
C |
396 |
C |
397 |
do J=2,JMR |
398 |
DTDX(j) = DT / ( DL*AE*COSP(J) ) |
399 |
|
400 |
DTDX5(j) = 0.5*DTDX(j) |
401 |
enddo |
402 |
C |
403 |
|
404 |
DTDY = DT /(AE*DP) |
405 |
DTDY5 = 0.5*DTDY |
406 |
C |
407 |
endif |
408 |
C |
409 |
C *********** End Initialization ********************** |
410 |
C |
411 |
C delp = pressure thickness: the psudo-density in a hydrostatic system. |
412 |
do k=1,NLAY |
413 |
do j=1,JNP |
414 |
do i=1,IMR |
415 |
delp1(i,j,k)=DAP(k)+DBK(k)*PS1(i,j) |
416 |
delp2(i,j,k)=DAP(k)+DBK(k)*PS2(i,j) |
417 |
enddo |
418 |
enddo |
419 |
enddo |
420 |
|
421 |
C |
422 |
if(j1.ne.2) then |
423 |
DO 40 IC=1,NC |
424 |
DO 40 L=1,NLAY |
425 |
DO 40 I=1,IMR |
426 |
Q(I, 2,L,IC) = Q(I, 1,L,IC) |
427 |
40 Q(I,JMR,L,IC) = Q(I,JNP,L,IC) |
428 |
endif |
429 |
C |
430 |
C Compute "tracer density" |
431 |
DO 550 IC=1,NC |
432 |
DO 44 k=1,NLAY |
433 |
DO 44 j=1,JNP |
434 |
DO 44 i=1,IMR |
435 |
44 DQ(i,j,k,IC) = Q(i,j,k,IC)*delp1(i,j,k) |
436 |
550 continue |
437 |
C |
438 |
do 1500 k=1,NLAY |
439 |
C |
440 |
if(IGD.eq.0) then |
441 |
C Convert winds on A-Grid to Courant number on C-Grid. |
442 |
call A2C(U(1,1,k),V(1,1,k),IMR,JMR,j1,j2,CRX,CRY,dtdx5,DTDY5) |
443 |
else |
444 |
C Convert winds on C-grid to Courant number |
445 |
do 45 j=j1,j2 |
446 |
do 45 i=2,IMR |
447 |
45 CRX(i,J) = dtdx(j)*U(i-1,j,k) |
448 |
|
449 |
C |
450 |
do 50 j=j1,j2 |
451 |
50 CRX(1,J) = dtdx(j)*U(IMR,j,k) |
452 |
C |
453 |
do 55 i=1,IMR*JMR |
454 |
55 CRY(i,2) = DTDY*V(i,1,k) |
455 |
endif |
456 |
C |
457 |
C Determine JS and JN |
458 |
JS = j1 |
459 |
JN = j2 |
460 |
C |
461 |
do j=JS0,j1+1,-1 |
462 |
do i=1,IMR |
463 |
if(abs(CRX(i,j)).GT.1.) then |
464 |
JS = j |
465 |
go to 2222 |
466 |
endif |
467 |
enddo |
468 |
enddo |
469 |
C |
470 |
2222 continue |
471 |
do j=JN0,j2-1 |
472 |
do i=1,IMR |
473 |
if(abs(CRX(i,j)).GT.1.) then |
474 |
JN = j |
475 |
go to 2233 |
476 |
endif |
477 |
enddo |
478 |
enddo |
479 |
2233 continue |
480 |
C |
481 |
if(j1.ne.2) then ! Enlarged polar cap. |
482 |
do i=1,IMR |
483 |
DPI(i, 2,k) = 0. |
484 |
DPI(i,JMR,k) = 0. |
485 |
enddo |
486 |
endif |
487 |
C |
488 |
C ******* Compute horizontal mass fluxes ************ |
489 |
C |
490 |
C N-S component |
491 |
do j=j1,j2+1 |
492 |
D5 = 0.5 * COSE(j) |
493 |
do i=1,IMR |
494 |
ymass(i,j) = CRY(i,j)*D5*(delp2(i,j,k) + delp2(i,j-1,k)) |
495 |
enddo |
496 |
enddo |
497 |
C |
498 |
do 95 j=j1,j2 |
499 |
DO 95 i=1,IMR |
500 |
95 DPI(i,j,k) = (ymass(i,j) - ymass(i,j+1)) * acosp(j) |
501 |
C |
502 |
C Poles |
503 |
sum1 = ymass(IMR,j1 ) |
504 |
sum2 = ymass(IMR,J2+1) |
505 |
do i=1,IMR-1 |
506 |
sum1 = sum1 + ymass(i,j1 ) |
507 |
sum2 = sum2 + ymass(i,J2+1) |
508 |
enddo |
509 |
C |
510 |
sum1 = - sum1 * RCAP |
511 |
sum2 = sum2 * RCAP |
512 |
do i=1,IMR |
513 |
DPI(i, 1,k) = sum1 |
514 |
DPI(i,JNP,k) = sum2 |
515 |
enddo |
516 |
C |
517 |
C E-W component |
518 |
C |
519 |
do j=j1,j2 |
520 |
do i=2,IMR |
521 |
PU(i,j) = 0.5 * (delp2(i,j,k) + delp2(i-1,j,k)) |
522 |
enddo |
523 |
enddo |
524 |
C |
525 |
do j=j1,j2 |
526 |
PU(1,j) = 0.5 * (delp2(1,j,k) + delp2(IMR,j,k)) |
527 |
enddo |
528 |
C |
529 |
do 110 j=j1,j2 |
530 |
DO 110 i=1,IMR |
531 |
110 xmass(i,j) = PU(i,j)*CRX(i,j) |
532 |
C |
533 |
DO 120 j=j1,j2 |
534 |
DO 120 i=1,IMR-1 |
535 |
120 DPI(i,j,k) = DPI(i,j,k) + xmass(i,j) - xmass(i+1,j) |
536 |
C |
537 |
DO 130 j=j1,j2 |
538 |
130 DPI(IMR,j,k) = DPI(IMR,j,k) + xmass(IMR,j) - xmass(1,j) |
539 |
C |
540 |
DO j=j1,j2 |
541 |
do i=1,IMR-1 |
542 |
UA(i,j) = 0.5 * (CRX(i,j)+CRX(i+1,j)) |
543 |
enddo |
544 |
enddo |
545 |
C |
546 |
DO j=j1,j2 |
547 |
UA(imr,j) = 0.5 * (CRX(imr,j)+CRX(1,j)) |
548 |
enddo |
549 |
ccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
550 |
c Rajouts pour LMDZ.3.3 |
551 |
ccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
552 |
do i=1,IMR |
553 |
do j=1,JNP |
554 |
VA(i,j)=0. |
555 |
enddo |
556 |
enddo |
557 |
|
558 |
do i=1,imr*(JMR-1) |
559 |
VA(i,2) = 0.5*(CRY(i,2)+CRY(i,3)) |
560 |
enddo |
561 |
C |
562 |
if(j1.eq.2) then |
563 |
IMH = IMR/2 |
564 |
do i=1,IMH |
565 |
VA(i, 1) = 0.5*(CRY(i,2)-CRY(i+IMH,2)) |
566 |
VA(i+IMH, 1) = -VA(i,1) |
567 |
VA(i, JNP) = 0.5*(CRY(i,JNP)-CRY(i+IMH,JMR)) |
568 |
VA(i+IMH,JNP) = -VA(i,JNP) |
569 |
enddo |
570 |
VA(IMR,1)=VA(1,1) |
571 |
VA(IMR,JNP)=VA(1,JNP) |
572 |
endif |
573 |
C |
574 |
C ****6***0*********0*********0*********0*********0*********0**********72 |
575 |
do 1000 IC=1,NC |
576 |
C |
577 |
do i=1,IMJM |
578 |
wk1(i,1,1) = 0. |
579 |
wk1(i,1,2) = 0. |
580 |
enddo |
581 |
C |
582 |
C E-W advective cross term |
583 |
do 250 j=J1,J2 |
584 |
if(J.GT.JS .and. J.LT.JN) GO TO 250 |
585 |
C |
586 |
do i=1,IMR |
587 |
qtmp(i) = q(i,j,k,IC) |
588 |
enddo |
589 |
C |
590 |
do i=-IML,0 |
591 |
qtmp(i) = q(IMR+i,j,k,IC) |
592 |
qtmp(IMR+1-i) = q(1-i,j,k,IC) |
593 |
enddo |
594 |
C |
595 |
DO 230 i=1,IMR |
596 |
iu = UA(i,j) |
597 |
ru = UA(i,j) - iu |
598 |
iiu = i-iu |
599 |
if(UA(i,j).GE.0.) then |
600 |
wk1(i,j,1) = qtmp(iiu)+ru*(qtmp(iiu-1)-qtmp(iiu)) |
601 |
else |
602 |
wk1(i,j,1) = qtmp(iiu)+ru*(qtmp(iiu)-qtmp(iiu+1)) |
603 |
endif |
604 |
wk1(i,j,1) = wk1(i,j,1) - qtmp(i) |
605 |
230 continue |
606 |
250 continue |
607 |
C |
608 |
if(JN.ne.0) then |
609 |
do j=JS+1,JN-1 |
610 |
C |
611 |
do i=1,IMR |
612 |
qtmp(i) = q(i,j,k,IC) |
613 |
enddo |
614 |
C |
615 |
qtmp(0) = q(IMR,J,k,IC) |
616 |
qtmp(IMR+1) = q( 1,J,k,IC) |
617 |
C |
618 |
do i=1,imr |
619 |
iu = i - UA(i,j) |
620 |
wk1(i,j,1) = UA(i,j)*(qtmp(iu) - qtmp(iu+1)) |
621 |
enddo |
622 |
enddo |
623 |
endif |
624 |
C ****6***0*********0*********0*********0*********0*********0**********72 |
625 |
C Contribution from the N-S advection |
626 |
do i=1,imr*(j2-j1+1) |
627 |
JT = float(J1) - VA(i,j1) |
628 |
wk1(i,j1,2) = VA(i,j1) * (q(i,jt,k,IC) - q(i,jt+1,k,IC)) |
629 |
enddo |
630 |
C |
631 |
do i=1,IMJM |
632 |
wk1(i,1,1) = q(i,1,k,IC) + 0.5*wk1(i,1,1) |
633 |
wk1(i,1,2) = q(i,1,k,IC) + 0.5*wk1(i,1,2) |
634 |
enddo |
635 |
C |
636 |
if(cross) then |
637 |
C Add cross terms in the vertical direction. |
638 |
if(IORD .GE. 2) then |
639 |
iad = 2 |
640 |
else |
641 |
iad = 1 |
642 |
endif |
643 |
C |
644 |
if(JORD .GE. 2) then |
645 |
jad = 2 |
646 |
else |
647 |
jad = 1 |
648 |
endif |
649 |
call xadv(IMR,JNP,j1,j2,wk1(1,1,2),UA,JS,JN,IML,DC2,iad) |
650 |
call yadv(IMR,JNP,j1,j2,wk1(1,1,1),VA,PV,W,jad) |
651 |
do j=1,JNP |
652 |
do i=1,IMR |
653 |
q(i,j,k,IC) = q(i,j,k,IC) + DC2(i,j) + PV(i,j) |
654 |
enddo |
655 |
enddo |
656 |
endif |
657 |
C |
658 |
call xtp(IMR,JNP,IML,j1,j2,JN,JS,PU,DQ(1,1,k,IC),wk1(1,1,2) |
659 |
& ,CRX,fx1,xmass,IORD) |
660 |
|
661 |
call ytp(IMR,JNP,j1,j2,acosp,RCAP,DQ(1,1,k,IC),wk1(1,1,1),CRY, |
662 |
& DC2,ymass,WK1(1,1,3),wk1(1,1,4),WK1(1,1,5),WK1(1,1,6),JORD) |
663 |
C |
664 |
1000 continue |
665 |
1500 continue |
666 |
C |
667 |
C ******* Compute vertical mass flux (same unit as PS) *********** |
668 |
C |
669 |
C 1st step: compute total column mass CONVERGENCE. |
670 |
C |
671 |
do 320 j=1,JNP |
672 |
do 320 i=1,IMR |
673 |
320 CRY(i,j) = DPI(i,j,1) |
674 |
C |
675 |
do 330 k=2,NLAY |
676 |
do 330 j=1,JNP |
677 |
do 330 i=1,IMR |
678 |
CRY(i,j) = CRY(i,j) + DPI(i,j,k) |
679 |
330 continue |
680 |
C |
681 |
do 360 j=1,JNP |
682 |
do 360 i=1,IMR |
683 |
C |
684 |
C 2nd step: compute PS2 (PS at n+1) using the hydrostatic assumption. |
685 |
C Changes (increases) to surface pressure = total column mass convergence |
686 |
C |
687 |
PS2(i,j) = PS1(i,j) + CRY(i,j) |
688 |
C |
689 |
C 3rd step: compute vertical mass flux from mass conservation principle. |
690 |
C |
691 |
W(i,j,1) = DPI(i,j,1) - DBK(1)*CRY(i,j) |
692 |
W(i,j,NLAY) = 0. |
693 |
360 continue |
694 |
C |
695 |
do 370 k=2,NLAY-1 |
696 |
do 370 j=1,JNP |
697 |
do 370 i=1,IMR |
698 |
W(i,j,k) = W(i,j,k-1) + DPI(i,j,k) - DBK(k)*CRY(i,j) |
699 |
370 continue |
700 |
C |
701 |
DO 380 k=1,NLAY |
702 |
DO 380 j=1,JNP |
703 |
DO 380 i=1,IMR |
704 |
delp2(i,j,k) = DAP(k) + DBK(k)*PS2(i,j) |
705 |
380 continue |
706 |
C |
707 |
KRD = max(3, KORD) |
708 |
do 4000 IC=1,NC |
709 |
C |
710 |
C****6***0*********0*********0*********0*********0*********0**********72 |
711 |
|
712 |
call FZPPM(IMR,JNP,NLAY,j1,DQ(1,1,1,IC),W,Q(1,1,1,IC),WK1,DPI, |
713 |
& DC2,CRX,CRY,PU,PV,xmass,ymass,delp1,KRD) |
714 |
C |
715 |
|
716 |
if(fill) call qckxyz(DQ(1,1,1,IC),DC2,IMR,JNP,NLAY,j1,j2, |
717 |
& cosp,acosp,.false.,IC,NSTEP) |
718 |
C |
719 |
C Recover tracer mixing ratio from "density" using predicted |
720 |
C "air density" (pressure thickness) at time-level n+1 |
721 |
C |
722 |
DO k=1,NLAY |
723 |
DO j=1,JNP |
724 |
DO i=1,IMR |
725 |
Q(i,j,k,IC) = DQ(i,j,k,IC) / delp2(i,j,k) |
726 |
enddo |
727 |
enddo |
728 |
enddo |
729 |
C |
730 |
if(j1.ne.2) then |
731 |
DO 400 k=1,NLAY |
732 |
DO 400 I=1,IMR |
733 |
c j=1 c'est le pôle Sud, j=JNP c'est le pôle Nord |
734 |
Q(I, 2,k,IC) = Q(I, 1,k,IC) |
735 |
Q(I,JMR,k,IC) = Q(I,JMP,k,IC) |
736 |
400 CONTINUE |
737 |
endif |
738 |
4000 continue |
739 |
C |
740 |
if(j1.ne.2) then |
741 |
DO 5000 k=1,NLAY |
742 |
DO 5000 i=1,IMR |
743 |
W(i, 2,k) = W(i, 1,k) |
744 |
W(i,JMR,k) = W(i,JNP,k) |
745 |
5000 continue |
746 |
endif |
747 |
C |
748 |
RETURN |
749 |
END |
750 |
C |
751 |
C****6***0*********0*********0*********0*********0*********0**********72 |
752 |
subroutine FZPPM(IMR,JNP,NLAY,j1,DQ,WZ,P,DC,DQDT,AR,AL,A6, |
753 |
& flux,wk1,wk2,wz2,delp,KORD) |
754 |
parameter ( kmax = 150 ) |
755 |
parameter ( R23 = 2./3., R3 = 1./3.) |
756 |
real WZ(IMR,JNP,NLAY),P(IMR,JNP,NLAY),DC(IMR,JNP,NLAY), |
757 |
& wk1(IMR,*),delp(IMR,JNP,NLAY),DQ(IMR,JNP,NLAY), |
758 |
& DQDT(IMR,JNP,NLAY) |
759 |
C Assuming JNP >= NLAY |
760 |
real AR(IMR,*),AL(IMR,*),A6(IMR,*),flux(IMR,*),wk2(IMR,*), |
761 |
& wz2(IMR,*) |
762 |
C |
763 |
JMR = JNP - 1 |
764 |
IMJM = IMR*JNP |
765 |
NLAYM1 = NLAY - 1 |
766 |
C |
767 |
LMT = KORD - 3 |
768 |
C |
769 |
C ****6***0*********0*********0*********0*********0*********0**********72 |
770 |
C Compute DC for PPM |
771 |
C ****6***0*********0*********0*********0*********0*********0**********72 |
772 |
C |
773 |
do 1000 k=1,NLAYM1 |
774 |
do 1000 i=1,IMJM |
775 |
DQDT(i,1,k) = P(i,1,k+1) - P(i,1,k) |
776 |
1000 continue |
777 |
C |
778 |
DO 1220 k=2,NLAYM1 |
779 |
DO 1220 I=1,IMJM |
780 |
c0 = delp(i,1,k) / (delp(i,1,k-1)+delp(i,1,k)+delp(i,1,k+1)) |
781 |
c1 = (delp(i,1,k-1)+0.5*delp(i,1,k))/(delp(i,1,k+1)+delp(i,1,k)) |
782 |
c2 = (delp(i,1,k+1)+0.5*delp(i,1,k))/(delp(i,1,k-1)+delp(i,1,k)) |
783 |
tmp = c0*(c1*DQDT(i,1,k) + c2*DQDT(i,1,k-1)) |
784 |
Qmax = max(P(i,1,k-1),P(i,1,k),P(i,1,k+1)) - P(i,1,k) |
785 |
Qmin = P(i,1,k) - min(P(i,1,k-1),P(i,1,k),P(i,1,k+1)) |
786 |
DC(i,1,k) = sign(min(abs(tmp),Qmax,Qmin), tmp) |
787 |
1220 CONTINUE |
788 |
|
789 |
C |
790 |
C ****6***0*********0*********0*********0*********0*********0**********72 |
791 |
C Loop over latitudes (to save memory) |
792 |
C ****6***0*********0*********0*********0*********0*********0**********72 |
793 |
C |
794 |
DO 2000 j=1,JNP |
795 |
if((j.eq.2 .or. j.eq.JMR) .and. j1.ne.2) goto 2000 |
796 |
C |
797 |
DO k=1,NLAY |
798 |
DO i=1,IMR |
799 |
wz2(i,k) = WZ(i,j,k) |
800 |
wk1(i,k) = P(i,j,k) |
801 |
wk2(i,k) = delp(i,j,k) |
802 |
flux(i,k) = DC(i,j,k) !this flux is actually the monotone slope |
803 |
enddo |
804 |
enddo |
805 |
C |
806 |
C****6***0*********0*********0*********0*********0*********0**********72 |
807 |
C Compute first guesses at cell interfaces |
808 |
C First guesses are required to be continuous. |
809 |
C ****6***0*********0*********0*********0*********0*********0**********72 |
810 |
C |
811 |
C three-cell parabolic subgrid distribution at model top |
812 |
C two-cell parabolic with zero gradient subgrid distribution |
813 |
C at the surface. |
814 |
C |
815 |
C First guess top edge value |
816 |
DO 10 i=1,IMR |
817 |
C three-cell PPM |
818 |
C Compute a,b, and c of q = aP**2 + bP + c using cell averages and delp |
819 |
a = 3.*( DQDT(i,j,2) - DQDT(i,j,1)*(wk2(i,2)+wk2(i,3))/ |
820 |
& (wk2(i,1)+wk2(i,2)) ) / |
821 |
& ( (wk2(i,2)+wk2(i,3))*(wk2(i,1)+wk2(i,2)+wk2(i,3)) ) |
822 |
b = 2.*DQDT(i,j,1)/(wk2(i,1)+wk2(i,2)) - |
823 |
& R23*a*(2.*wk2(i,1)+wk2(i,2)) |
824 |
AL(i,1) = wk1(i,1) - wk2(i,1)*(R3*a*wk2(i,1) + 0.5*b) |
825 |
AL(i,2) = wk2(i,1)*(a*wk2(i,1) + b) + AL(i,1) |
826 |
C |
827 |
C Check if change sign |
828 |
if(wk1(i,1)*AL(i,1).le.0.) then |
829 |
AL(i,1) = 0. |
830 |
flux(i,1) = 0. |
831 |
else |
832 |
flux(i,1) = wk1(i,1) - AL(i,1) |
833 |
endif |
834 |
10 continue |
835 |
C |
836 |
C Bottom |
837 |
DO 15 i=1,IMR |
838 |
C 2-cell PPM with zero gradient right at the surface |
839 |
C |
840 |
fct = DQDT(i,j,NLAYM1)*wk2(i,NLAY)**2 / |
841 |
& ( (wk2(i,NLAY)+wk2(i,NLAYM1))*(2.*wk2(i,NLAY)+wk2(i,NLAYM1))) |
842 |
AR(i,NLAY) = wk1(i,NLAY) + fct |
843 |
AL(i,NLAY) = wk1(i,NLAY) - (fct+fct) |
844 |
if(wk1(i,NLAY)*AR(i,NLAY).le.0.) AR(i,NLAY) = 0. |
845 |
flux(i,NLAY) = AR(i,NLAY) - wk1(i,NLAY) |
846 |
15 continue |
847 |
|
848 |
C |
849 |
C****6***0*********0*********0*********0*********0*********0**********72 |
850 |
C 4th order interpolation in the interior. |
851 |
C****6***0*********0*********0*********0*********0*********0**********72 |
852 |
C |
853 |
DO 14 k=3,NLAYM1 |
854 |
DO 12 i=1,IMR |
855 |
c1 = DQDT(i,j,k-1)*wk2(i,k-1) / (wk2(i,k-1)+wk2(i,k)) |
856 |
c2 = 2. / (wk2(i,k-2)+wk2(i,k-1)+wk2(i,k)+wk2(i,k+1)) |
857 |
A1 = (wk2(i,k-2)+wk2(i,k-1)) / (2.*wk2(i,k-1)+wk2(i,k)) |
858 |
A2 = (wk2(i,k )+wk2(i,k+1)) / (2.*wk2(i,k)+wk2(i,k-1)) |
859 |
AL(i,k) = wk1(i,k-1) + c1 + c2 * |
860 |
& ( wk2(i,k )*(c1*(A1 - A2)+A2*flux(i,k-1)) - |
861 |
& wk2(i,k-1)*A1*flux(i,k) ) |
862 |
12 CONTINUE |
863 |
14 continue |
864 |
C |
865 |
do 20 i=1,IMR*NLAYM1 |
866 |
AR(i,1) = AL(i,2) |
867 |
20 continue |
868 |
C |
869 |
do 30 i=1,IMR*NLAY |
870 |
A6(i,1) = 3.*(wk1(i,1)+wk1(i,1) - (AL(i,1)+AR(i,1))) |
871 |
30 continue |
872 |
C |
873 |
C****6***0*********0*********0*********0*********0*********0**********72 |
874 |
C Top & Bot always monotonic |
875 |
call lmtppm(flux(1,1),A6(1,1),AR(1,1),AL(1,1),wk1(1,1),IMR,0) |
876 |
call lmtppm(flux(1,NLAY),A6(1,NLAY),AR(1,NLAY),AL(1,NLAY), |
877 |
& wk1(1,NLAY),IMR,0) |
878 |
C |
879 |
C Interior depending on KORD |
880 |
if(LMT.LE.2) |
881 |
& call lmtppm(flux(1,2),A6(1,2),AR(1,2),AL(1,2),wk1(1,2), |
882 |
& IMR*(NLAY-2),LMT) |
883 |
C |
884 |
C****6***0*********0*********0*********0*********0*********0**********72 |
885 |
C |
886 |
DO 140 i=1,IMR*NLAYM1 |
887 |
IF(wz2(i,1).GT.0.) then |
888 |
CM = wz2(i,1) / wk2(i,1) |
889 |
flux(i,2) = AR(i,1)+0.5*CM*(AL(i,1)-AR(i,1)+A6(i,1)*(1.-R23*CM)) |
890 |
else |
891 |
CP= wz2(i,1) / wk2(i,2) |
892 |
flux(i,2) = AL(i,2)+0.5*CP*(AL(i,2)-AR(i,2)-A6(i,2)*(1.+R23*CP)) |
893 |
endif |
894 |
140 continue |
895 |
C |
896 |
DO 250 i=1,IMR*NLAYM1 |
897 |
flux(i,2) = wz2(i,1) * flux(i,2) |
898 |
250 continue |
899 |
C |
900 |
do 350 i=1,IMR |
901 |
DQ(i,j, 1) = DQ(i,j, 1) - flux(i, 2) |
902 |
DQ(i,j,NLAY) = DQ(i,j,NLAY) + flux(i,NLAY) |
903 |
350 continue |
904 |
C |
905 |
do 360 k=2,NLAYM1 |
906 |
do 360 i=1,IMR |
907 |
360 DQ(i,j,k) = DQ(i,j,k) + flux(i,k) - flux(i,k+1) |
908 |
2000 continue |
909 |
return |
910 |
end |
911 |
C |
912 |
subroutine xtp(IMR,JNP,IML,j1,j2,JN,JS,PU,DQ,Q,UC, |
913 |
& fx1,xmass,IORD) |
914 |
dimension UC(IMR,*),DC(-IML:IMR+IML+1),xmass(IMR,JNP) |
915 |
& ,fx1(IMR+1),DQ(IMR,JNP),qtmp(-IML:IMR+1+IML) |
916 |
dimension PU(IMR,JNP),Q(IMR,JNP),ISAVE(IMR) |
917 |
C |
918 |
IMP = IMR + 1 |
919 |
C |
920 |
C van Leer at high latitudes |
921 |
jvan = max(1,JNP/18) |
922 |
j1vl = j1+jvan |
923 |
j2vl = j2-jvan |
924 |
C |
925 |
do 1310 j=j1,j2 |
926 |
C |
927 |
do i=1,IMR |
928 |
qtmp(i) = q(i,j) |
929 |
enddo |
930 |
C |
931 |
if(j.ge.JN .or. j.le.JS) goto 2222 |
932 |
C ************* Eulerian ********** |
933 |
C |
934 |
qtmp(0) = q(IMR,J) |
935 |
qtmp(-1) = q(IMR-1,J) |
936 |
qtmp(IMP) = q(1,J) |
937 |
qtmp(IMP+1) = q(2,J) |
938 |
C |
939 |
IF(IORD.eq.1 .or. j.eq.j1. or. j.eq.j2) THEN |
940 |
DO 1406 i=1,IMR |
941 |
iu = float(i) - uc(i,j) |
942 |
1406 fx1(i) = qtmp(iu) |
943 |
ELSE |
944 |
call xmist(IMR,IML,Qtmp,DC) |
945 |
DC(0) = DC(IMR) |
946 |
C |
947 |
if(IORD.eq.2 .or. j.le.j1vl .or. j.ge.j2vl) then |
948 |
DO 1408 i=1,IMR |
949 |
iu = float(i) - uc(i,j) |
950 |
1408 fx1(i) = qtmp(iu) + DC(iu)*(sign(1.,uc(i,j))-uc(i,j)) |
951 |
else |
952 |
call fxppm(IMR,IML,UC(1,j),Qtmp,DC,fx1,IORD) |
953 |
endif |
954 |
C |
955 |
ENDIF |
956 |
C |
957 |
DO 1506 i=1,IMR |
958 |
1506 fx1(i) = fx1(i)*xmass(i,j) |
959 |
C |
960 |
goto 1309 |
961 |
C |
962 |
C ***** Conservative (flux-form) Semi-Lagrangian transport ***** |
963 |
C |
964 |
2222 continue |
965 |
C |
966 |
do i=-IML,0 |
967 |
qtmp(i) = q(IMR+i,j) |
968 |
qtmp(IMP-i) = q(1-i,j) |
969 |
enddo |
970 |
C |
971 |
IF(IORD.eq.1 .or. j.eq.j1. or. j.eq.j2) THEN |
972 |
DO 1306 i=1,IMR |
973 |
itmp = INT(uc(i,j)) |
974 |
ISAVE(i) = i - itmp |
975 |
iu = i - uc(i,j) |
976 |
1306 fx1(i) = (uc(i,j) - itmp)*qtmp(iu) |
977 |
ELSE |
978 |
call xmist(IMR,IML,Qtmp,DC) |
979 |
C |
980 |
do i=-IML,0 |
981 |
DC(i) = DC(IMR+i) |
982 |
DC(IMP-i) = DC(1-i) |
983 |
enddo |
984 |
C |
985 |
DO 1307 i=1,IMR |
986 |
itmp = INT(uc(i,j)) |
987 |
rut = uc(i,j) - itmp |
988 |
ISAVE(i) = i - itmp |
989 |
iu = i - uc(i,j) |
990 |
1307 fx1(i) = rut*(qtmp(iu) + DC(iu)*(sign(1.,rut) - rut)) |
991 |
ENDIF |
992 |
C |
993 |
do 1308 i=1,IMR |
994 |
IF(uc(i,j).GT.1.) then |
995 |
CDIR$ NOVECTOR |
996 |
do ist = ISAVE(i),i-1 |
997 |
fx1(i) = fx1(i) + qtmp(ist) |
998 |
enddo |
999 |
elseIF(uc(i,j).LT.-1.) then |
1000 |
do ist = i,ISAVE(i)-1 |
1001 |
fx1(i) = fx1(i) - qtmp(ist) |
1002 |
enddo |
1003 |
CDIR$ VECTOR |
1004 |
endif |
1005 |
1308 continue |
1006 |
do i=1,IMR |
1007 |
fx1(i) = PU(i,j)*fx1(i) |
1008 |
enddo |
1009 |
C |
1010 |
C *************************************** |
1011 |
C |
1012 |
1309 fx1(IMP) = fx1(1) |
1013 |
DO 1215 i=1,IMR |
1014 |
1215 DQ(i,j) = DQ(i,j) + fx1(i)-fx1(i+1) |
1015 |
C |
1016 |
C *************************************** |
1017 |
C |
1018 |
1310 continue |
1019 |
return |
1020 |
end |
1021 |
C |
1022 |
subroutine fxppm(IMR,IML,UT,P,DC,flux,IORD) |
1023 |
parameter ( R3 = 1./3., R23 = 2./3. ) |
1024 |
DIMENSION UT(*),flux(*),P(-IML:IMR+IML+1),DC(-IML:IMR+IML+1) |
1025 |
DIMENSION AR(0:IMR),AL(0:IMR),A6(0:IMR) |
1026 |
integer LMT |
1027 |
c logical first |
1028 |
c data first /.true./ |
1029 |
c SAVE LMT |
1030 |
c if(first) then |
1031 |
C |
1032 |
C correction calcul de LMT a chaque passage pour pouvoir choisir |
1033 |
c plusieurs schemas PPM pour differents traceurs |
1034 |
c IF (IORD.LE.0) then |
1035 |
c if(IMR.GE.144) then |
1036 |
c LMT = 0 |
1037 |
c elseif(IMR.GE.72) then |
1038 |
c LMT = 1 |
1039 |
c else |
1040 |
c LMT = 2 |
1041 |
c endif |
1042 |
c else |
1043 |
c LMT = IORD - 3 |
1044 |
c endif |
1045 |
C |
1046 |
LMT = IORD - 3 |
1047 |
|
1048 |
DO 10 i=1,IMR |
1049 |
10 AL(i) = 0.5*(p(i-1)+p(i)) + (DC(i-1) - DC(i))*R3 |
1050 |
C |
1051 |
do 20 i=1,IMR-1 |
1052 |
20 AR(i) = AL(i+1) |
1053 |
AR(IMR) = AL(1) |
1054 |
C |
1055 |
do 30 i=1,IMR |
1056 |
30 A6(i) = 3.*(p(i)+p(i) - (AL(i)+AR(i))) |
1057 |
C |
1058 |
if(LMT.LE.2) call lmtppm(DC(1),A6(1),AR(1),AL(1),P(1),IMR,LMT) |
1059 |
C |
1060 |
AL(0) = AL(IMR) |
1061 |
AR(0) = AR(IMR) |
1062 |
A6(0) = A6(IMR) |
1063 |
C |
1064 |
DO i=1,IMR |
1065 |
IF(UT(i).GT.0.) then |
1066 |
flux(i) = AR(i-1) + 0.5*UT(i)*(AL(i-1) - AR(i-1) + |
1067 |
& A6(i-1)*(1.-R23*UT(i)) ) |
1068 |
else |
1069 |
flux(i) = AL(i) - 0.5*UT(i)*(AR(i) - AL(i) + |
1070 |
& A6(i)*(1.+R23*UT(i))) |
1071 |
endif |
1072 |
enddo |
1073 |
return |
1074 |
end |
1075 |
C |
1076 |
subroutine xmist(IMR,IML,P,DC) |
1077 |
parameter( R24 = 1./24.) |
1078 |
dimension P(-IML:IMR+1+IML),DC(-IML:IMR+1+IML) |
1079 |
C |
1080 |
do 10 i=1,IMR |
1081 |
tmp = R24*(8.*(p(i+1) - p(i-1)) + p(i-2) - p(i+2)) |
1082 |
Pmax = max(P(i-1), p(i), p(i+1)) - p(i) |
1083 |
Pmin = p(i) - min(P(i-1), p(i), p(i+1)) |
1084 |
10 DC(i) = sign(min(abs(tmp),Pmax,Pmin), tmp) |
1085 |
return |
1086 |
end |
1087 |
C |
1088 |
subroutine ytp(IMR,JNP,j1,j2,acosp,RCAP,DQ,P,VC,DC2 |
1089 |
& ,ymass,fx,A6,AR,AL,JORD) |
1090 |
dimension P(IMR,JNP),VC(IMR,JNP),ymass(IMR,JNP) |
1091 |
& ,DC2(IMR,JNP),DQ(IMR,JNP),acosp(JNP) |
1092 |
C Work array |
1093 |
DIMENSION fx(IMR,JNP),AR(IMR,JNP),AL(IMR,JNP),A6(IMR,JNP) |
1094 |
C |
1095 |
JMR = JNP - 1 |
1096 |
len = IMR*(J2-J1+2) |
1097 |
C |
1098 |
if(JORD.eq.1) then |
1099 |
DO 1000 i=1,len |
1100 |
JT = float(J1) - VC(i,J1) |
1101 |
1000 fx(i,j1) = p(i,JT) |
1102 |
else |
1103 |
|
1104 |
call ymist(IMR,JNP,j1,P,DC2,4) |
1105 |
C |
1106 |
if(JORD.LE.0 .or. JORD.GE.3) then |
1107 |
|
1108 |
call fyppm(VC,P,DC2,fx,IMR,JNP,j1,j2,A6,AR,AL,JORD) |
1109 |
|
1110 |
else |
1111 |
DO 1200 i=1,len |
1112 |
JT = float(J1) - VC(i,J1) |
1113 |
1200 fx(i,j1) = p(i,JT) + (sign(1.,VC(i,j1))-VC(i,j1))*DC2(i,JT) |
1114 |
endif |
1115 |
endif |
1116 |
C |
1117 |
DO 1300 i=1,len |
1118 |
1300 fx(i,j1) = fx(i,j1)*ymass(i,j1) |
1119 |
C |
1120 |
DO 1400 j=j1,j2 |
1121 |
DO 1400 i=1,IMR |
1122 |
1400 DQ(i,j) = DQ(i,j) + (fx(i,j) - fx(i,j+1)) * acosp(j) |
1123 |
C |
1124 |
C Poles |
1125 |
sum1 = fx(IMR,j1 ) |
1126 |
sum2 = fx(IMR,J2+1) |
1127 |
do i=1,IMR-1 |
1128 |
sum1 = sum1 + fx(i,j1 ) |
1129 |
sum2 = sum2 + fx(i,J2+1) |
1130 |
enddo |
1131 |
C |
1132 |
sum1 = DQ(1, 1) - sum1 * RCAP |
1133 |
sum2 = DQ(1,JNP) + sum2 * RCAP |
1134 |
do i=1,IMR |
1135 |
DQ(i, 1) = sum1 |
1136 |
DQ(i,JNP) = sum2 |
1137 |
enddo |
1138 |
C |
1139 |
if(j1.ne.2) then |
1140 |
do i=1,IMR |
1141 |
DQ(i, 2) = sum1 |
1142 |
DQ(i,JMR) = sum2 |
1143 |
enddo |
1144 |
endif |
1145 |
C |
1146 |
return |
1147 |
end |
1148 |
C |
1149 |
subroutine ymist(IMR,JNP,j1,P,DC,ID) |
1150 |
parameter ( R24 = 1./24. ) |
1151 |
dimension P(IMR,JNP),DC(IMR,JNP) |
1152 |
C |
1153 |
IMH = IMR / 2 |
1154 |
JMR = JNP - 1 |
1155 |
IJM3 = IMR*(JMR-3) |
1156 |
C |
1157 |
IF(ID.EQ.2) THEN |
1158 |
do 10 i=1,IMR*(JMR-1) |
1159 |
tmp = 0.25*(p(i,3) - p(i,1)) |
1160 |
Pmax = max(p(i,1),p(i,2),p(i,3)) - p(i,2) |
1161 |
Pmin = p(i,2) - min(p(i,1),p(i,2),p(i,3)) |
1162 |
DC(i,2) = sign(min(abs(tmp),Pmin,Pmax),tmp) |
1163 |
10 CONTINUE |
1164 |
ELSE |
1165 |
do 12 i=1,IMH |
1166 |
C J=2 |
1167 |
tmp = (8.*(p(i,3) - p(i,1)) + p(i+IMH,2) - p(i,4))*R24 |
1168 |
Pmax = max(p(i,1),p(i,2),p(i,3)) - p(i,2) |
1169 |
Pmin = p(i,2) - min(p(i,1),p(i,2),p(i,3)) |
1170 |
DC(i,2) = sign(min(abs(tmp),Pmin,Pmax),tmp) |
1171 |
C J=JMR |
1172 |
tmp=(8.*(p(i,JNP)-p(i,JMR-1))+p(i,JMR-2)-p(i+IMH,JMR))*R24 |
1173 |
Pmax = max(p(i,JMR-1),p(i,JMR),p(i,JNP)) - p(i,JMR) |
1174 |
Pmin = p(i,JMR) - min(p(i,JMR-1),p(i,JMR),p(i,JNP)) |
1175 |
DC(i,JMR) = sign(min(abs(tmp),Pmin,Pmax),tmp) |
1176 |
12 CONTINUE |
1177 |
do 14 i=IMH+1,IMR |
1178 |
C J=2 |
1179 |
tmp = (8.*(p(i,3) - p(i,1)) + p(i-IMH,2) - p(i,4))*R24 |
1180 |
Pmax = max(p(i,1),p(i,2),p(i,3)) - p(i,2) |
1181 |
Pmin = p(i,2) - min(p(i,1),p(i,2),p(i,3)) |
1182 |
DC(i,2) = sign(min(abs(tmp),Pmin,Pmax),tmp) |
1183 |
C J=JMR |
1184 |
tmp=(8.*(p(i,JNP)-p(i,JMR-1))+p(i,JMR-2)-p(i-IMH,JMR))*R24 |
1185 |
Pmax = max(p(i,JMR-1),p(i,JMR),p(i,JNP)) - p(i,JMR) |
1186 |
Pmin = p(i,JMR) - min(p(i,JMR-1),p(i,JMR),p(i,JNP)) |
1187 |
DC(i,JMR) = sign(min(abs(tmp),Pmin,Pmax),tmp) |
1188 |
14 CONTINUE |
1189 |
C |
1190 |
do 15 i=1,IJM3 |
1191 |
tmp = (8.*(p(i,4) - p(i,2)) + p(i,1) - p(i,5))*R24 |
1192 |
Pmax = max(p(i,2),p(i,3),p(i,4)) - p(i,3) |
1193 |
Pmin = p(i,3) - min(p(i,2),p(i,3),p(i,4)) |
1194 |
DC(i,3) = sign(min(abs(tmp),Pmin,Pmax),tmp) |
1195 |
15 CONTINUE |
1196 |
ENDIF |
1197 |
C |
1198 |
if(j1.ne.2) then |
1199 |
do i=1,IMR |
1200 |
DC(i,1) = 0. |
1201 |
DC(i,JNP) = 0. |
1202 |
enddo |
1203 |
else |
1204 |
C Determine slopes in polar caps for scalars! |
1205 |
C |
1206 |
do 13 i=1,IMH |
1207 |
C South |
1208 |
tmp = 0.25*(p(i,2) - p(i+imh,2)) |
1209 |
Pmax = max(p(i,2),p(i,1), p(i+imh,2)) - p(i,1) |
1210 |
Pmin = p(i,1) - min(p(i,2),p(i,1), p(i+imh,2)) |
1211 |
DC(i,1)=sign(min(abs(tmp),Pmax,Pmin),tmp) |
1212 |
C North. |
1213 |
tmp = 0.25*(p(i+imh,JMR) - p(i,JMR)) |
1214 |
Pmax = max(p(i+imh,JMR),p(i,jnp), p(i,JMR)) - p(i,JNP) |
1215 |
Pmin = p(i,JNP) - min(p(i+imh,JMR),p(i,jnp), p(i,JMR)) |
1216 |
DC(i,JNP) = sign(min(abs(tmp),Pmax,pmin),tmp) |
1217 |
13 continue |
1218 |
C |
1219 |
do 25 i=imh+1,IMR |
1220 |
DC(i, 1) = - DC(i-imh, 1) |
1221 |
DC(i,JNP) = - DC(i-imh,JNP) |
1222 |
25 continue |
1223 |
endif |
1224 |
return |
1225 |
end |
1226 |
C |
1227 |
subroutine fyppm(VC,P,DC,flux,IMR,JNP,j1,j2,A6,AR,AL,JORD) |
1228 |
parameter ( R3 = 1./3., R23 = 2./3. ) |
1229 |
real VC(IMR,*),flux(IMR,*),P(IMR,*),DC(IMR,*) |
1230 |
C Local work arrays. |
1231 |
real AR(IMR,JNP),AL(IMR,JNP),A6(IMR,JNP) |
1232 |
integer LMT |
1233 |
c logical first |
1234 |
C data first /.true./ |
1235 |
C SAVE LMT |
1236 |
C |
1237 |
IMH = IMR / 2 |
1238 |
JMR = JNP - 1 |
1239 |
j11 = j1-1 |
1240 |
IMJM1 = IMR*(J2-J1+2) |
1241 |
len = IMR*(J2-J1+3) |
1242 |
C if(first) then |
1243 |
C IF(JORD.LE.0) then |
1244 |
C if(JMR.GE.90) then |
1245 |
C LMT = 0 |
1246 |
C elseif(JMR.GE.45) then |
1247 |
C LMT = 1 |
1248 |
C else |
1249 |
C LMT = 2 |
1250 |
C endif |
1251 |
C else |
1252 |
C LMT = JORD - 3 |
1253 |
C endif |
1254 |
C |
1255 |
C first = .false. |
1256 |
C endif |
1257 |
C |
1258 |
c modifs pour pouvoir choisir plusieurs schemas PPM |
1259 |
LMT = JORD - 3 |
1260 |
C |
1261 |
DO 10 i=1,IMR*JMR |
1262 |
AL(i,2) = 0.5*(p(i,1)+p(i,2)) + (DC(i,1) - DC(i,2))*R3 |
1263 |
AR(i,1) = AL(i,2) |
1264 |
10 CONTINUE |
1265 |
C |
1266 |
CPoles: |
1267 |
C |
1268 |
DO i=1,IMH |
1269 |
AL(i,1) = AL(i+IMH,2) |
1270 |
AL(i+IMH,1) = AL(i,2) |
1271 |
C |
1272 |
AR(i,JNP) = AR(i+IMH,JMR) |
1273 |
AR(i+IMH,JNP) = AR(i,JMR) |
1274 |
ENDDO |
1275 |
|
1276 |
ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
1277 |
c Rajout pour LMDZ.3.3 |
1278 |
cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
1279 |
AR(IMR,1)=AL(1,1) |
1280 |
AR(IMR,JNP)=AL(1,JNP) |
1281 |
ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
1282 |
|
1283 |
|
1284 |
do 30 i=1,len |
1285 |
30 A6(i,j11) = 3.*(p(i,j11)+p(i,j11) - (AL(i,j11)+AR(i,j11))) |
1286 |
C |
1287 |
if(LMT.le.2) call lmtppm(DC(1,j11),A6(1,j11),AR(1,j11) |
1288 |
& ,AL(1,j11),P(1,j11),len,LMT) |
1289 |
C |
1290 |
|
1291 |
DO 140 i=1,IMJM1 |
1292 |
IF(VC(i,j1).GT.0.) then |
1293 |
flux(i,j1) = AR(i,j11) + 0.5*VC(i,j1)*(AL(i,j11) - AR(i,j11) + |
1294 |
& A6(i,j11)*(1.-R23*VC(i,j1)) ) |
1295 |
else |
1296 |
flux(i,j1) = AL(i,j1) - 0.5*VC(i,j1)*(AR(i,j1) - AL(i,j1) + |
1297 |
& A6(i,j1)*(1.+R23*VC(i,j1))) |
1298 |
endif |
1299 |
140 continue |
1300 |
return |
1301 |
end |
1302 |
C |
1303 |
subroutine yadv(IMR,JNP,j1,j2,p,VA,ady,wk,IAD) |
1304 |
REAL p(IMR,JNP),ady(IMR,JNP),VA(IMR,JNP) |
1305 |
REAL WK(IMR,-1:JNP+2) |
1306 |
C |
1307 |
JMR = JNP-1 |
1308 |
IMH = IMR/2 |
1309 |
do j=1,JNP |
1310 |
do i=1,IMR |
1311 |
wk(i,j) = p(i,j) |
1312 |
enddo |
1313 |
enddo |
1314 |
C Poles: |
1315 |
do i=1,IMH |
1316 |
wk(i, -1) = p(i+IMH,3) |
1317 |
wk(i+IMH,-1) = p(i,3) |
1318 |
wk(i, 0) = p(i+IMH,2) |
1319 |
wk(i+IMH,0) = p(i,2) |
1320 |
wk(i,JNP+1) = p(i+IMH,JMR) |
1321 |
wk(i+IMH,JNP+1) = p(i,JMR) |
1322 |
wk(i,JNP+2) = p(i+IMH,JNP-2) |
1323 |
wk(i+IMH,JNP+2) = p(i,JNP-2) |
1324 |
enddo |
1325 |
|
1326 |
IF(IAD.eq.2) then |
1327 |
do j=j1-1,j2+1 |
1328 |
do i=1,IMR |
1329 |
JP = NINT(VA(i,j)) |
1330 |
rv = JP - VA(i,j) |
1331 |
JP = j - JP |
1332 |
a1 = 0.5*(wk(i,jp+1)+wk(i,jp-1)) - wk(i,jp) |
1333 |
b1 = 0.5*(wk(i,jp+1)-wk(i,jp-1)) |
1334 |
ady(i,j) = wk(i,jp) + rv*(a1*rv + b1) - wk(i,j) |
1335 |
enddo |
1336 |
enddo |
1337 |
|
1338 |
ELSEIF(IAD.eq.1) then |
1339 |
do j=j1-1,j2+1 |
1340 |
do i=1,imr |
1341 |
JP = float(j)-VA(i,j) |
1342 |
ady(i,j) = VA(i,j)*(wk(i,jp)-wk(i,jp+1)) |
1343 |
enddo |
1344 |
enddo |
1345 |
ENDIF |
1346 |
C |
1347 |
if(j1.ne.2) then |
1348 |
sum1 = 0. |
1349 |
sum2 = 0. |
1350 |
do i=1,imr |
1351 |
sum1 = sum1 + ady(i,2) |
1352 |
sum2 = sum2 + ady(i,JMR) |
1353 |
enddo |
1354 |
sum1 = sum1 / IMR |
1355 |
sum2 = sum2 / IMR |
1356 |
C |
1357 |
do i=1,imr |
1358 |
ady(i, 2) = sum1 |
1359 |
ady(i,JMR) = sum2 |
1360 |
ady(i, 1) = sum1 |
1361 |
ady(i,JNP) = sum2 |
1362 |
enddo |
1363 |
else |
1364 |
C Poles: |
1365 |
sum1 = 0. |
1366 |
sum2 = 0. |
1367 |
do i=1,imr |
1368 |
sum1 = sum1 + ady(i,1) |
1369 |
sum2 = sum2 + ady(i,JNP) |
1370 |
enddo |
1371 |
sum1 = sum1 / IMR |
1372 |
sum2 = sum2 / IMR |
1373 |
C |
1374 |
do i=1,imr |
1375 |
ady(i, 1) = sum1 |
1376 |
ady(i,JNP) = sum2 |
1377 |
enddo |
1378 |
endif |
1379 |
C |
1380 |
return |
1381 |
end |
1382 |
C |
1383 |
subroutine xadv(IMR,JNP,j1,j2,p,UA,JS,JN,IML,adx,IAD) |
1384 |
REAL p(IMR,JNP),adx(IMR,JNP),qtmp(-IMR:IMR+IMR),UA(IMR,JNP) |
1385 |
C |
1386 |
JMR = JNP-1 |
1387 |
do 1309 j=j1,j2 |
1388 |
if(J.GT.JS .and. J.LT.JN) GO TO 1309 |
1389 |
C |
1390 |
do i=1,IMR |
1391 |
qtmp(i) = p(i,j) |
1392 |
enddo |
1393 |
C |
1394 |
do i=-IML,0 |
1395 |
qtmp(i) = p(IMR+i,j) |
1396 |
qtmp(IMR+1-i) = p(1-i,j) |
1397 |
enddo |
1398 |
C |
1399 |
IF(IAD.eq.2) THEN |
1400 |
DO i=1,IMR |
1401 |
IP = NINT(UA(i,j)) |
1402 |
ru = IP - UA(i,j) |
1403 |
IP = i - IP |
1404 |
a1 = 0.5*(qtmp(ip+1)+qtmp(ip-1)) - qtmp(ip) |
1405 |
b1 = 0.5*(qtmp(ip+1)-qtmp(ip-1)) |
1406 |
adx(i,j) = qtmp(ip) + ru*(a1*ru + b1) |
1407 |
enddo |
1408 |
ELSEIF(IAD.eq.1) then |
1409 |
DO i=1,IMR |
1410 |
iu = UA(i,j) |
1411 |
ru = UA(i,j) - iu |
1412 |
iiu = i-iu |
1413 |
if(UA(i,j).GE.0.) then |
1414 |
adx(i,j) = qtmp(iiu)+ru*(qtmp(iiu-1)-qtmp(iiu)) |
1415 |
else |
1416 |
adx(i,j) = qtmp(iiu)+ru*(qtmp(iiu)-qtmp(iiu+1)) |
1417 |
endif |
1418 |
enddo |
1419 |
ENDIF |
1420 |
C |
1421 |
do i=1,IMR |
1422 |
adx(i,j) = adx(i,j) - p(i,j) |
1423 |
enddo |
1424 |
1309 continue |
1425 |
C |
1426 |
C Eulerian upwind |
1427 |
C |
1428 |
do j=JS+1,JN-1 |
1429 |
C |
1430 |
do i=1,IMR |
1431 |
qtmp(i) = p(i,j) |
1432 |
enddo |
1433 |
C |
1434 |
qtmp(0) = p(IMR,J) |
1435 |
qtmp(IMR+1) = p(1,J) |
1436 |
C |
1437 |
IF(IAD.eq.2) THEN |
1438 |
qtmp(-1) = p(IMR-1,J) |
1439 |
qtmp(IMR+2) = p(2,J) |
1440 |
do i=1,imr |
1441 |
IP = NINT(UA(i,j)) |
1442 |
ru = IP - UA(i,j) |
1443 |
IP = i - IP |
1444 |
a1 = 0.5*(qtmp(ip+1)+qtmp(ip-1)) - qtmp(ip) |
1445 |
b1 = 0.5*(qtmp(ip+1)-qtmp(ip-1)) |
1446 |
adx(i,j) = qtmp(ip)- p(i,j) + ru*(a1*ru + b1) |
1447 |
enddo |
1448 |
ELSEIF(IAD.eq.1) then |
1449 |
C 1st order |
1450 |
DO i=1,IMR |
1451 |
IP = i - UA(i,j) |
1452 |
adx(i,j) = UA(i,j)*(qtmp(ip)-qtmp(ip+1)) |
1453 |
enddo |
1454 |
ENDIF |
1455 |
enddo |
1456 |
C |
1457 |
if(j1.ne.2) then |
1458 |
do i=1,IMR |
1459 |
adx(i, 2) = 0. |
1460 |
adx(i,JMR) = 0. |
1461 |
enddo |
1462 |
endif |
1463 |
C set cross term due to x-adv at the poles to zero. |
1464 |
do i=1,IMR |
1465 |
adx(i, 1) = 0. |
1466 |
adx(i,JNP) = 0. |
1467 |
enddo |
1468 |
return |
1469 |
end |
1470 |
C |
1471 |
subroutine lmtppm(DC,A6,AR,AL,P,IM,LMT) |
1472 |
C |
1473 |
C A6 = CURVATURE OF THE TEST PARABOLA |
1474 |
C AR = RIGHT EDGE VALUE OF THE TEST PARABOLA |
1475 |
C AL = LEFT EDGE VALUE OF THE TEST PARABOLA |
1476 |
C DC = 0.5 * MISMATCH |
1477 |
C P = CELL-AVERAGED VALUE |
1478 |
C IM = VECTOR LENGTH |
1479 |
C |
1480 |
C OPTIONS: |
1481 |
C |
1482 |
C LMT = 0: FULL MONOTONICITY |
1483 |
C LMT = 1: SEMI-MONOTONIC CONSTRAINT (NO UNDERSHOOTS) |
1484 |
C LMT = 2: POSITIVE-DEFINITE CONSTRAINT |
1485 |
C |
1486 |
parameter ( R12 = 1./12. ) |
1487 |
dimension A6(IM),AR(IM),AL(IM),P(IM),DC(IM) |
1488 |
C |
1489 |
if(LMT.eq.0) then |
1490 |
C Full constraint |
1491 |
do 100 i=1,IM |
1492 |
if(DC(i).eq.0.) then |
1493 |
AR(i) = p(i) |
1494 |
AL(i) = p(i) |
1495 |
A6(i) = 0. |
1496 |
else |
1497 |
da1 = AR(i) - AL(i) |
1498 |
da2 = da1**2 |
1499 |
A6DA = A6(i)*da1 |
1500 |
if(A6DA .lt. -da2) then |
1501 |
A6(i) = 3.*(AL(i)-p(i)) |
1502 |
AR(i) = AL(i) - A6(i) |
1503 |
elseif(A6DA .gt. da2) then |
1504 |
A6(i) = 3.*(AR(i)-p(i)) |
1505 |
AL(i) = AR(i) - A6(i) |
1506 |
endif |
1507 |
endif |
1508 |
100 continue |
1509 |
elseif(LMT.eq.1) then |
1510 |
C Semi-monotonic constraint |
1511 |
do 150 i=1,IM |
1512 |
if(abs(AR(i)-AL(i)) .GE. -A6(i)) go to 150 |
1513 |
if(p(i).lt.AR(i) .and. p(i).lt.AL(i)) then |
1514 |
AR(i) = p(i) |
1515 |
AL(i) = p(i) |
1516 |
A6(i) = 0. |
1517 |
elseif(AR(i) .gt. AL(i)) then |
1518 |
A6(i) = 3.*(AL(i)-p(i)) |
1519 |
AR(i) = AL(i) - A6(i) |
1520 |
else |
1521 |
A6(i) = 3.*(AR(i)-p(i)) |
1522 |
AL(i) = AR(i) - A6(i) |
1523 |
endif |
1524 |
150 continue |
1525 |
elseif(LMT.eq.2) then |
1526 |
do 250 i=1,IM |
1527 |
if(abs(AR(i)-AL(i)) .GE. -A6(i)) go to 250 |
1528 |
fmin = p(i) + 0.25*(AR(i)-AL(i))**2/A6(i) + A6(i)*R12 |
1529 |
if(fmin.ge.0.) go to 250 |
1530 |
if(p(i).lt.AR(i) .and. p(i).lt.AL(i)) then |
1531 |
AR(i) = p(i) |
1532 |
AL(i) = p(i) |
1533 |
A6(i) = 0. |
1534 |
elseif(AR(i) .gt. AL(i)) then |
1535 |
A6(i) = 3.*(AL(i)-p(i)) |
1536 |
AR(i) = AL(i) - A6(i) |
1537 |
else |
1538 |
A6(i) = 3.*(AR(i)-p(i)) |
1539 |
AL(i) = AR(i) - A6(i) |
1540 |
endif |
1541 |
250 continue |
1542 |
endif |
1543 |
return |
1544 |
end |
1545 |
C |
1546 |
subroutine A2C(U,V,IMR,JMR,j1,j2,CRX,CRY,dtdx5,DTDY5) |
1547 |
dimension U(IMR,*),V(IMR,*),CRX(IMR,*),CRY(IMR,*),DTDX5(*) |
1548 |
C |
1549 |
do 35 j=j1,j2 |
1550 |
do 35 i=2,IMR |
1551 |
35 CRX(i,J) = dtdx5(j)*(U(i,j)+U(i-1,j)) |
1552 |
C |
1553 |
do 45 j=j1,j2 |
1554 |
45 CRX(1,J) = dtdx5(j)*(U(1,j)+U(IMR,j)) |
1555 |
C |
1556 |
do 55 i=1,IMR*JMR |
1557 |
55 CRY(i,2) = DTDY5*(V(i,2)+V(i,1)) |
1558 |
return |
1559 |
end |
1560 |
C |
1561 |
subroutine cosa(cosp,cose,JNP,PI,DP) |
1562 |
dimension cosp(*),cose(*) |
1563 |
JMR = JNP-1 |
1564 |
do 55 j=2,JNP |
1565 |
ph5 = -0.5*PI + (FLOAT(J-1)-0.5)*DP |
1566 |
55 cose(j) = cos(ph5) |
1567 |
C |
1568 |
JEQ = (JNP+1) / 2 |
1569 |
if(JMR .eq. 2*(JMR/2) ) then |
1570 |
do j=JNP, JEQ+1, -1 |
1571 |
cose(j) = cose(JNP+2-j) |
1572 |
enddo |
1573 |
else |
1574 |
C cell edge at equator. |
1575 |
cose(JEQ+1) = 1. |
1576 |
do j=JNP, JEQ+2, -1 |
1577 |
cose(j) = cose(JNP+2-j) |
1578 |
enddo |
1579 |
endif |
1580 |
C |
1581 |
do 66 j=2,JMR |
1582 |
66 cosp(j) = 0.5*(cose(j)+cose(j+1)) |
1583 |
cosp(1) = 0. |
1584 |
cosp(JNP) = 0. |
1585 |
return |
1586 |
end |
1587 |
C |
1588 |
subroutine cosc(cosp,cose,JNP,PI,DP) |
1589 |
dimension cosp(*),cose(*) |
1590 |
C |
1591 |
phi = -0.5*PI |
1592 |
do 55 j=2,JNP-1 |
1593 |
phi = phi + DP |
1594 |
55 cosp(j) = cos(phi) |
1595 |
cosp( 1) = 0. |
1596 |
cosp(JNP) = 0. |
1597 |
C |
1598 |
do 66 j=2,JNP |
1599 |
cose(j) = 0.5*(cosp(j)+cosp(j-1)) |
1600 |
66 CONTINUE |
1601 |
C |
1602 |
do 77 j=2,JNP-1 |
1603 |
cosp(j) = 0.5*(cose(j)+cose(j+1)) |
1604 |
77 CONTINUE |
1605 |
return |
1606 |
end |
1607 |
C |
1608 |
SUBROUTINE qckxyz (Q,qtmp,IMR,JNP,NLAY,j1,j2,cosp,acosp, |
1609 |
& cross,IC,NSTEP) |
1610 |
C |
1611 |
parameter( tiny = 1.E-60 ) |
1612 |
DIMENSION Q(IMR,JNP,NLAY),qtmp(IMR,JNP),cosp(*),acosp(*) |
1613 |
logical cross |
1614 |
C |
1615 |
NLAYM1 = NLAY-1 |
1616 |
len = IMR*(j2-j1+1) |
1617 |
ip = 0 |
1618 |
C |
1619 |
C Top layer |
1620 |
L = 1 |
1621 |
icr = 1 |
1622 |
call filns(q(1,1,L),IMR,JNP,j1,j2,cosp,acosp,ipy,tiny) |
1623 |
if(ipy.eq.0) goto 50 |
1624 |
call filew(q(1,1,L),qtmp,IMR,JNP,j1,j2,ipx,tiny) |
1625 |
if(ipx.eq.0) goto 50 |
1626 |
C |
1627 |
if(cross) then |
1628 |
call filcr(q(1,1,L),IMR,JNP,j1,j2,cosp,acosp,icr,tiny) |
1629 |
endif |
1630 |
if(icr.eq.0) goto 50 |
1631 |
C |
1632 |
C Vertical filling... |
1633 |
do i=1,len |
1634 |
IF( Q(i,j1,1).LT.0.) THEN |
1635 |
ip = ip + 1 |
1636 |
Q(i,j1,2) = Q(i,j1,2) + Q(i,j1,1) |
1637 |
Q(i,j1,1) = 0. |
1638 |
endif |
1639 |
enddo |
1640 |
C |
1641 |
50 continue |
1642 |
DO 225 L = 2,NLAYM1 |
1643 |
icr = 1 |
1644 |
C |
1645 |
call filns(q(1,1,L),IMR,JNP,j1,j2,cosp,acosp,ipy,tiny) |
1646 |
if(ipy.eq.0) goto 225 |
1647 |
call filew(q(1,1,L),qtmp,IMR,JNP,j1,j2,ipx,tiny) |
1648 |
if(ipx.eq.0) go to 225 |
1649 |
if(cross) then |
1650 |
call filcr(q(1,1,L),IMR,JNP,j1,j2,cosp,acosp,icr,tiny) |
1651 |
endif |
1652 |
if(icr.eq.0) goto 225 |
1653 |
C |
1654 |
do i=1,len |
1655 |
IF( Q(I,j1,L).LT.0.) THEN |
1656 |
C |
1657 |
ip = ip + 1 |
1658 |
C From above |
1659 |
qup = Q(I,j1,L-1) |
1660 |
qly = -Q(I,j1,L) |
1661 |
dup = min(qly,qup) |
1662 |
Q(I,j1,L-1) = qup - dup |
1663 |
Q(I,j1,L ) = dup-qly |
1664 |
C Below |
1665 |
Q(I,j1,L+1) = Q(I,j1,L+1) + Q(I,j1,L) |
1666 |
Q(I,j1,L) = 0. |
1667 |
ENDIF |
1668 |
ENDDO |
1669 |
225 CONTINUE |
1670 |
C |
1671 |
C BOTTOM LAYER |
1672 |
sum = 0. |
1673 |
L = NLAY |
1674 |
C |
1675 |
call filns(q(1,1,L),IMR,JNP,j1,j2,cosp,acosp,ipy,tiny) |
1676 |
if(ipy.eq.0) goto 911 |
1677 |
call filew(q(1,1,L),qtmp,IMR,JNP,j1,j2,ipx,tiny) |
1678 |
if(ipx.eq.0) goto 911 |
1679 |
C |
1680 |
call filcr(q(1,1,L),IMR,JNP,j1,j2,cosp,acosp,icr,tiny) |
1681 |
if(icr.eq.0) goto 911 |
1682 |
C |
1683 |
DO I=1,len |
1684 |
IF( Q(I,j1,L).LT.0.) THEN |
1685 |
ip = ip + 1 |
1686 |
c |
1687 |
C From above |
1688 |
C |
1689 |
qup = Q(I,j1,NLAYM1) |
1690 |
qly = -Q(I,j1,L) |
1691 |
dup = min(qly,qup) |
1692 |
Q(I,j1,NLAYM1) = qup - dup |
1693 |
C From "below" the surface. |
1694 |
sum = sum + qly-dup |
1695 |
Q(I,j1,L) = 0. |
1696 |
ENDIF |
1697 |
ENDDO |
1698 |
C |
1699 |
911 continue |
1700 |
C |
1701 |
if(ip.gt.IMR) then |
1702 |
write(6,*) 'IC=',IC,' STEP=',NSTEP, |
1703 |
& ' Vertical filling pts=',ip |
1704 |
endif |
1705 |
C |
1706 |
if(sum.gt.1.e-25) then |
1707 |
write(6,*) IC,NSTEP,' Mass source from the ground=',sum |
1708 |
endif |
1709 |
RETURN |
1710 |
END |
1711 |
C |
1712 |
subroutine filcr(q,IMR,JNP,j1,j2,cosp,acosp,icr,tiny) |
1713 |
dimension q(IMR,*),cosp(*),acosp(*) |
1714 |
icr = 0 |
1715 |
do 65 j=j1+1,j2-1 |
1716 |
DO 50 i=1,IMR-1 |
1717 |
IF(q(i,j).LT.0.) THEN |
1718 |
icr = 1 |
1719 |
dq = - q(i,j)*cosp(j) |
1720 |
C N-E |
1721 |
dn = q(i+1,j+1)*cosp(j+1) |
1722 |
d0 = max(0.,dn) |
1723 |
d1 = min(dq,d0) |
1724 |
q(i+1,j+1) = (dn - d1)*acosp(j+1) |
1725 |
dq = dq - d1 |
1726 |
C S-E |
1727 |
ds = q(i+1,j-1)*cosp(j-1) |
1728 |
d0 = max(0.,ds) |
1729 |
d2 = min(dq,d0) |
1730 |
q(i+1,j-1) = (ds - d2)*acosp(j-1) |
1731 |
q(i,j) = (d2 - dq)*acosp(j) + tiny |
1732 |
endif |
1733 |
50 continue |
1734 |
if(icr.eq.0 .and. q(IMR,j).ge.0.) goto 65 |
1735 |
DO 55 i=2,IMR |
1736 |
IF(q(i,j).LT.0.) THEN |
1737 |
icr = 1 |
1738 |
dq = - q(i,j)*cosp(j) |
1739 |
C N-W |
1740 |
dn = q(i-1,j+1)*cosp(j+1) |
1741 |
d0 = max(0.,dn) |
1742 |
d1 = min(dq,d0) |
1743 |
q(i-1,j+1) = (dn - d1)*acosp(j+1) |
1744 |
dq = dq - d1 |
1745 |
C S-W |
1746 |
ds = q(i-1,j-1)*cosp(j-1) |
1747 |
d0 = max(0.,ds) |
1748 |
d2 = min(dq,d0) |
1749 |
q(i-1,j-1) = (ds - d2)*acosp(j-1) |
1750 |
q(i,j) = (d2 - dq)*acosp(j) + tiny |
1751 |
endif |
1752 |
55 continue |
1753 |
C ***************************************** |
1754 |
C i=1 |
1755 |
i=1 |
1756 |
IF(q(i,j).LT.0.) THEN |
1757 |
icr = 1 |
1758 |
dq = - q(i,j)*cosp(j) |
1759 |
C N-W |
1760 |
dn = q(IMR,j+1)*cosp(j+1) |
1761 |
d0 = max(0.,dn) |
1762 |
d1 = min(dq,d0) |
1763 |
q(IMR,j+1) = (dn - d1)*acosp(j+1) |
1764 |
dq = dq - d1 |
1765 |
C S-W |
1766 |
ds = q(IMR,j-1)*cosp(j-1) |
1767 |
d0 = max(0.,ds) |
1768 |
d2 = min(dq,d0) |
1769 |
q(IMR,j-1) = (ds - d2)*acosp(j-1) |
1770 |
q(i,j) = (d2 - dq)*acosp(j) + tiny |
1771 |
endif |
1772 |
C ***************************************** |
1773 |
C i=IMR |
1774 |
i=IMR |
1775 |
IF(q(i,j).LT.0.) THEN |
1776 |
icr = 1 |
1777 |
dq = - q(i,j)*cosp(j) |
1778 |
C N-E |
1779 |
dn = q(1,j+1)*cosp(j+1) |
1780 |
d0 = max(0.,dn) |
1781 |
d1 = min(dq,d0) |
1782 |
q(1,j+1) = (dn - d1)*acosp(j+1) |
1783 |
dq = dq - d1 |
1784 |
C S-E |
1785 |
ds = q(1,j-1)*cosp(j-1) |
1786 |
d0 = max(0.,ds) |
1787 |
d2 = min(dq,d0) |
1788 |
q(1,j-1) = (ds - d2)*acosp(j-1) |
1789 |
q(i,j) = (d2 - dq)*acosp(j) + tiny |
1790 |
endif |
1791 |
C ***************************************** |
1792 |
65 continue |
1793 |
C |
1794 |
do i=1,IMR |
1795 |
if(q(i,j1).lt.0. .or. q(i,j2).lt.0.) then |
1796 |
icr = 1 |
1797 |
goto 80 |
1798 |
endif |
1799 |
enddo |
1800 |
C |
1801 |
80 continue |
1802 |
C |
1803 |
if(q(1,1).lt.0. .or. q(1,jnp).lt.0.) then |
1804 |
icr = 1 |
1805 |
endif |
1806 |
C |
1807 |
return |
1808 |
end |
1809 |
C |
1810 |
subroutine filns(q,IMR,JNP,j1,j2,cosp,acosp,ipy,tiny) |
1811 |
dimension q(IMR,*),cosp(*),acosp(*) |
1812 |
c logical first |
1813 |
c data first /.true./ |
1814 |
c save cap1 |
1815 |
C |
1816 |
c if(first) then |
1817 |
DP = 4.*ATAN(1.)/float(JNP-1) |
1818 |
CAP1 = IMR*(1.-COS((j1-1.5)*DP))/DP |
1819 |
c first = .false. |
1820 |
c endif |
1821 |
C |
1822 |
ipy = 0 |
1823 |
do 55 j=j1+1,j2-1 |
1824 |
DO 55 i=1,IMR |
1825 |
IF(q(i,j).LT.0.) THEN |
1826 |
ipy = 1 |
1827 |
dq = - q(i,j)*cosp(j) |
1828 |
C North |
1829 |
dn = q(i,j+1)*cosp(j+1) |
1830 |
d0 = max(0.,dn) |
1831 |
d1 = min(dq,d0) |
1832 |
q(i,j+1) = (dn - d1)*acosp(j+1) |
1833 |
dq = dq - d1 |
1834 |
C South |
1835 |
ds = q(i,j-1)*cosp(j-1) |
1836 |
d0 = max(0.,ds) |
1837 |
d2 = min(dq,d0) |
1838 |
q(i,j-1) = (ds - d2)*acosp(j-1) |
1839 |
q(i,j) = (d2 - dq)*acosp(j) + tiny |
1840 |
endif |
1841 |
55 continue |
1842 |
C |
1843 |
do i=1,imr |
1844 |
IF(q(i,j1).LT.0.) THEN |
1845 |
ipy = 1 |
1846 |
dq = - q(i,j1)*cosp(j1) |
1847 |
C North |
1848 |
dn = q(i,j1+1)*cosp(j1+1) |
1849 |
d0 = max(0.,dn) |
1850 |
d1 = min(dq,d0) |
1851 |
q(i,j1+1) = (dn - d1)*acosp(j1+1) |
1852 |
q(i,j1) = (d1 - dq)*acosp(j1) + tiny |
1853 |
endif |
1854 |
enddo |
1855 |
C |
1856 |
j = j2 |
1857 |
do i=1,imr |
1858 |
IF(q(i,j).LT.0.) THEN |
1859 |
ipy = 1 |
1860 |
dq = - q(i,j)*cosp(j) |
1861 |
C South |
1862 |
ds = q(i,j-1)*cosp(j-1) |
1863 |
d0 = max(0.,ds) |
1864 |
d2 = min(dq,d0) |
1865 |
q(i,j-1) = (ds - d2)*acosp(j-1) |
1866 |
q(i,j) = (d2 - dq)*acosp(j) + tiny |
1867 |
endif |
1868 |
enddo |
1869 |
C |
1870 |
C Check Poles. |
1871 |
if(q(1,1).lt.0.) then |
1872 |
dq = q(1,1)*cap1/float(IMR)*acosp(j1) |
1873 |
do i=1,imr |
1874 |
q(i,1) = 0. |
1875 |
q(i,j1) = q(i,j1) + dq |
1876 |
if(q(i,j1).lt.0.) ipy = 1 |
1877 |
enddo |
1878 |
endif |
1879 |
C |
1880 |
if(q(1,JNP).lt.0.) then |
1881 |
dq = q(1,JNP)*cap1/float(IMR)*acosp(j2) |
1882 |
do i=1,imr |
1883 |
q(i,JNP) = 0. |
1884 |
q(i,j2) = q(i,j2) + dq |
1885 |
if(q(i,j2).lt.0.) ipy = 1 |
1886 |
enddo |
1887 |
endif |
1888 |
C |
1889 |
return |
1890 |
end |
1891 |
C |
1892 |
subroutine filew(q,qtmp,IMR,JNP,j1,j2,ipx,tiny) |
1893 |
dimension q(IMR,*),qtmp(JNP,IMR) |
1894 |
C |
1895 |
ipx = 0 |
1896 |
C Copy & swap direction for vectorization. |
1897 |
do 25 i=1,imr |
1898 |
do 25 j=j1,j2 |
1899 |
25 qtmp(j,i) = q(i,j) |
1900 |
C |
1901 |
do 55 i=2,imr-1 |
1902 |
do 55 j=j1,j2 |
1903 |
if(qtmp(j,i).lt.0.) then |
1904 |
ipx = 1 |
1905 |
c west |
1906 |
d0 = max(0.,qtmp(j,i-1)) |
1907 |
d1 = min(-qtmp(j,i),d0) |
1908 |
qtmp(j,i-1) = qtmp(j,i-1) - d1 |
1909 |
qtmp(j,i) = qtmp(j,i) + d1 |
1910 |
c east |
1911 |
d0 = max(0.,qtmp(j,i+1)) |
1912 |
d2 = min(-qtmp(j,i),d0) |
1913 |
qtmp(j,i+1) = qtmp(j,i+1) - d2 |
1914 |
qtmp(j,i) = qtmp(j,i) + d2 + tiny |
1915 |
endif |
1916 |
55 continue |
1917 |
c |
1918 |
i=1 |
1919 |
do 65 j=j1,j2 |
1920 |
if(qtmp(j,i).lt.0.) then |
1921 |
ipx = 1 |
1922 |
c west |
1923 |
d0 = max(0.,qtmp(j,imr)) |
1924 |
d1 = min(-qtmp(j,i),d0) |
1925 |
qtmp(j,imr) = qtmp(j,imr) - d1 |
1926 |
qtmp(j,i) = qtmp(j,i) + d1 |
1927 |
c east |
1928 |
d0 = max(0.,qtmp(j,i+1)) |
1929 |
d2 = min(-qtmp(j,i),d0) |
1930 |
qtmp(j,i+1) = qtmp(j,i+1) - d2 |
1931 |
c |
1932 |
qtmp(j,i) = qtmp(j,i) + d2 + tiny |
1933 |
endif |
1934 |
65 continue |
1935 |
i=IMR |
1936 |
do 75 j=j1,j2 |
1937 |
if(qtmp(j,i).lt.0.) then |
1938 |
ipx = 1 |
1939 |
c west |
1940 |
d0 = max(0.,qtmp(j,i-1)) |
1941 |
d1 = min(-qtmp(j,i),d0) |
1942 |
qtmp(j,i-1) = qtmp(j,i-1) - d1 |
1943 |
qtmp(j,i) = qtmp(j,i) + d1 |
1944 |
c east |
1945 |
d0 = max(0.,qtmp(j,1)) |
1946 |
d2 = min(-qtmp(j,i),d0) |
1947 |
qtmp(j,1) = qtmp(j,1) - d2 |
1948 |
c |
1949 |
qtmp(j,i) = qtmp(j,i) + d2 + tiny |
1950 |
endif |
1951 |
75 continue |
1952 |
C |
1953 |
if(ipx.ne.0) then |
1954 |
do 85 j=j1,j2 |
1955 |
do 85 i=1,imr |
1956 |
85 q(i,j) = qtmp(j,i) |
1957 |
else |
1958 |
C |
1959 |
C Poles. |
1960 |
if(q(1,1).lt.0. or. q(1,JNP).lt.0.) ipx = 1 |
1961 |
endif |
1962 |
return |
1963 |
end |