1 |
! |
2 |
! $Header: /home/cvsroot/LMDZ4/libf/phylmd/convect3.F,v 1.1.1.1 2004/05/19 12:53:09 lmdzadmin Exp $ |
3 |
! |
4 |
SUBROUTINE CONVECT3 |
5 |
* (DTIME,EPMAX,ok_adj, |
6 |
* T1, R1, RS, U, V, TRA, P, PH, |
7 |
* ND, NDP1, NL, NTRA, DELT, IFLAG, |
8 |
* FT, FR, FU, FV, FTRA, PRECIP, |
9 |
* icb,inb, upwd,dnwd,dnwd0,SIG, W0,Mike,Mke, |
10 |
* Ma,MENTS,QENTS,TPS,TLS,SIGIJ,CAPE,TVP,PBASE,BUOYBASE, |
11 |
cccc * DTVPDT1,DTVPDQ1,DPLCLDT,DPLCLDR) |
12 |
* DTVPDT1,DTVPDQ1,DPLCLDT,DPLCLDR, ! sbl |
13 |
* FT2,FR2,FU2,FV2,WD,QCOND,QCONDC) ! sbl |
14 |
C |
15 |
C *** THE PARAMETER NA SHOULD IN GENERAL EQUAL ND *** |
16 |
C |
17 |
|
18 |
cFleur Introduction des traceurs dans convect3 le 6 juin 200 |
19 |
|
20 |
use dimens_m |
21 |
use dimphy |
22 |
use YOMCST |
23 |
PARAMETER (NA=60) |
24 |
|
25 |
integer NTRAC |
26 |
PARAMETER (NTRAC=nqmx-2) |
27 |
REAL DELTAC ! cld |
28 |
PARAMETER (DELTAC=0.01) ! cld |
29 |
|
30 |
INTEGER NENT(NA) |
31 |
REAL T1(ND),R1(ND),RS(ND),U(ND),V(ND),TRA(ND,NTRA) |
32 |
REAL P(ND),PH(NDP1) |
33 |
REAL FT(ND),FR(ND),FU(ND),FV(ND),FTRA(ND,NTRA) |
34 |
REAL SIG(ND),W0(ND) |
35 |
REAL UENT(NA,NA),VENT(NA,NA),TRAENT(NA,NA,NTRAC),TRATM(NA) |
36 |
REAL UP(NA),VP(NA),TRAP(NA,NTRAC) |
37 |
REAL M(NA),MP(NA),MENT(NA,NA),QENT(NA,NA),ELIJ(NA,NA) |
38 |
REAL SIJ(NA,NA),TVP(NA),TV(NA),WATER(NA) |
39 |
REAL RP(NA),EP(NA),TH(NA),WT(NA),EVAP(NA),CLW(NA) |
40 |
REAL SIGP(NA),B(NA),TP(NA),CPN(NA) |
41 |
REAL LV(NA),LVCP(NA),H(NA),HP(NA),GZ(NA) |
42 |
REAL T(NA),RR(NA) |
43 |
C |
44 |
REAL FT2(ND),FR2(ND),FU2(ND),FV2(ND) ! added sbl |
45 |
REAL U1(ND),V1(ND) ! added sbl |
46 |
C |
47 |
REAL BUOY(NA) ! Lifted parcel buoyancy |
48 |
REAL DTVPDT1(ND),DTVPDQ1(ND) ! Derivatives of parcel virtual |
49 |
C temperature wrt T1 and Q1 |
50 |
REAL CLW_NEW(NA),QI(NA) |
51 |
C |
52 |
REAL WD, BETAD ! for gust factor (sb) |
53 |
REAL QCONDC(ND) ! interface cld param (sb) |
54 |
REAL QCOND(ND),NQCOND(NA),WA(NA),MAA(NA),SIGA(NA),AXC(NA) ! cld |
55 |
C |
56 |
LOGICAL ICE_CONV,ok_adj |
57 |
PARAMETER (ICE_CONV=.TRUE.) |
58 |
|
59 |
cccccccccccccccccccccccccccccccccccccccccccccc |
60 |
c declaration des variables a sortir |
61 |
ccccccccccccccccccccccccccccccccccccccccccccc |
62 |
real Mke(nd) |
63 |
real Mike(nd) |
64 |
real Ma(nd) |
65 |
real TPS(ND) !temperature dans les ascendances non diluees |
66 |
real TLS(ND) !temperature potentielle |
67 |
real MENTS(nd,nd) |
68 |
real QENTS(nd,nd) |
69 |
REAL SIGIJ(KLEV,KLEV) |
70 |
REAL PBASE ! pressure at the cloud base level |
71 |
REAL BUOYBASE ! buoyancy at the cloud base level |
72 |
cccccccccccccccccccccccccccccccccccccccccccccc |
73 |
|
74 |
|
75 |
|
76 |
c |
77 |
real dnwd0(nd) ! precipitation driven unsaturated downdraft flux |
78 |
real dnwd(nd), dn1 ! in-cloud saturated downdraft mass flux |
79 |
real upwd(nd), up1 ! in-cloud saturated updraft mass flux |
80 |
C |
81 |
C *** ASSIGN VALUES OF THERMODYNAMIC CONSTANTS *** |
82 |
C *** THESE SHOULD BE CONSISTENT WITH *** |
83 |
C *** THOSE USED IN CALLING PROGRAM *** |
84 |
C *** NOTE: THESE ARE ALSO SPECIFIED IN SUBROUTINE TLIFT *** |
85 |
C |
86 |
c sb CPD=1005.7 |
87 |
c sb CPV=1870.0 |
88 |
c sb CL=4190.0 |
89 |
c sb CPVMCL=CL-CPV |
90 |
c sb RV=461.5 |
91 |
c sb RD=287.04 |
92 |
c sb EPS=RD/RV |
93 |
c sb ALV0=2.501E6 |
94 |
ccccccccccccccccccccccc |
95 |
c constantes coherentes avec le modele du Centre Europeen |
96 |
c sb RD = 1000.0 * 1.380658E-23 * 6.0221367E+23 / 28.9644 |
97 |
c sb RV = 1000.0 * 1.380658E-23 * 6.0221367E+23 / 18.0153 |
98 |
c sb CPD = 3.5 * RD |
99 |
c sb CPV = 4.0 * RV |
100 |
c sb CL = 4218.0 |
101 |
c sb CPVMCL=CL-CPV |
102 |
c sb EPS=RD/RV |
103 |
c sb ALV0=2.5008E+06 |
104 |
cccccccccccccccccccccc |
105 |
c on utilise les constantes thermo du Centre Europeen: (SB) |
106 |
c |
107 |
c |
108 |
CPD = RCPD |
109 |
CPV = RCPV |
110 |
CL = RCW |
111 |
CPVMCL = CL-CPV |
112 |
EPS = RD/RV |
113 |
ALV0 = RLVTT |
114 |
c |
115 |
NK = 1 ! origin level of the lifted parcel |
116 |
c |
117 |
cccccccccccccccccccccc |
118 |
C |
119 |
C *** INITIALIZE OUTPUT ARRAYS AND PARAMETERS *** |
120 |
C |
121 |
DO 5 I=1,ND |
122 |
FT(I)=0.0 |
123 |
FR(I)=0.0 |
124 |
FU(I)=0.0 |
125 |
FV(I)=0.0 |
126 |
|
127 |
FT2(I)=0.0 |
128 |
FR2(I)=0.0 |
129 |
FU2(I)=0.0 |
130 |
FV2(I)=0.0 |
131 |
|
132 |
DO 4 J=1,NTRA |
133 |
FTRA(I,J)=0.0 |
134 |
4 CONTINUE |
135 |
|
136 |
QCONDC(I)=0.0 ! cld |
137 |
QCOND(I)=0.0 ! cld |
138 |
NQCOND(I)=0.0 ! cld |
139 |
|
140 |
T(I)=T1(I) |
141 |
RR(I)=R1(I) |
142 |
U1(I)=U(I) ! added sbl |
143 |
V1(I)=V(I) ! added sbl |
144 |
5 CONTINUE |
145 |
DO 7 I=1,NL |
146 |
RDCP=(RD*(1.-RR(I))+RR(I)*RV)/ (CPD*(1.-RR(I))+RR(I)*CPV) |
147 |
TH(I)=T(I)*(1000.0/P(I))**RDCP |
148 |
7 CONTINUE |
149 |
C |
150 |
************************************************************* |
151 |
** CALCUL DES TEMPERATURES POTENTIELLES A SORTIR |
152 |
************************************************************* |
153 |
do i=1,ND |
154 |
RDCP=(RD*(1.-RR(I))+RR(I)*RV)/ (CPD*(1.-RR(I))+RR(I)*CPV) |
155 |
|
156 |
TLS(i)=T(I)*(1000.0/P(I))**RDCP |
157 |
enddo |
158 |
|
159 |
|
160 |
|
161 |
|
162 |
************************************************************ |
163 |
|
164 |
|
165 |
PRECIP=0.0 |
166 |
WD=0.0 ! sb |
167 |
IFLAG=1 |
168 |
C |
169 |
C *** SPECIFY PARAMETERS *** |
170 |
C *** PBCRIT IS THE CRITICAL CLOUD DEPTH (MB) BENEATH WHICH THE *** |
171 |
C *** PRECIPITATION EFFICIENCY IS ASSUMED TO BE ZERO *** |
172 |
C *** PTCRIT IS THE CLOUD DEPTH (MB) ABOVE WHICH THE PRECIP. *** |
173 |
C *** EFFICIENCY IS ASSUMED TO BE UNITY *** |
174 |
C *** SIGD IS THE FRACTIONAL AREA COVERED BY UNSATURATED DNDRAFT *** |
175 |
C *** SPFAC IS THE FRACTION OF PRECIPITATION FALLING OUTSIDE *** |
176 |
C *** OF CLOUD *** |
177 |
C *** ALPHA AND BETA ARE PARAMETERS THAT CONTROL THE RATE OF *** |
178 |
C *** APPROACH TO QUASI-EQUILIBRIUM *** |
179 |
C *** (THEIR STANDARD VALUES ARE 1.0 AND 0.96, RESPECTIVELY) *** |
180 |
C *** (BETA MUST BE LESS THAN OR EQUAL TO 1) *** |
181 |
C *** DTCRIT IS THE CRITICAL BUOYANCY (K) USED TO ADJUST THE *** |
182 |
C *** APPROACH TO QUASI-EQUILIBRIUM *** |
183 |
C *** IT MUST BE LESS THAN 0 *** |
184 |
C |
185 |
PBCRIT=150.0 |
186 |
PTCRIT=500.0 |
187 |
SIGD=0.01 |
188 |
SPFAC=0.15 |
189 |
c sb: |
190 |
c EPMAX=0.993 ! precip efficiency less than unity |
191 |
c EPMAX=1. ! precip efficiency less than unity |
192 |
C |
193 |
Cjyg |
194 |
CCC BETA=0.96 |
195 |
C Beta is now expressed as a function of the characteristic time |
196 |
C of the convective process. |
197 |
CCC Old value : TAU = 15000. !(for dtime = 600.s) |
198 |
CCC Other value (inducing little change) :TAU = 8000. |
199 |
TAU = 8000. |
200 |
BETA = 1.-DTIME/TAU |
201 |
Cjyg |
202 |
CCC ALPHA=1.0 |
203 |
ALPHA=1.5E-3*DTIME/TAU |
204 |
C Increase alpha in order to compensate W decrease |
205 |
ALPHA = ALPHA*1.5 |
206 |
C |
207 |
Cjyg (voir CONVECT 3) |
208 |
CCC DTCRIT=-0.2 |
209 |
DTCRIT=-2. |
210 |
Cgf&jyg |
211 |
CCC DT pour l'overshoot. |
212 |
DTOVSH = -0.2 |
213 |
|
214 |
C |
215 |
C *** INCREMENT THE COUNTER *** |
216 |
C |
217 |
SIG(ND)=SIG(ND)+1.0 |
218 |
SIG(ND)=AMIN1(SIG(ND),12.1) |
219 |
C |
220 |
C *** IF NOPT IS AN INTEGER OTHER THAN 0, CONVECT *** |
221 |
C *** RETURNS ARRAYS T AND R THAT MAY HAVE BEEN *** |
222 |
C *** ALTERED BY DRY ADIABATIC ADJUSTMENT; OTHERWISE *** |
223 |
C *** THE RETURNED ARRAYS ARE UNALTERED. *** |
224 |
C |
225 |
NOPT=0 |
226 |
c! NOPT=1 ! sbl |
227 |
C |
228 |
C *** PERFORM DRY ADIABATIC ADJUSTMENT *** |
229 |
C |
230 |
C *** DO NOT BYPASS THIS EVEN IF THE CALLING PROGRAM HAS A *** |
231 |
C *** BOUNDARY LAYER SCHEME !!! *** |
232 |
C |
233 |
IF (ok_adj) THEN ! added sbl |
234 |
|
235 |
DO 30 I=NL-1,1,-1 |
236 |
JN=0 |
237 |
DO 10 J=I+1,NL |
238 |
10 IF(TH(J).LT.TH(I))JN=J |
239 |
IF(JN.EQ.0)GOTO 30 |
240 |
AHM=0.0 |
241 |
RM=0.0 |
242 |
UM=0.0 |
243 |
VM=0.0 |
244 |
DO K=1,NTRA |
245 |
TRATM(K)=0.0 |
246 |
END DO |
247 |
DO 15 J=I,JN |
248 |
AHM=AHM+(CPD*(1.-RR(J))+RR(J)*CPV)*T(J)*(PH(J)-PH(J+1)) |
249 |
RM=RM+RR(J)*(PH(J)-PH(J+1)) |
250 |
UM=UM+U(J)*(PH(J)-PH(J+1)) |
251 |
VM=VM+V(J)*(PH(J)-PH(J+1)) |
252 |
DO K=1,NTRA |
253 |
TRATM(K)=TRATM(K)+TRA(J,K)*(PH(J)-PH(J+1)) |
254 |
END DO |
255 |
15 CONTINUE |
256 |
DPHINV=1./(PH(I)-PH(JN+1)) |
257 |
RM=RM*DPHINV |
258 |
UM=UM*DPHINV |
259 |
VM=VM*DPHINV |
260 |
DO K=1,NTRA |
261 |
TRATM(K)=TRATM(K)*DPHINV |
262 |
END DO |
263 |
A2=0.0 |
264 |
DO 20 J=I,JN |
265 |
RR(J)=RM |
266 |
U(J)=UM |
267 |
V(J)=VM |
268 |
DO K=1,NTRA |
269 |
TRA(J,K)=TRATM(K) |
270 |
END DO |
271 |
RDCP=(RD*(1.-RR(J))+RR(J)*RV)/ (CPD*(1.-RR(J))+RR(J)*CPV) |
272 |
X=(0.001*P(J))**RDCP |
273 |
T(J)=X |
274 |
A2=A2+(CPD*(1.-RR(J))+RR(J)*CPV)*X*(PH(J)-PH(J+1)) |
275 |
20 CONTINUE |
276 |
DO 25 J=I,JN |
277 |
TH(J)=AHM/A2 |
278 |
T(J)=T(J)*TH(J) |
279 |
25 CONTINUE |
280 |
30 CONTINUE |
281 |
|
282 |
ENDIF ! added sbl |
283 |
C |
284 |
C *** RESET INPUT ARRAYS IF ok_adj 0 *** |
285 |
C |
286 |
IF (ok_adj)THEN |
287 |
DO 35 I=1,ND |
288 |
|
289 |
FT2(I)=(T(I)-T1(I))/DELT ! sbl |
290 |
FR2(I)=(RR(I)-R1(I))/DELT ! sbl |
291 |
FU2(I)=(U(I)-U1(I))/DELT ! sbl |
292 |
FV2(I)=(V(I)-V1(I))/DELT ! sbl |
293 |
|
294 |
c! T1(I)=T(I) ! commente sbl |
295 |
c! R1(I)=RR(I) ! commente sbl |
296 |
35 CONTINUE |
297 |
END IF |
298 |
C |
299 |
C *** CALCULATE ARRAYS OF GEOPOTENTIAL, HEAT CAPACITY AND STATIC ENERGY |
300 |
C |
301 |
GZ(1)=0.0 |
302 |
CPN(1)=CPD*(1.-RR(1))+RR(1)*CPV |
303 |
H(1)=T(1)*CPN(1) |
304 |
DO 40 I=2,NL |
305 |
TVX=T(I)*(1.+RR(I)/EPS-RR(I)) |
306 |
TVY=T(I-1)*(1.+RR(I-1)/EPS-RR(I-1)) |
307 |
GZ(I)=GZ(I-1)+0.5*RD*(TVX+TVY)*(P(I-1)-P(I))/PH(I) |
308 |
CPN(I)=CPD*(1.-RR(I))+CPV*RR(I) |
309 |
H(I)=T(I)*CPN(I)+GZ(I) |
310 |
40 CONTINUE |
311 |
C |
312 |
C *** CALCULATE LIFTED CONDENSATION LEVEL OF AIR AT LOWEST MODEL LEVEL *** |
313 |
C *** (WITHIN 0.2% OF FORMULA OF BOLTON, MON. WEA. REV.,1980) *** |
314 |
C |
315 |
IF (T(1).LT.250.0.OR.RR(1).LE.0.0)THEN |
316 |
IFLAG=0 |
317 |
c sb3d print*,'je suis passe par 366' |
318 |
RETURN |
319 |
END IF |
320 |
|
321 |
cjyg1 Utilisation de la subroutine CLIFT |
322 |
CC RH=RR(1)/RS(1) |
323 |
CC CHI=T(1)/(1669.0-122.0*RH-T(1)) |
324 |
CC PLCL=P(1)*(RH**CHI) |
325 |
CALL CLIFT(P(1),T(1),RR(1),RS(1),PLCL,DPLCLDT,DPLCLDR) |
326 |
cjyg2 |
327 |
c sb3d PRINT *,' em_plcl,p1,t1,r1,rs1,rh ' |
328 |
c sb3d $ ,PLCL,P(1),T(1),RR(1),RS(1),RH |
329 |
c |
330 |
IF (PLCL.LT.200.0.OR.PLCL.GE.2000.0)THEN |
331 |
IFLAG=2 |
332 |
RETURN |
333 |
END IF |
334 |
Cjyg1 |
335 |
C Essais de modification de ICB |
336 |
C |
337 |
C *** CALCULATE FIRST LEVEL ABOVE LCL (=ICB) *** |
338 |
C |
339 |
CC ICB=NL-1 |
340 |
CC DO 50 I=2,NL-1 |
341 |
CC IF(P(I).LT.PLCL)THEN |
342 |
CC ICB=MIN(ICB,I) ! ICB sup ou egal a 2 |
343 |
CC END IF |
344 |
CC 50 CONTINUE |
345 |
CC IF(ICB.EQ.(NL-1))THEN |
346 |
CC IFLAG=3 |
347 |
CC RETURN |
348 |
CC END IF |
349 |
C |
350 |
C *** CALCULATE LAYER CONTAINING LCL (=ICB) *** |
351 |
C |
352 |
ICB=NL-1 |
353 |
c sb DO 50 I=2,NL-1 |
354 |
DO 50 I=3,NL-1 ! modif sb pour que ICB soit sup/egal a 2 |
355 |
C la modification consiste a comparer PLCL a PH et non a P: |
356 |
C ICB est defini par : PH(ICB)<PLCL<PH(ICB-!) |
357 |
IF(PH(I).LT.PLCL)THEN |
358 |
ICB=MIN(ICB,I) |
359 |
END IF |
360 |
50 CONTINUE |
361 |
IF(ICB.EQ.(NL-1))THEN |
362 |
IFLAG=3 |
363 |
RETURN |
364 |
END IF |
365 |
ICB = ICB-1 ! ICB sup ou egal a 2 |
366 |
Cjyg2 |
367 |
C |
368 |
C |
369 |
|
370 |
C *** SUBROUTINE TLIFT CALCULATES PART OF THE LIFTED PARCEL VIRTUAL *** |
371 |
C *** TEMPERATURE, THE ACTUAL TEMPERATURE AND THE ADIABATIC *** |
372 |
C *** LIQUID WATER CONTENT *** |
373 |
C |
374 |
|
375 |
cjyg1 |
376 |
c make sure that "Cloud base" seen by TLIFT is actually the |
377 |
c fisrt level where adiabatic ascent is saturated |
378 |
IF (PLCL .GT. P(ICB)) THEN |
379 |
c sb CALL TLIFT(P,T,RR,RS,GZ,PLCL,ICB,TVP,TP,CLW,ND,NL) |
380 |
CALL TLIFT(P,T,RR,RS,GZ,PLCL,ICB,NK,TVP,TP,CLW,ND,NL |
381 |
: ,DTVPDT1,DTVPDQ1) |
382 |
ELSE |
383 |
c sb CALL TLIFT(P,T,RR,RS,GZ,PLCL,ICB+1,TVP,TP,CLW,ND,NL) |
384 |
CALL TLIFT(P,T,RR,RS,GZ,PLCL,ICB+1,NK,TVP,TP,CLW,ND,NL |
385 |
: ,DTVPDT1,DTVPDQ1) |
386 |
ENDIF |
387 |
cjyg2 |
388 |
|
389 |
****************************************************************************** |
390 |
**** SORTIE DE LA TEMPERATURE DE L ASCENDANCE NON DILUE |
391 |
****************************************************************************** |
392 |
do i=1,ND |
393 |
TPS(i)=TP(i) |
394 |
enddo |
395 |
|
396 |
|
397 |
****************************************************************************** |
398 |
|
399 |
C |
400 |
C *** SET THE PRECIPITATION EFFICIENCIES AND THE FRACTION OF *** |
401 |
C *** PRECIPITATION FALLING OUTSIDE OF CLOUD *** |
402 |
C *** THESE MAY BE FUNCTIONS OF TP(I), P(I) AND CLW(I) *** |
403 |
C |
404 |
DO 55 I=1,NL |
405 |
PDEN=PTCRIT-PBCRIT |
406 |
c |
407 |
cjyg |
408 |
ccc EP(I)=(P(ICB)-P(I)-PBCRIT)/PDEN |
409 |
c sb EP(I)=(PLCL-P(I)-PBCRIT)/PDEN |
410 |
EP(I)=(PLCL-P(I)-PBCRIT)/PDEN * EPMAX ! sb |
411 |
c |
412 |
EP(I)=AMAX1(EP(I),0.0) |
413 |
c sb EP(I)=AMIN1(EP(I),1.0) |
414 |
EP(I)=AMIN1(EP(I),EPMAX) ! sb |
415 |
SIGP(I)=SPFAC |
416 |
C |
417 |
C *** CALCULATE VIRTUAL TEMPERATURE AND LIFTED PARCEL *** |
418 |
C *** VIRTUAL TEMPERATURE *** |
419 |
C |
420 |
TV(I)=T(I)*(1.+RR(I)/EPS-RR(I)) |
421 |
Ccd1 |
422 |
C . Keep all liquid water in lifted parcel (-> adiabatic CAPE) |
423 |
C |
424 |
ccc TVP(I)=TVP(I)-TP(I)*(RR(1)-EP(I)*CLW(I)) |
425 |
c!!! sb TVP(I)=TVP(I)-TP(I)*RR(1) ! calcule dans tlift |
426 |
Ccd2 |
427 |
C |
428 |
C *** Calculate first estimate of buoyancy |
429 |
C |
430 |
BUOY(I) = TVP(I) - TV(I) |
431 |
55 CONTINUE |
432 |
C |
433 |
C *** Set Cloud Base Buoyancy at (Plcl+DPbase) level buoyancy |
434 |
C |
435 |
DPBASE = -40. !That is 400m above LCL |
436 |
PBASE = PLCL + DPBASE |
437 |
TVPBASE = TVP(ICB )*(PBASE -P(ICB+1))/(P(ICB)-P(ICB+1)) |
438 |
$ +TVP(ICB+1)*(P(ICB)-PBASE )/(P(ICB)-P(ICB+1)) |
439 |
TVBASE = TV(ICB )*(PBASE -P(ICB+1))/(P(ICB)-P(ICB+1)) |
440 |
$ +TV(ICB+1)*(P(ICB)-PBASE )/(P(ICB)-P(ICB+1)) |
441 |
C |
442 |
c test sb: |
443 |
c@ write(*,*) '++++++++++++++++++++++++++++++++++++++++' |
444 |
c@ write(*,*) 'plcl,dpbas,tvpbas,tvbas,tvp(icb),tvp(icb1) |
445 |
c@ : ,tv(icb),tv(icb1)' |
446 |
c@ write(*,*) plcl,dpbase,tvpbase,tvbase,tvp(icb) |
447 |
c@ L ,tvp(icb+1),tv(icb),tv(icb+1) |
448 |
c@ write(*,*) '++++++++++++++++++++++++++++++++++++++++' |
449 |
c fin test sb |
450 |
BUOYBASE = TVPBASE - TVBASE |
451 |
C |
452 |
CC Set buoyancy = BUOYBASE for all levels below BASE. |
453 |
CC For safety, set : BUOY(ICB) = BUOYBASE |
454 |
DO I = ICB,NL |
455 |
IF (P(I) .GE. PBASE) THEN |
456 |
BUOY(I) = BUOYBASE |
457 |
ENDIF |
458 |
ENDDO |
459 |
BUOY(ICB) = BUOYBASE |
460 |
C |
461 |
c sb3d print *,'buoybase,tvp_tv,tvpbase,tvbase,pbase,plcl' |
462 |
c sb3d $, buoybase,tvp(icb)-tv(icb),tvpbase,tvbase,pbase,plcl |
463 |
c sb3d print *,'TVP ',(tvp(i),i=1,nl) |
464 |
c sb3d print *,'TV ',(tv(i),i=1,nl) |
465 |
c sb3d print *, 'P ',(p(i),i=1,nl) |
466 |
c sb3d print *,'ICB ',icb |
467 |
c test sb: |
468 |
c@ write(*,*) '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@' |
469 |
c@ write(*,*) 'icb,icbs,inb,buoybase:' |
470 |
c@ write(*,*) icb,icb+1,inb,buoybase |
471 |
c@ write(*,*) 'k,tvp,tv,tp,buoy,ep: ' |
472 |
c@ do k=1,nl |
473 |
c@ write(*,*) k,tvp(k),tv(k),tp(k),buoy(k),ep(k) |
474 |
c@ enddo |
475 |
c@ write(*,*) '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@' |
476 |
c fin test sb |
477 |
|
478 |
|
479 |
C |
480 |
C *** MAKE SURE THAT COLUMN IS DRY ADIABATIC BETWEEN THE SURFACE *** |
481 |
C *** AND CLOUD BASE, AND THAT LIFTED AIR IS POSITIVELY BUOYANT *** |
482 |
C *** AT CLOUD BASE *** |
483 |
C *** IF NOT, RETURN TO CALLING PROGRAM AFTER RESETTING *** |
484 |
C *** SIG(I) AND W0(I) *** |
485 |
C |
486 |
Cjyg |
487 |
CCC TDIF=TVP(ICB)-TV(ICB) |
488 |
TDIF = BUOY(ICB) |
489 |
ATH1=TH(1) |
490 |
Cjyg |
491 |
CCC ATH=TH(ICB-1)-1.0 |
492 |
ATH=TH(ICB-1)-5.0 |
493 |
c ATH=0. ! ajout sb |
494 |
c IF (ICB.GT.1) ATH=TH(ICB-1)-5.0 ! modif sb |
495 |
IF(TDIF.LT.DTCRIT.OR.ATH.GT.ATH1)THEN |
496 |
DO 60 I=1,NL |
497 |
SIG(I)=BETA*SIG(I)-2.*ALPHA*TDIF*TDIF |
498 |
SIG(I)=AMAX1(SIG(I),0.0) |
499 |
W0(I)=BETA*W0(I) |
500 |
60 CONTINUE |
501 |
IFLAG=0 |
502 |
RETURN |
503 |
END IF |
504 |
C |
505 |
|
506 |
|
507 |
C *** IF THIS POINT IS REACHED, MOIST CONVECTIVE ADJUSTMENT IS NECESSARY *** |
508 |
C *** NOW INITIALIZE VARIOUS ARRAYS USED IN THE COMPUTATIONS *** |
509 |
C |
510 |
DO 70 I=1,NL |
511 |
HP(I)=H(I) |
512 |
WT(I)=0.001 |
513 |
RP(I)=RR(I) |
514 |
UP(I)=U(I) |
515 |
VP(I)=V(I) |
516 |
DO 71 J=1,NTRA |
517 |
TRAP(I,J)=TRA(I,J) |
518 |
71 CONTINUE |
519 |
NENT(I)=0 |
520 |
WATER(I)=0.0 |
521 |
EVAP(I)=0.0 |
522 |
B(I)=0.0 |
523 |
MP(I)=0.0 |
524 |
M(I)=0.0 |
525 |
LV(I)=ALV0-CPVMCL*(T(I)-273.15) |
526 |
LVCP(I)=LV(I)/CPN(I) |
527 |
DO 70 J=1,NL |
528 |
QENT(I,J)=RR(J) |
529 |
ELIJ(I,J)=0.0 |
530 |
MENT(I,J)=0.0 |
531 |
SIJ(I,J)=0.0 |
532 |
UENT(I,J)=U(J) |
533 |
VENT(I,J)=V(J) |
534 |
DO 70 K=1,NTRA |
535 |
TRAENT(I,J,K)=TRA(J,K) |
536 |
70 CONTINUE |
537 |
|
538 |
DELTI=1.0/DELT |
539 |
C |
540 |
C *** FIND THE FIRST MODEL LEVEL (INB) ABOVE THE PARCEL'S *** |
541 |
C *** LEVEL OF NEUTRAL BUOYANCY *** |
542 |
C |
543 |
INB=NL-1 |
544 |
DO 80 I=ICB,NL-1 |
545 |
Cjyg |
546 |
CCC IF((TVP(I)-TV(I)).LT.DTCRIT)THEN |
547 |
IF(BUOY(I).LT.DTOVSH)THEN |
548 |
INB=MIN(INB,I) |
549 |
END IF |
550 |
80 CONTINUE |
551 |
|
552 |
|
553 |
|
554 |
|
555 |
C |
556 |
C *** RESET SIG(I) AND W0(I) FOR I>INB AND I<ICB *** |
557 |
C |
558 |
IF(INB.LT.(NL-1))THEN |
559 |
DO 85 I=INB+1,NL-1 |
560 |
Cjyg |
561 |
CCC SIG(I)=BETA*SIG(I)-2.0E-4*ALPHA*(TV(INB)-TVP(INB))* |
562 |
CCC 1 ABS(TV(INB)-TVP(INB)) |
563 |
SIG(I)=BETA*SIG(I)+2.*ALPHA*BUOY(INB)* |
564 |
1 ABS(BUOY(INB)) |
565 |
SIG(I)=AMAX1(SIG(I),0.0) |
566 |
W0(I)=BETA*W0(I) |
567 |
85 CONTINUE |
568 |
END IF |
569 |
DO 87 I=1,ICB |
570 |
Cjyg |
571 |
CCC SIG(I)=BETA*SIG(I)-2.0E-4*ALPHA*(TV(ICB)-TVP(ICB))* |
572 |
CCC 1 (TV(ICB)-TVP(ICB)) |
573 |
SIG(I)=BETA*SIG(I)-2.*ALPHA*BUOY(ICB)*BUOY(ICB) |
574 |
SIG(I)=AMAX1(SIG(I),0.0) |
575 |
W0(I)=BETA*W0(I) |
576 |
87 CONTINUE |
577 |
C |
578 |
C *** RESET FRACTIONAL AREAS OF UPDRAFTS AND W0 AT INITIAL TIME *** |
579 |
C *** AND AFTER 10 TIME STEPS OF NO CONVECTION *** |
580 |
C |
581 |
|
582 |
IF(SIG(ND).LT.1.5.OR.SIG(ND).GT.12.0)THEN |
583 |
DO 90 I=1,NL-1 |
584 |
SIG(I)=0.0 |
585 |
W0(I)=0.0 |
586 |
90 CONTINUE |
587 |
END IF |
588 |
C |
589 |
C *** CALCULATE LIQUID WATER STATIC ENERGY OF LIFTED PARCEL *** |
590 |
C |
591 |
DO 95 I=ICB,INB |
592 |
HP(I)=H(1)+(LV(I)+(CPD-CPV)*T(I))*EP(I)*CLW(I) |
593 |
95 CONTINUE |
594 |
C |
595 |
C *** CALCULATE CONVECTIVE AVAILABLE POTENTIAL ENERGY (CAPE), *** |
596 |
C *** VERTICAL VELOCITY (W), FRACTIONAL AREA COVERED BY *** |
597 |
C *** UNDILUTE UPDRAFT (SIG), AND UPDRAFT MASS FLUX (M) *** |
598 |
C |
599 |
CAPE=0.0 |
600 |
C |
601 |
DO 98 I=ICB+1,INB |
602 |
Cjyg1 |
603 |
CCC CAPE=CAPE+RD*(TVP(I-1)-TV(I-1))*(PH(I-1)-PH(I))/P(I-1) |
604 |
CCC DCAPE=RD*BUOY(I-1)*(PH(I-1)-PH(I))/P(I-1) |
605 |
CCC DLNP=(PH(I-1)-PH(I))/P(I-1) |
606 |
C The interval on which CAPE is computed starts at PBASE : |
607 |
DELTAP = MIN(PBASE,PH(I-1))-MIN(PBASE,PH(I)) |
608 |
CAPE=CAPE+RD*BUOY(I-1)*DELTAP/P(I-1) |
609 |
DCAPE=RD*BUOY(I-1)*DELTAP/P(I-1) |
610 |
DLNP=DELTAP/P(I-1) |
611 |
|
612 |
CAPE=AMAX1(0.0,CAPE) |
613 |
C |
614 |
SIGOLD=SIG(I) |
615 |
DTMIN=100.0 |
616 |
DO 97 J=ICB,I-1 |
617 |
Cjyg |
618 |
CCC DTMIN=AMIN1(DTMIN,(TVP(J)-TV(J))) |
619 |
DTMIN=AMIN1(DTMIN,BUOY(J)) |
620 |
97 CONTINUE |
621 |
c sb3d print *, 'DTMIN, BETA, ALPHA, SIG = ',DTMIN,BETA,ALPHA,SIG(I) |
622 |
SIG(I)=BETA*SIG(I)+ALPHA*DTMIN*ABS(DTMIN) |
623 |
SIG(I)=AMAX1(SIG(I),0.0) |
624 |
SIG(I)=AMIN1(SIG(I),0.01) |
625 |
FAC=AMIN1(((DTCRIT-DTMIN)/DTCRIT),1.0) |
626 |
Cjyg |
627 |
CC Essais de reduction de la vitesse |
628 |
CC FAC = FAC*.5 |
629 |
C |
630 |
W=(1.-BETA)*FAC*SQRT(CAPE)+BETA*W0(I) |
631 |
AMU=0.5*(SIG(I)+SIGOLD)*W |
632 |
M(I)=AMU*0.007*P(I)*(PH(I)-PH(I+1))/TV(I) |
633 |
|
634 |
W0(I)=W |
635 |
98 CONTINUE |
636 |
W0(ICB)=0.5*W0(ICB+1) |
637 |
M(ICB)=0.5*M(ICB+1)*(PH(ICB)-PH(ICB+1))/(PH(ICB+1)-PH(ICB+2)) |
638 |
SIG(ICB)=SIG(ICB+1) |
639 |
SIG(ICB-1)=SIG(ICB) |
640 |
cjyg1 |
641 |
c sb3d print *, 'Cloud base, c. top, CAPE',ICB,INB,cape |
642 |
c sb3d print *, 'SIG ',(sig(i),i=1,inb) |
643 |
c sb3d print *, 'W ',(w0(i),i=1,inb) |
644 |
c sb3d print *, 'M ',(m(i), i=1,inb) |
645 |
c sb3d print *, 'Dt1 ',(tvp(i)-tv(i),i=1,inb) |
646 |
c sb3d print *, 'Dt_vrai ',(buoy(i),i=1,inb) |
647 |
Cjyg2 |
648 |
C |
649 |
C *** CALCULATE ENTRAINED AIR MASS FLUX (MENT), TOTAL WATER MIXING *** |
650 |
C *** RATIO (QENT), TOTAL CONDENSED WATER (ELIJ), AND MIXING *** |
651 |
C *** FRACTION (SIJ) *** |
652 |
C |
653 |
|
654 |
|
655 |
DO 170 I=ICB,INB |
656 |
RTI=RR(1)-EP(I)*CLW(I) |
657 |
DO 160 J=ICB-1,INB |
658 |
BF2=1.+LV(J)*LV(J)*RS(J)/(RV*T(J)*T(J)*CPD) |
659 |
ANUM=H(J)-HP(I)+(CPV-CPD)*T(J)*(RTI-RR(J)) |
660 |
DENOM=H(I)-HP(I)+(CPD-CPV)*(RR(I)-RTI)*T(J) |
661 |
DEI=DENOM |
662 |
IF(ABS(DEI).LT.0.01)DEI=0.01 |
663 |
SIJ(I,J)=ANUM/DEI |
664 |
SIJ(I,I)=1.0 |
665 |
ALTEM=SIJ(I,J)*RR(I)+(1.-SIJ(I,J))*RTI-RS(J) |
666 |
ALTEM=ALTEM/BF2 |
667 |
CWAT=CLW(J)*(1.-EP(J)) |
668 |
STEMP=SIJ(I,J) |
669 |
IF((STEMP.LT.0.0.OR.STEMP.GT.1.0.OR. |
670 |
1 ALTEM.GT.CWAT).AND.J.GT.I)THEN |
671 |
ANUM=ANUM-LV(J)*(RTI-RS(J)-CWAT*BF2) |
672 |
DENOM=DENOM+LV(J)*(RR(I)-RTI) |
673 |
IF(ABS(DENOM).LT.0.01)DENOM=0.01 |
674 |
SIJ(I,J)=ANUM/DENOM |
675 |
ALTEM=SIJ(I,J)*RR(I)+(1.-SIJ(I,J))*RTI-RS(J) |
676 |
ALTEM=ALTEM-(BF2-1.)*CWAT |
677 |
END IF |
678 |
|
679 |
|
680 |
IF(SIJ(I,J).GT.0.0.AND.SIJ(I,J).LT.0.95)THEN |
681 |
QENT(I,J)=SIJ(I,J)*RR(I)+(1.-SIJ(I,J))*RTI |
682 |
UENT(I,J)=SIJ(I,J)*U(I)+(1.-SIJ(I,J))*U(NK) |
683 |
VENT(I,J)=SIJ(I,J)*V(I)+(1.-SIJ(I,J))*V(NK) |
684 |
DO K=1,NTRA |
685 |
TRAENT(I,J,K)=SIJ(I,J)*TRA(I,K)+(1.-SIJ(I,J))* |
686 |
1 TRA(NK,K) |
687 |
END DO |
688 |
ELIJ(I,J)=ALTEM |
689 |
ELIJ(I,J)=AMAX1(0.0,ELIJ(I,J)) |
690 |
MENT(I,J)=M(I)/(1.-SIJ(I,J)) |
691 |
NENT(I)=NENT(I)+1 |
692 |
END IF |
693 |
SIJ(I,J)=AMAX1(0.0,SIJ(I,J)) |
694 |
SIJ(I,J)=AMIN1(1.0,SIJ(I,J)) |
695 |
160 CONTINUE |
696 |
C |
697 |
C *** IF NO AIR CAN ENTRAIN AT LEVEL I ASSUME THAT UPDRAFT DETRAINS *** |
698 |
C *** AT THAT LEVEL AND CALCULATE DETRAINED AIR FLUX AND PROPERTIES *** |
699 |
C |
700 |
IF(NENT(I).EQ.0)THEN |
701 |
MENT(I,I)=M(I) |
702 |
QENT(I,I)=RR(NK)-EP(I)*CLW(I) |
703 |
UENT(I,I)=U(NK) |
704 |
VENT(I,I)=V(NK) |
705 |
DO J=1,NTRA |
706 |
TRAENT(I,I,J)=TRA(NK,J) |
707 |
END DO |
708 |
ELIJ(I,I)=CLW(I) |
709 |
SIJ(I,I)=1.0 |
710 |
END IF |
711 |
C |
712 |
DO J = ICB-1,INB |
713 |
SIGIJ(I,J) = SIJ(I,J) |
714 |
ENDDO |
715 |
C |
716 |
170 CONTINUE |
717 |
C |
718 |
C *** NORMALIZE ENTRAINED AIR MASS FLUXES TO REPRESENT EQUAL *** |
719 |
C *** PROBABILITIES OF MIXING *** |
720 |
C |
721 |
|
722 |
DO 200 I=ICB,INB |
723 |
IF(NENT(I).NE.0)THEN |
724 |
QP=RR(1)-EP(I)*CLW(I) |
725 |
ANUM=H(I)-HP(I)-LV(I)*(QP-RS(I))+(CPV-CPD)*T(I)* |
726 |
1 (QP-RR(I)) |
727 |
DENOM=H(I)-HP(I)+LV(I)*(RR(I)-QP)+ |
728 |
1 (CPD-CPV)*T(I)*(RR(I)-QP) |
729 |
IF(ABS(DENOM).LT.0.01)DENOM=0.01 |
730 |
SCRIT=ANUM/DENOM |
731 |
ALT=QP-RS(I)+SCRIT*(RR(I)-QP) |
732 |
IF(SCRIT.LE.0.0.OR.ALT.LE.0.0)SCRIT=1.0 |
733 |
SMAX=0.0 |
734 |
ASIJ=0.0 |
735 |
DO 175 J=INB,ICB-1,-1 |
736 |
IF(SIJ(I,J).GT.1.0E-16.AND.SIJ(I,J).LT.0.95)THEN |
737 |
WGH=1.0 |
738 |
IF(J.GT.I)THEN |
739 |
SJMAX=AMAX1(SIJ(I,J+1),SMAX) |
740 |
SJMAX=AMIN1(SJMAX,SCRIT) |
741 |
SMAX=AMAX1(SIJ(I,J),SMAX) |
742 |
SJMIN=AMAX1(SIJ(I,J-1),SMAX) |
743 |
SJMIN=AMIN1(SJMIN,SCRIT) |
744 |
IF(SIJ(I,J).LT.(SMAX-1.0E-16))WGH=0.0 |
745 |
SMID=AMIN1(SIJ(I,J),SCRIT) |
746 |
ELSE |
747 |
SJMAX=AMAX1(SIJ(I,J+1),SCRIT) |
748 |
SMID=AMAX1(SIJ(I,J),SCRIT) |
749 |
SJMIN=0.0 |
750 |
IF(J.GT.1)SJMIN=SIJ(I,J-1) |
751 |
SJMIN=AMAX1(SJMIN,SCRIT) |
752 |
END IF |
753 |
DELP=ABS(SJMAX-SMID) |
754 |
DELM=ABS(SJMIN-SMID) |
755 |
ASIJ=ASIJ+WGH*(DELP+DELM) |
756 |
MENT(I,J)=MENT(I,J)*(DELP+DELM)*WGH |
757 |
END IF |
758 |
175 CONTINUE |
759 |
ASIJ=AMAX1(1.0E-16,ASIJ) |
760 |
ASIJ=1.0/ASIJ |
761 |
DO 180 J=ICB-1,INB |
762 |
MENT(I,J)=MENT(I,J)*ASIJ |
763 |
180 CONTINUE |
764 |
ASUM=0.0 |
765 |
BSUM=0.0 |
766 |
DO 190 J=ICB-1,INB |
767 |
ASUM=ASUM+MENT(I,J) |
768 |
MENT(I,J)=MENT(I,J)*SIG(J) |
769 |
BSUM=BSUM+MENT(I,J) |
770 |
190 CONTINUE |
771 |
BSUM=AMAX1(BSUM,1.0E-16) |
772 |
BSUM=1.0/BSUM |
773 |
DO 195 J=ICB-1,INB |
774 |
MENT(I,J)=MENT(I,J)*ASUM*BSUM |
775 |
195 CONTINUE |
776 |
CSUM=0.0 |
777 |
DO 197 J=ICB-1,INB |
778 |
CSUM=CSUM+MENT(I,J) |
779 |
197 CONTINUE |
780 |
|
781 |
IF(CSUM.LT.M(I))THEN |
782 |
NENT(I)=0 |
783 |
MENT(I,I)=M(I) |
784 |
QENT(I,I)=RR(1)-EP(I)*CLW(I) |
785 |
UENT(I,I)=U(NK) |
786 |
VENT(I,I)=V(NK) |
787 |
DO J=1,NTRA |
788 |
TRAENT(I,I,J)=TRA(NK,J) |
789 |
END DO |
790 |
ELIJ(I,I)=CLW(I) |
791 |
SIJ(I,I)=1.0 |
792 |
END IF |
793 |
END IF |
794 |
200 CONTINUE |
795 |
|
796 |
|
797 |
|
798 |
*************************************************************** |
799 |
** CALCUL DES MENTS(I,J) ET DES QENTS(I,J) |
800 |
************************************************************** |
801 |
|
802 |
DO im=1,nd |
803 |
do jm=1,nd |
804 |
|
805 |
QENTS(im,jm)=QENT(im,jm) |
806 |
MENTS(im,jm)=MENT(im,jm) |
807 |
enddo |
808 |
enddo |
809 |
|
810 |
*********************************************************** |
811 |
c--- test sb: |
812 |
c@ write(*,*) '^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^' |
813 |
c@ write(*,*) 'inb,m(inb),ment(inb,inb),sigij(inb,inb):' |
814 |
c@ write(*,*) inb,m(inb),ment(inb,inb),sigij(inb,inb) |
815 |
c@ write(*,*) '^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^' |
816 |
c--- |
817 |
|
818 |
|
819 |
|
820 |
|
821 |
C |
822 |
C *** CHECK WHETHER EP(INB)=0, IF SO, SKIP PRECIPITATING *** |
823 |
C *** DOWNDRAFT CALCULATION *** |
824 |
C |
825 |
IF(EP(INB).LT.0.0001)GOTO 405 |
826 |
C |
827 |
C *** INTEGRATE LIQUID WATER EQUATION TO FIND CONDENSED WATER *** |
828 |
C *** AND CONDENSED WATER FLUX *** |
829 |
C |
830 |
WFLUX=0.0 |
831 |
TINV=1./3. |
832 |
C |
833 |
C *** BEGIN DOWNDRAFT LOOP *** |
834 |
C |
835 |
DO 400 I=INB,1,-1 |
836 |
C |
837 |
C *** CALCULATE DETRAINED PRECIPITATION *** |
838 |
C |
839 |
|
840 |
|
841 |
WDTRAIN=10.0*EP(I)*M(I)*CLW(I) |
842 |
IF(I.GT.1)THEN |
843 |
DO 320 J=1,I-1 |
844 |
AWAT=ELIJ(J,I)-(1.-EP(I))*CLW(I) |
845 |
AWAT=AMAX1(AWAT,0.0) |
846 |
320 WDTRAIN=WDTRAIN+10.0*AWAT*MENT(J,I) |
847 |
END IF |
848 |
C |
849 |
C *** FIND RAIN WATER AND EVAPORATION USING PROVISIONAL *** |
850 |
C *** ESTIMATES OF RP(I)AND RP(I-1) *** |
851 |
C |
852 |
|
853 |
|
854 |
WT(I)=45.0 |
855 |
IF(I.LT.INB)THEN |
856 |
RP(I)=RP(I+1)+(CPD*(T(I+1)-T(I))+GZ(I+1)-GZ(I))/LV(I) |
857 |
RP(I)=0.5*(RP(I)+RR(I)) |
858 |
END IF |
859 |
RP(I)=AMAX1(RP(I),0.0) |
860 |
RP(I)=AMIN1(RP(I),RS(I)) |
861 |
RP(INB)=RR(INB) |
862 |
IF(I.EQ.1)THEN |
863 |
AFAC=P(1)*(RS(1)-RP(1))/(1.0E4+2000.0*P(1)*RS(1)) |
864 |
ELSE |
865 |
RP(I-1)=RP(I)+(CPD*(T(I)-T(I-1))+GZ(I)-GZ(I-1))/LV(I) |
866 |
RP(I-1)=0.5*(RP(I-1)+RR(I-1)) |
867 |
RP(I-1)=AMIN1(RP(I-1),RS(I-1)) |
868 |
RP(I-1)=AMAX1(RP(I-1),0.0) |
869 |
AFAC1=P(I)*(RS(I)-RP(I))/(1.0E4+2000.0*P(I)*RS(I)) |
870 |
AFAC2=P(I-1)*(RS(I-1)-RP(I-1))/(1.0E4+ |
871 |
1 2000.0*P(I-1)*RS(I-1)) |
872 |
AFAC=0.5*(AFAC1+AFAC2) |
873 |
END IF |
874 |
IF(I.EQ.INB)AFAC=0.0 |
875 |
AFAC=AMAX1(AFAC,0.0) |
876 |
BFAC=1./(SIGD*WT(I)) |
877 |
C |
878 |
Cjyg1 |
879 |
CCC SIGT=1.0 |
880 |
CCC IF(I.GE.ICB)SIGT=SIGP(I) |
881 |
C Prise en compte de la variation progressive de SIGT dans |
882 |
C les couches ICB et ICB-1: |
883 |
C Pour PLCL<PH(I+1), PR1=0 & PR2=1 |
884 |
C Pour PLCL>PH(I), PR1=1 & PR2=0 |
885 |
C Pour PH(I+1)<PLCL<PH(I), PR1 est la proportion a cheval |
886 |
C sur le nuage, et PR2 est la proportion sous la base du |
887 |
C nuage. |
888 |
PR1 =(PLCL-PH(I+1))/(PH(I)-PH(I+1)) |
889 |
PR1 = MAX(0.,MIN(1.,PR1)) |
890 |
PR2 = (PH(I)-PLCL)/(PH(I)-PH(I+1)) |
891 |
PR2 = MAX(0.,MIN(1.,PR2)) |
892 |
SIGT = SIGP(I)*PR1 + PR2 |
893 |
c sb3d print *,'i,sigt,pr1,pr2', i,sigt,pr1,pr2 |
894 |
Cjyg2 |
895 |
C |
896 |
B6=BFAC*50.*SIGD*(PH(I)-PH(I+1))*SIGT*AFAC |
897 |
C6=WATER(I+1)+BFAC*WDTRAIN-50.*SIGD*BFAC* |
898 |
1 (PH(I)-PH(I+1))*EVAP(I+1) |
899 |
IF(C6.GT.0.0)THEN |
900 |
REVAP=0.5*(-B6+SQRT(B6*B6+4.*C6)) |
901 |
EVAP(I)=SIGT*AFAC*REVAP |
902 |
WATER(I)=REVAP*REVAP |
903 |
ELSE |
904 |
EVAP(I)=-EVAP(I+1)+0.02*(WDTRAIN+SIGD*WT(I)* |
905 |
1 WATER(I+1))/(SIGD*(PH(I)-PH(I+1))) |
906 |
END IF |
907 |
|
908 |
|
909 |
C |
910 |
C *** CALCULATE PRECIPITATING DOWNDRAFT MASS FLUX UNDER *** |
911 |
C *** HYDROSTATIC APPROXIMATION *** |
912 |
C |
913 |
IF(I.EQ.1)GOTO 360 |
914 |
TEVAP=AMAX1(0.0,EVAP(I)) |
915 |
DELTH=AMAX1(0.001,(TH(I)-TH(I-1))) |
916 |
MP(I)=10.*LVCP(I)*SIGD*TEVAP*(P(I-1)-P(I))/DELTH |
917 |
C |
918 |
C *** IF HYDROSTATIC ASSUMPTION FAILS, *** |
919 |
C *** SOLVE CUBIC DIFFERENCE EQUATION FOR DOWNDRAFT THETA *** |
920 |
C *** AND MASS FLUX FROM TWO SIMULTANEOUS DIFFERENTIAL EQNS *** |
921 |
C |
922 |
AMFAC=SIGD*SIGD*70.0*PH(I)*(P(I-1)-P(I))* |
923 |
1 (TH(I)-TH(I-1))/(TV(I)*TH(I)) |
924 |
AMP2=ABS(MP(I+1)*MP(I+1)-MP(I)*MP(I)) |
925 |
IF(AMP2.GT.(0.1*AMFAC))THEN |
926 |
XF=100.0*SIGD*SIGD*SIGD*(PH(I)-PH(I+1)) |
927 |
TF=B(I)-5.0*(TH(I)-TH(I-1))*T(I)/(LVCP(I)*SIGD*TH(I)) |
928 |
AF=XF*TF+MP(I+1)*MP(I+1)*TINV |
929 |
BF=2.*(TINV*MP(I+1))**3+TINV*MP(I+1)*XF*TF+50.* |
930 |
1 (P(I-1)-P(I))*XF*TEVAP |
931 |
FAC2=1.0 |
932 |
IF(BF.LT.0.0)FAC2=-1.0 |
933 |
BF=ABS(BF) |
934 |
UR=0.25*BF*BF-AF*AF*AF*TINV*TINV*TINV |
935 |
IF(UR.GE.0.0)THEN |
936 |
SRU=SQRT(UR) |
937 |
FAC=1.0 |
938 |
IF((0.5*BF-SRU).LT.0.0)FAC=-1.0 |
939 |
MP(I)=MP(I+1)*TINV+(0.5*BF+SRU)**TINV+ |
940 |
1 FAC*(ABS(0.5*BF-SRU))**TINV |
941 |
ELSE |
942 |
D=ATAN(2.*SQRT(-UR)/(BF+1.0E-28)) |
943 |
IF(FAC2.LT.0.0)D=3.14159-D |
944 |
MP(I)=MP(I+1)*TINV+2.*SQRT(AF*TINV)*COS(D*TINV) |
945 |
END IF |
946 |
MP(I)=AMAX1(0.0,MP(I)) |
947 |
B(I-1)=B(I)+100.0*(P(I-1)-P(I))*TEVAP/(MP(I)+SIGD*0.1)- |
948 |
1 10.0*(TH(I)-TH(I-1))*T(I)/(LVCP(I)*SIGD*TH(I)) |
949 |
B(I-1)=AMAX1(B(I-1),0.0) |
950 |
END IF |
951 |
|
952 |
|
953 |
C |
954 |
C *** LIMIT MAGNITUDE OF MP(I) TO MEET CFL CONDITION *** |
955 |
C |
956 |
AMPMAX=2.0*(PH(I)-PH(I+1))*DELTI |
957 |
AMP2=2.0*(PH(I-1)-PH(I))*DELTI |
958 |
AMPMAX=AMIN1(AMPMAX,AMP2) |
959 |
MP(I)=AMIN1(MP(I),AMPMAX) |
960 |
C |
961 |
C *** FORCE MP TO DECREASE LINEARLY TO ZERO *** |
962 |
C *** BETWEEN CLOUD BASE AND THE SURFACE *** |
963 |
C |
964 |
IF(P(I).GT.P(ICB))THEN |
965 |
MP(I)=MP(ICB)*(P(1)-P(I))/(P(1)-P(ICB)) |
966 |
END IF |
967 |
360 CONTINUE |
968 |
C |
969 |
C *** FIND MIXING RATIO OF PRECIPITATING DOWNDRAFT *** |
970 |
C |
971 |
IF(I.EQ.INB)GOTO 400 |
972 |
RP(I)=RR(I) |
973 |
IF(MP(I).GT.MP(I+1))THEN |
974 |
RP(I)=RP(I+1)*MP(I+1)+RR(I)*(MP(I)-MP(I+1))+ |
975 |
1 5.*SIGD*(PH(I)-PH(I+1))*(EVAP(I+1)+EVAP(I)) |
976 |
RP(I)=RP(I)/MP(I) |
977 |
UP(I)=UP(I+1)*MP(I+1)+U(I)*(MP(I)-MP(I+1)) |
978 |
UP(I)=UP(I)/MP(I) |
979 |
VP(I)=VP(I+1)*MP(I+1)+V(I)*(MP(I)-MP(I+1)) |
980 |
VP(I)=VP(I)/MP(I) |
981 |
DO J=1,NTRA |
982 |
TRAP(I,J)=TRAP(I+1,J)*MP(I+1)+ |
983 |
s TRAP(I,J)*(MP(I)-MP(I+1)) |
984 |
TRAP(I,J)=TRAP(I,J)/MP(I) |
985 |
END DO |
986 |
ELSE |
987 |
IF(MP(I+1).GT.1.0E-16)THEN |
988 |
RP(I)=RP(I+1)+5.0*SIGD*(PH(I)-PH(I+1))*(EVAP(I+1)+ |
989 |
1 EVAP(I))/MP(I+1) |
990 |
UP(I)=UP(I+1) |
991 |
VP(I)=VP(I+1) |
992 |
DO J=1,NTRA |
993 |
TRAP(I,J)=TRAP(I+1,J) |
994 |
END DO |
995 |
END IF |
996 |
END IF |
997 |
RP(I)=AMIN1(RP(I),RS(I)) |
998 |
RP(I)=AMAX1(RP(I),0.0) |
999 |
400 CONTINUE |
1000 |
C |
1001 |
C *** CALCULATE SURFACE PRECIPITATION IN MM/DAY *** |
1002 |
C |
1003 |
PRECIP=WT(1)*SIGD*WATER(1)*8640.0 |
1004 |
|
1005 |
c sb *** Calculate downdraft velocity scale and surface temperature and *** |
1006 |
c sb *** water vapor fluctuations *** |
1007 |
c sb (inspire de convect 4.3) |
1008 |
|
1009 |
c BETAD=10.0 |
1010 |
BETAD=5.0 |
1011 |
WD=BETAD*ABS(MP(ICB))*0.01*RD*T(ICB)/(SIGD*P(ICB)) |
1012 |
|
1013 |
405 CONTINUE |
1014 |
C |
1015 |
C *** CALCULATE TENDENCIES OF LOWEST LEVEL POTENTIAL TEMPERATURE *** |
1016 |
C *** AND MIXING RATIO *** |
1017 |
C |
1018 |
DPINV=1.0/(PH(1)-PH(2)) |
1019 |
AM=0.0 |
1020 |
DO 410 K=2,INB |
1021 |
410 AM=AM+M(K) |
1022 |
IF((0.1*DPINV*AM).GE.DELTI)IFLAG=4 |
1023 |
FT(1)=0.1*DPINV*AM*(T(2)-T(1)+(GZ(2)-GZ(1))/CPN(1)) |
1024 |
FT(1)=FT(1)-0.5*LVCP(1)*SIGD*(EVAP(1)+EVAP(2)) |
1025 |
FT(1)=FT(1)-0.09*SIGD*MP(2)*T(1)*B(1)*DPINV |
1026 |
FT(1)=FT(1)+0.01*SIGD*WT(1)*(CL-CPD)*WATER(2)*(T(2)- |
1027 |
1 T(1))*DPINV/CPN(1) |
1028 |
FR(1)=0.1*MP(2)*(RP(2)-RR(1))* |
1029 |
Ccorrection bug conservation eau |
1030 |
C 1 DPINV+SIGD*0.5*(EVAP(1)+EVAP(2)) |
1031 |
1 DPINV+SIGD*0.5*(EVAP(1)+EVAP(2)) |
1032 |
cIM cf. SBL |
1033 |
C 1 DPINV+SIGD*EVAP(1) |
1034 |
FR(1)=FR(1)+0.1*AM*(RR(2)-RR(1))*DPINV |
1035 |
FU(1)=FU(1)+0.1*DPINV*(MP(2)*(UP(2)-U(1))+AM*(U(2)-U(1))) |
1036 |
FV(1)=FV(1)+0.1*DPINV*(MP(2)*(VP(2)-V(1))+AM*(V(2)-V(1))) |
1037 |
DO J=1,NTRA |
1038 |
FTRA(1,J)=FTRA(1,J)+0.1*DPINV*(MP(2)*(TRAP(2,J)-TRA(1,J))+ |
1039 |
1 AM*(TRA(2,J)-TRA(1,J))) |
1040 |
END DO |
1041 |
AMDE=0.0 |
1042 |
DO 415 J=2,INB |
1043 |
FR(1)=FR(1)+0.1*DPINV*MENT(J,1)*(QENT(J,1)-RR(1)) |
1044 |
FU(1)=FU(1)+0.1*DPINV*MENT(J,1)*(UENT(J,1)-U(1)) |
1045 |
FV(1)=FV(1)+0.1*DPINV*MENT(J,1)*(VENT(J,1)-V(1)) |
1046 |
DO K=1,NTRA |
1047 |
FTRA(1,K)=FTRA(1,K)+0.1*DPINV*MENT(J,1)*(TRAENT(J,1,K)- |
1048 |
1 TRA(1,K)) |
1049 |
END DO |
1050 |
415 CONTINUE |
1051 |
C |
1052 |
C *** CALCULATE TENDENCIES OF POTENTIAL TEMPERATURE AND MIXING RATIO *** |
1053 |
C *** AT LEVELS ABOVE THE LOWEST LEVEL *** |
1054 |
C |
1055 |
C *** FIRST FIND THE NET SATURATED UPDRAFT AND DOWNDRAFT MASS FLUXES *** |
1056 |
C *** THROUGH EACH LEVEL *** |
1057 |
C |
1058 |
|
1059 |
|
1060 |
DO 500 I=2,INB |
1061 |
DPINV=1.0/(PH(I)-PH(I+1)) |
1062 |
CPINV=1.0/CPN(I) |
1063 |
AMP1=0.0 |
1064 |
DO 440 K=I+1,INB+1 |
1065 |
440 AMP1=AMP1+M(K) |
1066 |
DO 450 K=1,I |
1067 |
DO 450 J=I+1,INB+1 |
1068 |
AMP1=AMP1+MENT(K,J) |
1069 |
450 CONTINUE |
1070 |
IF((0.1*DPINV*AMP1).GE.DELTI)IFLAG=4 |
1071 |
AD=0.0 |
1072 |
DO 470 K=1,I-1 |
1073 |
DO 470 J=I,INB |
1074 |
470 AD=AD+MENT(J,K) |
1075 |
FT(I)=0.1*DPINV*(AMP1*(T(I+1)-T(I)+(GZ(I+1)-GZ(I))* |
1076 |
1 CPINV)-AD*(T(I)-T(I-1)+(GZ(I)-GZ(I-1))*CPINV)) |
1077 |
2 -0.5*SIGD*LVCP(I)*(EVAP(I)+EVAP(I+1)) |
1078 |
RAT=CPN(I-1)*CPINV |
1079 |
FT(I)=FT(I)-0.09*SIGD*(MP(I+1)*T(I)* |
1080 |
1 B(I)-MP(I)*T(I-1)*RAT*B(I-1))*DPINV |
1081 |
FT(I)=FT(I)+0.1*DPINV*MENT(I,I)*(HP(I)-H(I)+ |
1082 |
1 T(I)*(CPV-CPD)*(RR(I)-QENT(I,I)))*CPINV |
1083 |
FT(I)=FT(I)+0.01*SIGD*WT(I)*(CL-CPD)*WATER(I+1)* |
1084 |
1 (T(I+1)-T(I))*DPINV*CPINV |
1085 |
FR(I)=0.1*DPINV*(AMP1*(RR(I+1)-RR(I))- |
1086 |
1 AD*(RR(I)-RR(I-1))) |
1087 |
FU(I)=FU(I)+0.1*DPINV*(AMP1*(U(I+1)-U(I))- |
1088 |
1 AD*(U(I)-U(I-1))) |
1089 |
FV(I)=FV(I)+0.1*DPINV*(AMP1*(V(I+1)-V(I))- |
1090 |
1 AD*(V(I)-V(I-1))) |
1091 |
DO K=1,NTRA |
1092 |
FTRA(I,K)=FTRA(I,K)+0.1*DPINV*(AMP1*(TRA(I+1,K)- |
1093 |
1 TRA(I,K))-AD*(TRA(I,K)-TRA(I-1,K))) |
1094 |
END DO |
1095 |
DO 480 K=1,I-1 |
1096 |
AWAT=ELIJ(K,I)-(1.-EP(I))*CLW(I) |
1097 |
AWAT=AMAX1(AWAT,0.0) |
1098 |
FR(I)=FR(I)+0.1*DPINV*MENT(K,I)*(QENT(K,I)-AWAT |
1099 |
1 -RR(I)) |
1100 |
FU(I)=FU(I)+0.1*DPINV*MENT(K,I)*(UENT(K,I)-U(I)) |
1101 |
FV(I)=FV(I)+0.1*DPINV*MENT(K,I)*(VENT(K,I)-V(I)) |
1102 |
C (saturated updrafts resulting from mixing) ! cld |
1103 |
QCOND(I)=QCOND(I)+(ELIJ(K,I)-AWAT) ! cld |
1104 |
NQCOND(I)=NQCOND(I)+1. ! cld |
1105 |
DO J=1,NTRA |
1106 |
FTRA(I,J)=FTRA(I,J)+0.1*DPINV*MENT(K,I)*(TRAENT(K,I,J)- |
1107 |
1 TRA(I,J)) |
1108 |
END DO |
1109 |
480 CONTINUE |
1110 |
DO 490 K=I,INB |
1111 |
FR(I)=FR(I)+0.1*DPINV*MENT(K,I)*(QENT(K,I)-RR(I)) |
1112 |
FU(I)=FU(I)+0.1*DPINV*MENT(K,I)*(UENT(K,I)-U(I)) |
1113 |
FV(I)=FV(I)+0.1*DPINV*MENT(K,I)*(VENT(K,I)-V(I)) |
1114 |
DO J=1,NTRA |
1115 |
FTRA(I,J)=FTRA(I,J)+0.1*DPINV*MENT(K,I)*(TRAENT(K,I,J)- |
1116 |
1 TRA(I,J)) |
1117 |
END DO |
1118 |
490 CONTINUE |
1119 |
FR(I)=FR(I)+0.5*SIGD*(EVAP(I)+EVAP(I+1))+0.1*(MP(I+1)* |
1120 |
1 (RP(I+1)-RR(I))-MP(I)*(RP(I)-RR(I-1)))*DPINV |
1121 |
FU(I)=FU(I)+0.1*(MP(I+1)*(UP(I+1)-U(I))-MP(I)* |
1122 |
1 (UP(I)-U(I-1)))*DPINV |
1123 |
FV(I)=FV(I)+0.1*(MP(I+1)*(VP(I+1)-V(I))-MP(I)* |
1124 |
1 (VP(I)-V(I-1)))*DPINV |
1125 |
DO J=1,NTRA |
1126 |
FTRA(I,J)=FTRA(I,J)+0.1*DPINV*(MP(I+1)*(TRAP(I+1,J)-TRA(I,J))- |
1127 |
1 MP(I)*(TRAP(I,J)-TRAP(I-1,J))) |
1128 |
END DO |
1129 |
C (saturated downdrafts resulting from mixing) ! cld |
1130 |
DO K=I+1,INB ! cld |
1131 |
QCOND(I)=QCOND(I)+ELIJ(K,I) ! cld |
1132 |
NQCOND(I)=NQCOND(I)+1. ! cld |
1133 |
ENDDO ! cld |
1134 |
C (particular case: no detraining level is found) ! cld |
1135 |
IF (NENT(I).EQ.0) THEN ! cld |
1136 |
QCOND(I)=QCOND(I)+(1-EP(I))*CLW(I) ! cld |
1137 |
NQCOND(I)=NQCOND(I)+1. ! cld |
1138 |
ENDIF ! cld |
1139 |
IF (NQCOND(I).NE.0.) THEN ! cld |
1140 |
QCOND(I)=QCOND(I)/NQCOND(I) ! cld |
1141 |
ENDIF ! cld |
1142 |
500 CONTINUE |
1143 |
|
1144 |
|
1145 |
|
1146 |
C |
1147 |
C *** MOVE THE DETRAINMENT AT LEVEL INB DOWN TO LEVEL INB-1 *** |
1148 |
C *** IN SUCH A WAY AS TO PRESERVE THE VERTICALLY *** |
1149 |
C *** INTEGRATED ENTHALPY AND WATER TENDENCIES *** |
1150 |
C |
1151 |
c test sb: |
1152 |
c@ write(*,*) '--------------------------------------------' |
1153 |
c@ write(*,*) 'inb,ft,hp,h,t,rr,qent,ment,water,waterp,wt,mp,b' |
1154 |
c@ write(*,*) inb,ft(inb),hp(inb),h(inb) |
1155 |
c@ : ,t(inb),rr(inb),qent(inb,inb) |
1156 |
c@ : ,ment(inb,inb),water(inb) |
1157 |
c@ : ,water(inb+1),wt(inb),mp(inb),b(inb) |
1158 |
c@ write(*,*) '--------------------------------------------' |
1159 |
c fin test sb: |
1160 |
|
1161 |
AX=0.1*MENT(INB,INB)*(HP(INB)-H(INB)+T(INB)* |
1162 |
1 (CPV-CPD)*(RR(INB)-QENT(INB,INB)))/(CPN(INB)* |
1163 |
2 (PH(INB)-PH(INB+1))) |
1164 |
FT(INB)=FT(INB)-AX |
1165 |
FT(INB-1)=FT(INB-1)+AX*CPN(INB)*(PH(INB)-PH(INB+1))/ |
1166 |
1 (CPN(INB-1)*(PH(INB-1)-PH(INB))) |
1167 |
BX=0.1*MENT(INB,INB)*(QENT(INB,INB)-RR(INB))/ |
1168 |
1 (PH(INB)-PH(INB+1)) |
1169 |
FR(INB)=FR(INB)-BX |
1170 |
FR(INB-1)=FR(INB-1)+BX*(PH(INB)-PH(INB+1))/ |
1171 |
1 (PH(INB-1)-PH(INB)) |
1172 |
CX=0.1*MENT(INB,INB)*(UENT(INB,INB)-U(INB))/ |
1173 |
1 (PH(INB)-PH(INB+1)) |
1174 |
FU(INB)=FU(INB)-CX |
1175 |
FU(INB-1)=FU(INB-1)+CX*(PH(INB)-PH(INB+1))/ |
1176 |
1 (PH(INB-1)-PH(INB)) |
1177 |
DX=0.1*MENT(INB,INB)*(VENT(INB,INB)-V(INB))/ |
1178 |
1 (PH(INB)-PH(INB+1)) |
1179 |
FV(INB)=FV(INB)-DX |
1180 |
FV(INB-1)=FV(INB-1)+DX*(PH(INB)-PH(INB+1))/ |
1181 |
1 (PH(INB-1)-PH(INB)) |
1182 |
DO J=1,NTRA |
1183 |
EX=0.1*MENT(INB,INB)*(TRAENT(INB,INB,J) |
1184 |
1 -TRA(INB,J))/(PH(INB)-PH(INB+1)) |
1185 |
FTRA(INB,J)=FTRA(INB,J)-EX |
1186 |
FTRA(INB-1,J)=FTRA(INB-1,J)+EX* |
1187 |
1 (PH(INB)-PH(INB+1))/(PH(INB-1)-PH(INB)) |
1188 |
ENDDO |
1189 |
C |
1190 |
C *** HOMOGINIZE TENDENCIES BELOW CLOUD BASE *** |
1191 |
C |
1192 |
ASUM=0.0 |
1193 |
BSUM=0.0 |
1194 |
CSUM=0.0 |
1195 |
DSUM=0.0 |
1196 |
DO 650 I=1,ICB-1 |
1197 |
ASUM=ASUM+FT(I)*(PH(I)-PH(I+1)) |
1198 |
BSUM=BSUM+FR(I)*(LV(I)+(CL-CPD)*(T(I)-T(1)))* |
1199 |
1 (PH(I)-PH(I+1)) |
1200 |
CSUM=CSUM+(LV(I)+(CL-CPD)*(T(I)-T(1)))*(PH(I)-PH(I+1)) |
1201 |
DSUM=DSUM+T(I)*(PH(I)-PH(I+1))/TH(I) |
1202 |
650 CONTINUE |
1203 |
DO 700 I=1,ICB-1 |
1204 |
FT(I)=ASUM*T(I)/(TH(I)*DSUM) |
1205 |
FR(I)=BSUM/CSUM |
1206 |
700 CONTINUE |
1207 |
C |
1208 |
C *** RESET COUNTER AND RETURN *** |
1209 |
C |
1210 |
SIG(ND)=2.0 |
1211 |
c |
1212 |
c |
1213 |
do i = 1, nd |
1214 |
upwd(i) = 0.0 |
1215 |
dnwd(i) = 0.0 |
1216 |
c sb dnwd0(i) = - mp(i) |
1217 |
enddo |
1218 |
c |
1219 |
do i = 1, nl |
1220 |
dnwd0(i) = - mp(i) |
1221 |
enddo |
1222 |
do i = nl+1, nd |
1223 |
dnwd0(i) = 0. |
1224 |
enddo |
1225 |
c |
1226 |
do i = icb, inb |
1227 |
upwd(i) = 0.0 |
1228 |
dnwd(i) = 0.0 |
1229 |
|
1230 |
do k =i, inb |
1231 |
up1=0.0 |
1232 |
dn1=0.0 |
1233 |
do n = 1, i-1 |
1234 |
up1 = up1 + ment(n,k) |
1235 |
dn1 = dn1 - ment(k,n) |
1236 |
enddo |
1237 |
upwd(i) = upwd(i)+ m(k) + up1 |
1238 |
dnwd(i) = dnwd(i) + dn1 |
1239 |
enddo |
1240 |
enddo |
1241 |
|
1242 |
cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
1243 |
c DETERMINATION DE LA VARIATION DE FLUX ASCENDANT ENTRE |
1244 |
C DEUX NIVEAU NON DILUE Mike |
1245 |
cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
1246 |
|
1247 |
|
1248 |
c sb do i=1,ND |
1249 |
c sb Mike(i)=M(i) |
1250 |
c sb enddo |
1251 |
|
1252 |
do i = 1, NL |
1253 |
Mike(i) = M(i) |
1254 |
enddo |
1255 |
do i = NL+1, ND |
1256 |
Mike(i) = 0. |
1257 |
enddo |
1258 |
|
1259 |
do i=1,nd |
1260 |
Ma(i)=0 |
1261 |
enddo |
1262 |
|
1263 |
c sb do i=1,nd |
1264 |
c sb do j=i,nd |
1265 |
c sb Ma(i)=Ma(i)+M(j) |
1266 |
c sb enddo |
1267 |
c sb enddo |
1268 |
|
1269 |
do i = 1, NL |
1270 |
do j = i, NL |
1271 |
Ma(i) = Ma(i) + M(j) |
1272 |
enddo |
1273 |
enddo |
1274 |
c |
1275 |
do i = NL+1, ND |
1276 |
Ma(i) = 0. |
1277 |
enddo |
1278 |
c |
1279 |
do i=1,ICB-1 |
1280 |
Ma(i)=0 |
1281 |
enddo |
1282 |
|
1283 |
|
1284 |
|
1285 |
ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
1286 |
C ICB REPRESENTE DE NIVEAU OU SE TROUVE LA |
1287 |
c BASE DU NUAGE , ET INB LE TOP DU NUAGE |
1288 |
cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
1289 |
|
1290 |
|
1291 |
do i=1,ND |
1292 |
Mke(i)=upwd(i)+dnwd(i) |
1293 |
enddo |
1294 |
|
1295 |
C |
1296 |
C *** Diagnose the in-cloud mixing ratio *** ! cld |
1297 |
C *** of condensed water *** ! cld |
1298 |
C ! cld |
1299 |
DO I=1,ND ! cld |
1300 |
MAA(I)=0.0 ! cld |
1301 |
WA(I)=0.0 ! cld |
1302 |
SIGA(I)=0.0 ! cld |
1303 |
ENDDO ! cld |
1304 |
DO I=NK,INB ! cld |
1305 |
DO K=I+1,INB+1 ! cld |
1306 |
MAA(I)=MAA(I)+M(K) ! cld |
1307 |
ENDDO ! cld |
1308 |
ENDDO ! cld |
1309 |
DO I=ICB,INB-1 ! cld |
1310 |
AXC(I)=0. ! cld |
1311 |
DO J=ICB,I ! cld |
1312 |
AXC(I)=AXC(I)+RD*(TVP(J)-TV(J))*(PH(J)-PH(J+1))/P(J) ! cld |
1313 |
ENDDO ! cld |
1314 |
IF (AXC(I).GT.0.0) THEN ! cld |
1315 |
WA(I)=SQRT(2.*AXC(I)) ! cld |
1316 |
ENDIF ! cld |
1317 |
ENDDO ! cld |
1318 |
DO I=1,NL ! cld |
1319 |
IF (WA(I).GT.0.0) ! cld |
1320 |
1 SIGA(I)=MAA(I)/WA(I)*RD*TVP(I)/P(I)/100./DELTAC ! cld |
1321 |
SIGA(I) = MIN(SIGA(I),1.0) ! cld |
1322 |
QCONDC(I)=SIGA(I)*CLW(I)*(1.-EP(I)) ! cld |
1323 |
1 + (1.-SIGA(I))*QCOND(I) ! cld |
1324 |
ENDDO ! cld |
1325 |
|
1326 |
|
1327 |
c$$$cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
1328 |
c$$$ call writeg1d(1,klev,ma,'ma ','ma ') |
1329 |
c$$$ call writeg1d(1,klev,upwd,'upwd ','upwd ') |
1330 |
c$$$ call writeg1d(1,klev,dnwd,'dnwd ','dnwd ') |
1331 |
c$$$ call writeg1d(1,klev,dnwd0,'dnwd0 ','dnwd0 ') |
1332 |
c$$$ call writeg1d(1,klev,tvp,'tvp ','tvp ') |
1333 |
c$$$ call writeg1d(1,klev,tra(1:klev,3),'tra3 ','tra3 ') |
1334 |
c$$$ call writeg1d(1,klev,tra(1:klev,4),'tra4 ','tra4 ') |
1335 |
c$$$ call writeg1d(1,klev,tra(1:klev,5),'tra5 ','tra5 ') |
1336 |
c$$$ call writeg1d(1,klev,tra(1:klev,6),'tra6 ','tra6 ') |
1337 |
c$$$ call writeg1d(1,klev,tra(1:klev,7),'tra7 ','tra7 ') |
1338 |
c$$$ call writeg1d(1,klev,tra(1:klev,8),'tra8 ','tra8 ') |
1339 |
c$$$ call writeg1d(1,klev,tra(1:klev,9),'tra9 ','tra9 ') |
1340 |
c$$$ call writeg1d(1,klev,tra(1:klev,10),'tra10','tra10') |
1341 |
c$$$ call writeg1d(1,klev,tra(1:klev,11),'tra11','tra11') |
1342 |
c$$$ call writeg1d(1,klev,tra(1:klev,12),'tra12','tra12') |
1343 |
c$$$ call writeg1d(1,klev,tra(1:klev,13),'tra13','tra13') |
1344 |
c$$$ call writeg1d(1,klev,tra(1:klev,14),'tra14','tra14') |
1345 |
c$$$ call writeg1d(1,klev,tra(1:klev,15),'tra15','tra15') |
1346 |
c$$$ call writeg1d(1,klev,tra(1:klev,16),'tra16','tra16') |
1347 |
c$$$ call writeg1d(1,klev,tra(1:klev,17),'tra17','tra17') |
1348 |
c$$$ call writeg1d(1,klev,tra(1:klev,18),'tra18','tra18') |
1349 |
c$$$ call writeg1d(1,klev,tra(1:klev,19),'tra19','tra19') |
1350 |
c$$$ call writeg1d(1,klev,tra(1:klev,20),'tra20','tra20 ') |
1351 |
c$$$ call writeg1d(1,klev,trap(1:klev,1),'trp1','trp1') |
1352 |
c$$$ call writeg1d(1,klev,trap(1:klev,2),'trp2','trp2') |
1353 |
c$$$ call writeg1d(1,klev,trap(1:klev,3),'trp3','trp3') |
1354 |
c$$$ call writeg1d(1,klev,trap(1:klev,4),'trp4','trp4') |
1355 |
c$$$ call writeg1d(1,klev,trap(1:klev,5),'trp5','trp5') |
1356 |
c$$$ call writeg1d(1,klev,trap(1:klev,10),'trp10','trp10') |
1357 |
c$$$ call writeg1d(1,klev,trap(1:klev,12),'trp12','trp12') |
1358 |
c$$$ call writeg1d(1,klev,trap(1:klev,15),'trp15','trp15') |
1359 |
c$$$ call writeg1d(1,klev,trap(1:klev,20),'trp20','trp20') |
1360 |
c$$$ call writeg1d(1,klev,ftra(1:klev,1),'ftr1 ','ftr1 ') |
1361 |
c$$$ call writeg1d(1,klev,ftra(1:klev,2),'ftr2 ','ftr2 ') |
1362 |
c$$$ call writeg1d(1,klev,ftra(1:klev,3),'ftr3 ','ftr3 ') |
1363 |
c$$$ call writeg1d(1,klev,ftra(1:klev,4),'ftr4 ','ftr4 ') |
1364 |
c$$$ call writeg1d(1,klev,ftra(1:klev,5),'ftr5 ','ftr5 ') |
1365 |
c$$$ call writeg1d(1,klev,ftra(1:klev,6),'ftr6 ','ftr6 ') |
1366 |
c$$$ call writeg1d(1,klev,ftra(1:klev,7),'ftr7 ','ftr7 ') |
1367 |
c$$$ call writeg1d(1,klev,ftra(1:klev,8),'ftr8 ','ftr8 ') |
1368 |
c$$$ call writeg1d(1,klev,ftra(1:klev,9),'ftr9 ','ftr9 ') |
1369 |
c$$$ call writeg1d(1,klev,ftra(1:klev,10),'ftr10','ftr10') |
1370 |
c$$$ call writeg1d(1,klev,ftra(1:klev,11),'ftr11','ftr11') |
1371 |
c$$$ call writeg1d(1,klev,ftra(1:klev,12),'ftr12','ftr12') |
1372 |
c$$$ call writeg1d(1,klev,ftra(1:klev,13),'ftr13','ftr13') |
1373 |
c$$$ call writeg1d(1,klev,ftra(1:klev,14),'ftr14','ftr14') |
1374 |
c$$$ call writeg1d(1,klev,ftra(1:klev,15),'ftr15','ftr15') |
1375 |
c$$$ call writeg1d(1,klev,ftra(1:klev,16),'ftr16','ftr16') |
1376 |
c$$$ call writeg1d(1,klev,ftra(1:klev,17),'ftr17','ftr17') |
1377 |
c$$$ call writeg1d(1,klev,ftra(1:klev,18),'ftr18','ftr18') |
1378 |
c$$$ call writeg1d(1,klev,ftra(1:klev,19),'ftr19','ftr19') |
1379 |
c$$$ call writeg1d(1,klev,ftra(1:klev,20),'ftr20','ftr20 ') |
1380 |
c$$$ call writeg1d(1,klev,mp,'mp ','mp ') |
1381 |
c$$$ call writeg1d(1,klev,Mke,'Mke ','Mke ') |
1382 |
|
1383 |
|
1384 |
|
1385 |
ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
1386 |
c |
1387 |
c |
1388 |
RETURN |
1389 |
END |
1390 |
C --------------------------------------------------------------------------- |