/[MITgcm]/MITgcm/pkg/fizhi/fizhi_moist.F
ViewVC logotype

Contents of /MITgcm/pkg/fizhi/fizhi_moist.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.8 - (show annotations) (download)
Tue Jul 13 23:44:43 2004 UTC (19 years, 11 months ago) by molod
Branch: MAIN
Changes since 1.7: +42 -69 lines
Debugging

1 C $Header: /u/gcmpack/MITgcm/pkg/fizhi/fizhi_moist.F,v 1.7 2004/07/13 21:18:41 molod Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6 subroutine moistio (ndmoist,istrip,npcs,
7 . lowlevel,midlevel,nltop,nsubmin,nsubmax,Lup,
8 . pz,plz,plze,dpres,pkht,pkl,tz,qz,bi,bj,ntracer,ptracer,
9 . qqz,dumoist,dvmoist,dtmoist,dqmoist,
10 . im,jm,lm,ptop,
11 . iras,rainlsp,rainconv,snowfall,
12 . nswcld,cldtot_sw,cldras_sw,cldlsp_sw,nswlz,swlz,
13 . nlwcld,cldtot_lw,cldras_lw,cldlsp_lw,nlwlz,lwlz,
14 . lpnt,myid)
15
16 implicit none
17
18 #ifdef ALLOW_DIAGNOSTICS
19 #include "SIZE.h"
20 #include "diagnostics_SIZE.h"
21 #include "diagnostics.h"
22 #endif
23
24 c Input Variables
25 c ---------------
26 integer im,jm,lm
27 integer ndmoist,istrip,npcs
28 integer bi,bj,ntracer,ptracer
29 integer lowlevel,midlevel,nltop,nsubmin,nsubmax,Lup
30 real pz(im,jm),plz(im,jm,lm),plze(im,jm,lm+1),dpres(im,jm,lm)
31 real pkht(im,jm,lm+1),pkl(im,jm,lm)
32 real tz(im,jm,lm),qz(im,jm,lm,ntracer)
33 real qqz(im,jm,lm)
34 real dumoist(im,jm,lm),dvmoist(im,jm,lm)
35 real dtmoist(im,jm,lm),dqmoist(im,jm,lm,ntracer)
36 real ptop
37 integer iras
38 real rainlsp(im,jm),rainconv(im,jm),snowfall(im,jm)
39 integer nswcld,nswlz
40 real cldlsp_sw(im,jm,lm),cldras_sw(im,jm,lm)
41 real cldtot_sw(im,jm,lm),swlz(im,jm,lm)
42 integer nlwcld,nlwlz
43 real cldlsp_lw(im,jm,lm),cldras_lw(im,jm,lm)
44 real cldtot_lw(im,jm,lm),lwlz(im,jm,lm)
45 logical lpnt
46 integer myid
47
48 c Local Variables
49 c ---------------
50 integer ncrnd,nsecf
51
52 real fracqq, rh,temp1,temp2,dum
53 integer snowcrit
54 parameter (fracqq = 0.1)
55
56 real cldsr(im,jm,lm)
57 real srcld(istrip,lm)
58
59 real plev
60 real cldnow,cldlsp_mem,cldlsp,cldras_mem,cldras
61 real watnow,watmin,cldmin
62 real cldprs(im,jm),cldtmp(im,jm)
63 real cldhi (im,jm),cldlow(im,jm)
64 real cldmid(im,jm),totcld(im,jm)
65
66 real CLDLS(im,jm,lm) , CPEN(im,jm,lm)
67 real tmpimjm(im,jm)
68 real lsp_new(im,jm)
69 real conv_new(im,jm)
70 real snow_new(im,jm)
71
72 real qqcolmin(im,jm)
73 real qqcolmax(im,jm)
74 integer levpbl(im,jm)
75
76 c Gathered Arrays for Variable Cloud Base
77 c ---------------------------------------
78 real raincgath(im*jm)
79 real pigather(im*jm)
80 real thgather(im*jm,lm)
81 real shgather(im*jm,lm)
82 real pkzgather(im*jm,lm)
83 real pkegather(im*jm,lm+1)
84 real plzgather(im*jm,lm)
85 real plegather(im*jm,lm+1)
86 real dpgather(im*jm,lm)
87 real tmpgather(im*jm,lm)
88 real deltgather(im*jm,lm)
89 real delqgather(im*jm,lm)
90 real ugather(im*jm,lm,ntracer)
91 real delugather(im*jm,lm,ntracer)
92 real deltrnev(im*jm,lm)
93 real delqrnev(im*jm,lm)
94
95 integer nindeces(lm)
96 integer pblindex(im*jm)
97 integer levgather(im*jm)
98
99 c Stripped Arrays
100 c ---------------
101 real saveth (istrip,lm)
102 real saveq (istrip,lm)
103 real saveu (istrip,lm,ntracer)
104 real usubcl (istrip, ntracer)
105
106 real ple(istrip,lm+1), gam(istrip,lm)
107 real dp(istrip,lm)
108 real TL(ISTRIP,lm) , SHL(ISTRIP,lm)
109 real PL(ISTRIP,lm) , PLK(ISTRIP,lm)
110 real PLKE(ISTRIP,lm+1)
111 real TH(ISTRIP,lm) ,CVTH(ISTRIP,lm)
112 real SHSAT(ISTRIP,lm) , CVQ(ISTRIP,lm)
113 real UL(ISTRIP,lm,ntracer)
114 real cvu(istrip,lm,ntracer)
115 real CLMAXO(ISTRIP,lm),CLBOTH(ISTRIP,lm)
116 real CLSBTH(ISTRIP,lm)
117 real TMP1(ISTRIP,lm), TMP2(ISTRIP,lm)
118 real TMP3(ISTRIP,lm), TMP4(ISTRIP,lm+1)
119 real TMP5(ISTRIP,lm+1)
120 integer ITMP1(ISTRIP,lm), ITMP2(ISTRIP,lm)
121 integer ITMP3(ISTRIP,lm)
122
123 real PRECIP(ISTRIP), PCMID(ISTRIP), PCNET(ISTRIP)
124 real PCLOW (ISTRIP), SP(ISTRIP), PREP(ISTRIP)
125 real PCPEN (ISTRIP,lm)
126 integer pbl(istrip),depths(lm)
127
128 real cldlz(istrip,lm), cldwater(im,jm,lm)
129 real rhfrac(istrip), rhmin, pup, ppbl, rhcrit(istrip,lm)
130 real offset, alpha, rasmax
131
132 logical first
133 logical lras
134 real clfrac (istrip,lm)
135 real cldmas (istrip,lm)
136 real detrain(istrip,lm)
137 real psubcld (istrip), psubcldg (im,jm)
138 real psubcld_cnt(istrip), psubcldgc(im,jm)
139 real rnd(lm/2)
140 DATA FIRST /.TRUE./
141
142 integer imstp,nsubcl,nlras
143 integer i,j,iloop,index,l,nn,num,numdeps,nt
144 real tmstp,tminv,sday,grav,alhl,cp,elocp,gamfac
145 real rkappa,p0kappa,p0kinv,ptopkap,pcheck
146 real tice,getcon,pi
147
148 C **********************************************************************
149 C **** INITIALIZATION ****
150 C **********************************************************************
151
152 IMSTP = nsecf(NDMOIST)
153 TMSTP = FLOAT(IMSTP)
154 TMINV = 1. / TMSTP
155
156 C Minimum Large-Scale Cloud Fraction at rhcrit
157 alpha = 0.80
158 C Difference in fraction between SR and LS Threshold
159 offset = 0.10
160 C Large-Scale Relative Humidity Threshold in PBL
161 rhmin = 0.90
162 C Maximum Cloud Fraction associated with RAS
163 rasmax = 1.00
164
165 nn = 3*3600.0/tmstp + 1
166 C Threshold for Cloud Fraction Memory
167 cldmin = rasmax*(1.0-tmstp/3600.)**nn
168 C Threshold for Cloud Liquid Water Memory
169 watmin = 1.0e-8
170
171 SDAY = GETCON('SDAY')
172 GRAV = GETCON('GRAVITY')
173 ALHL = GETCON('LATENT HEAT COND')
174 CP = GETCON('CP')
175 ELOCP = GETCON('LATENT HEAT COND') / GETCON('CP')
176 GAMFAC = GETCON('LATENT HEAT COND') * GETCON('EPS') * ELOCP
177 . / GETCON('RGAS')
178 RKAPPA = GETCON('KAPPA')
179 P0KAPPA = 1000.0**RKAPPA
180 P0KINV = 1. / P0KAPPA
181 PTOPKAP = PTOP**RKAPPA
182 tice = getcon('FREEZING-POINT')
183 PI = 4.*atan(1.)
184
185 c Determine Total number of Random Clouds to Check
186 c ---------------------------------------------
187 ncrnd = (lm-nltop+1)/2
188
189 if(first .and. myid.eq.0) then
190 print *
191 print *,'Top Level Allowed for Convection : ',nltop
192 print *,' Highest Sub-Cloud Level: ',nsubmax
193 print *,' Lowest Sub-Cloud Level: ',nsubmin
194 print *,' Total Number of Random Clouds: ',ncrnd
195 print *
196 first = .false.
197 endif
198
199 c And now find PBL depth - the level where qq = fracqq * qq at surface
200 c --------------------------------------------------------------------
201 do j = 1,jm
202 do i = 1,im
203 qqcolmin(i,j) = qqz(i,j,lm)*fracqq
204 qqcolmax(i,j) = qqz(i,j,lm)
205 levpbl(i,j) = lm+1
206 enddo
207 enddo
208
209 DO L = lm-1,1,-1
210 DO j = 1,jm
211 DO i = 1,im
212 IF((qqz(i,j,l).gt.qqcolmax(i,j))
213 1 .and.(levpbl(i,j).eq.lm+1))then
214 qqcolmax(i,j) = qqz(i,j,l)
215 qqcolmin(i,j) = fracqq*qqcolmax(i,j)
216 endif
217 if((qqz(i,j,l).lt.qqcolmin(i,j))
218 1 .and.(levpbl(i,j).eq.lm+1))levpbl(i,j)=L+1
219 enddo
220 enddo
221 enddo
222
223 do j = 1,jm
224 do i = 1,im
225 if(levpbl(i,j).gt.nsubmin) levpbl(i,j) = nsubmin
226 if(levpbl(i,j).lt.nsubmax) levpbl(i,j) = nsubmax
227 enddo
228 enddo
229
230
231 c Set up the array of indeces of subcloud levels for the gathering
232 c ----------------------------------------------------------------
233 index = 0
234 do L = nsubmin,nltop,-1
235 do j = 1,jm
236 do i = 1,im
237 if(levpbl(i,j).eq.L) then
238 index = index + 1
239 pblindex(index) = (j-1)*im + i
240 endif
241 enddo
242 enddo
243 enddo
244
245 do index = 1,im*jm
246 levgather(index) = levpbl(pblindex(index),1)
247 pigather(index) = pz(pblindex(index),1)
248 pkegather(index,lm+1) = pkht(pblindex(index),1,lm+1)
249 plegather(index,lm+1) = plze(pblindex(index),1,lm+1)
250 enddo
251
252 do L = 1,lm
253 do index = 1,im*jm
254 thgather(index,L) = tz(pblindex(index),1,L)
255 shgather(index,L) = qz(pblindex(index),1,L,1)
256 pkegather(index,L) = pkht(pblindex(index),1,L)
257 pkzgather(index,L) = pkl(pblindex(index),1,L)
258 plegather(index,L) = plze(pblindex(index),1,L)
259 plzgather(index,L) = plz(pblindex(index),1,L)
260 dpgather(index,L) = dpres(pblindex(index),1,L)
261 enddo
262 enddo
263 do nt = 1,ntracer-ptracer
264 do L = 1,lm
265 do index = 1,im*jm
266 ugather(index,L,nt) = qz(pblindex(index),1,L,nt+ptracer)
267 enddo
268 enddo
269 enddo
270
271 c bump the counter for number of calls to convection
272 c --------------------------------------------------
273 iras = iras + 1
274 if( iras.ge.1e9 ) iras = 1
275
276 c select the 'random' cloud detrainment levels for RAS
277 c ----------------------------------------------------
278 call rndcloud(iras,ncrnd,rnd,myid)
279
280 do l=1,lm
281 do j=1,jm
282 do i=1,im
283 dtmoist(i,j,l) = 0.
284 do nt = 1,ntracer
285 dqmoist(i,j,l,nt) = 0.
286 enddo
287 enddo
288 enddo
289 enddo
290
291 C***********************************************************************
292 C **** LOOP OVER NPCS PEICES ****
293 C **********************************************************************
294
295 DO 1000 NN = 1,NPCS
296
297 C **********************************************************************
298 C **** VARIABLE INITIALIZATION ****
299 C **********************************************************************
300
301 CALL STRIP ( pigather, SP ,im*jm,ISTRIP,1 ,NN )
302 CALL STRIP ( pkzgather, PLK ,im*jm,ISTRIP,lm,NN )
303 CALL STRIP ( pkegather, PLKE ,im*jm,ISTRIP,lm+1,NN )
304 CALL STRIP ( plzgather, PL ,im*jm,ISTRIP,lm,NN )
305 CALL STRIP ( plegather, PLE ,im*jm,ISTRIP,lm+1,NN )
306 CALL STRIP ( dpgather, dp ,im*jm,ISTRIP,lm,NN )
307 CALL STRIP ( thgather, TH ,im*jm,ISTRIP,lm,NN )
308 CALL STRIP ( shgather, SHL ,im*jm,ISTRIP,lm,NN )
309 CALL STRINT( levgather, pbl ,im*jm,ISTRIP,1 ,NN )
310
311 do nt = 1,ntracer-ptracer
312 call strip ( ugather(1,1,nt), ul(1,1,nt),im*jm,istrip,lm,nn )
313 enddo
314
315 C **********************************************************************
316 C **** SETUP FOR RAS CUMULUS PARAMETERIZATION ****
317 C **********************************************************************
318
319 DO L = 1,lm
320 DO I = 1,ISTRIP
321 TH(I,L) = TH(I,L) * P0KAPPA
322 CLMAXO(I,L) = 0.
323 CLBOTH(I,L) = 0.
324 cldmas(I,L) = 0.
325 detrain(I,L) = 0.
326 ENDDO
327 ENDDO
328
329 do L = 1,lm
330 depths(L) = 0
331 enddo
332
333 numdeps = 0
334 do L = nsubmin,nltop,-1
335 nindeces(L) = 0
336 do i = 1,istrip
337 if(pbl(i).eq.L) nindeces(L) = nindeces(L) + 1
338 enddo
339 if(nindeces(L).gt.0) then
340 numdeps = numdeps + 1
341 depths(numdeps) = L
342 endif
343 enddo
344
345
346 C Initiate a do-loop around RAS for the number of different
347 C sub-cloud layer depths in this strip
348 C --If all subcloud depths are the same, execute loop once
349 C Otherwise loop over different subcloud layer depths
350
351 num = 1
352 DO iloop = 1,numdeps
353
354 nsubcl = depths(iloop)
355
356 c Compute sub-cloud values for Temperature and Spec.Hum.
357 c ------------------------------------------------------
358 DO 600 I=num,num+nindeces(nsubcl)-1
359 TMP1(I,2) = 0.
360 TMP1(I,3) = 0.
361 600 CONTINUE
362
363 NLRAS = NSUBCL - NLTOP + 1
364 DO 601 L=NSUBCL,lm
365 DO 602 I=num,num+nindeces(nsubcl)-1
366 TMP1(I,2) = TMP1(I,2) + (PLE(I,L+1)-PLE(I,L))*TH (I,L)/sp(i)
367 TMP1(I,3) = TMP1(I,3) + (PLE(I,L+1)-PLE(I,L))*SHL(I,L)/sp(i)
368 602 CONTINUE
369 601 CONTINUE
370 DO 603 I=num,num+nindeces(nsubcl)-1
371 TMP1(I,4) = 1. / ( (PLE(I,lm+1)-PLE(I,NSUBCL))/sp(I) )
372 TH(I,NSUBCL) = TMP1(I,2)*TMP1(I,4)
373 SHL(I,NSUBCL) = TMP1(I,3)*TMP1(I,4)
374 603 CONTINUE
375
376 c Save initial value of tracers and compute sub-cloud value
377 c ---------------------------------------------------------
378 DO NT = 1,ntracer-ptracer
379 do L = 1,lm
380 do i = num,num+nindeces(nsubcl)-1
381 saveu(i,L,nt) = ul(i,L,nt)
382 enddo
383 enddo
384 DO I=num,num+nindeces(nsubcl)-1
385 TMP1(I,2) = 0.
386 ENDDO
387 DO L=NSUBCL,lm
388 DO I=num,num+nindeces(nsubcl)-1
389 TMP1(I,2) = TMP1(I,2)+(PLE(I,L+1)-PLE(I,L))*UL(I,L,NT)/sp(i)
390 ENDDO
391 ENDDO
392 DO I=num,num+nindeces(nsubcl)-1
393 UL(I,NSUBCL,NT) = TMP1(I,2)*TMP1(I,4)
394 usubcl(i,nt) = ul(i,nsubcl,nt)
395 ENDDO
396 ENDDO
397
398 c Compute Pressure Arrays for RAS
399 c -------------------------------
400 DO 111 L=1,lm
401 DO 112 I=num,num+nindeces(nsubcl)-1
402 TMP4(I,L) = PLE(I,L)
403 112 CONTINUE
404 111 CONTINUE
405 DO I=num,num+nindeces(nsubcl)-1
406 TMP5(I,1) = PTOPKAP / P0KAPPA
407 ENDDO
408 DO L=2,lm
409 DO I=num,num+nindeces(nsubcl)-1
410 TMP5(I,L) = PLKE(I,L)*P0KINV
411 ENDDO
412 ENDDO
413 DO I=num,num+nindeces(nsubcl)-1
414 TMP4(I,lm+1) = PLE (I,lm+1)
415 TMP5(I,lm+1) = PLKE(I,lm+1)*P0KINV
416 ENDDO
417 DO 113 I=num,num+nindeces(nsubcl)-1
418 TMP4(I,NSUBCL+1) = PLE (I,lm+1)
419 TMP5(I,NSUBCL+1) = PLKE(I,lm+1)*P0KINV
420 113 CONTINUE
421
422 do i=num,num+nindeces(nsubcl)-1
423 C Temperature at top of sub-cloud layer
424 tmp2(i,1) = TH(i,NSUBCL) * PLKE(i,NSUBCL)/P0KAPPA
425 C Pressure at top of sub-cloud layer
426 tmp2(i,2) = tmp4(i,nsubcl)
427 enddo
428
429 C CHANGED THIS: no RH requirement for RAS
430 c call vqsat ( tmp2(num,1),tmp2(num,2),tmp2(num,3),
431 c . dum,.false.,nindeces(nsubcl) )
432 c do i=num,num+nindeces(nsubcl)-1
433 c rh = SHL(I,NSUBCL) / tmp2(i,3)
434 c if (rh .le. 0.85) then
435 c rhfrac(i) = 0.
436 c else if (rh .ge. 0.95) then
437 c rhfrac(i) = 1.
438 c else
439 c rhfrac(i) = (rh-0.85)*10.
440 c endif
441 c enddo
442 do i=num,num+nindeces(nsubcl)-1
443 rhfrac(i) = 1.
444 enddo
445
446 C Compute RH threshold for Large-scale condensation
447 C Used in Slingo-Ritter clouds as well - define offset between SR and LS
448
449 C Top level of atan func above this rh_threshold = rhmin
450 pup = 600.
451 do i=num,num+nindeces(nsubcl)-1
452 do L = nsubcl, lm
453 rhcrit(i,L) = 1.
454 enddo
455 do L = 1, nsubcl-1
456 pcheck = pl(i,L)
457 if (pcheck .le. pup) then
458 rhcrit(i,L) = rhmin
459 else
460 ppbl = pl(i,nsubcl)
461 rhcrit(i,L) = rhmin + (1.-rhmin)/(19.) *
462 . ((atan( (2.*(pcheck-pup)/(ppbl-pup)-1.) *
463 . tan(20.*pi/21.-0.5*pi) )
464 . + 0.5*pi) * 21./pi - 1.)
465 endif
466 enddo
467 enddo
468
469 c Save Initial Values of Temperature and Specific Humidity
470 c --------------------------------------------------------
471 do L = 1,lm
472 do i = num,num+nindeces(nsubcl)-1
473 saveth(i,L) = th (i,L)
474 saveq (i,L) = shl(i,L)
475 PCPEN (i,L) = 0.
476 CLFRAC(i,L) = 0.
477 enddo
478 enddo
479
480 CALL RAS ( NN,istrip,nindeces(nsubcl),NLRAS,NLTOP,lm,TMSTP
481 1, UL(num,1,1),ntracer-ptracer,TH(num,NLTOP),SHL(num,NLTOP)
482 2, TMP4(num,NLTOP), TMP5(num,NLTOP),rnd, ncrnd, PCPEN(num,NLTOP)
483 3, CLBOTH(num,NLTOP), CLFRAC(num,NLTOP)
484 4, cldmas(num,nltop), detrain(num,nltop)
485 8, cp,grav,rkappa,alhl,rhfrac(num),rasmax )
486
487 c Compute Diagnostic CLDMAS in RAS Subcloud Layers
488 c ------------------------------------------------
489 do L=nsubcl,lm
490 do I=num,num+nindeces(nsubcl)-1
491 dum = dp(i,L)/(ple(i,lm+1)-ple(i,nsubcl))
492 cldmas(i,L) = cldmas(i,L-1) - dum*cldmas(i,nsubcl-1)
493 enddo
494 enddo
495
496 c Update Theta and Moisture due to RAS
497 c ------------------------------------
498 DO L=1,nsubcl
499 DO I=num,num+nindeces(nsubcl)-1
500 CVTH(I,L) = (TH (I,L) - saveth(i,l))
501 CVQ (I,L) = (SHL(I,L) - saveq (i,l))
502 ENDDO
503 ENDDO
504 DO L=nsubcl+1,lm
505 DO I=num,num+nindeces(nsubcl)-1
506 CVTH(I,L) = cvth(i,nsubcl)
507 CVQ (I,L) = cvq (i,nsubcl)
508 ENDDO
509 ENDDO
510
511 DO L=nsubcl+1,lm
512 DO I=num,num+nindeces(nsubcl)-1
513 TH (I,L) = saveth(i,l) + cvth(i,l)
514 SHL(I,L) = saveq (i,l) + cvq (i,l)
515 ENDDO
516 ENDDO
517 DO L=1,lm
518 DO I=num,num+nindeces(nsubcl)-1
519 CVTH(I,L) = CVTH(I,L) *P0KINV*SP(I)*tminv
520 CVQ (I,L) = CVQ (I,L) *SP(I)*tminv
521 ENDDO
522 ENDDO
523
524 c Compute Tracer Tendency due to RAS
525 c ----------------------------------
526 do nt = 1,ntracer-ptracer
527 DO L=1,nsubcl-1
528 DO I=num,num+nindeces(nsubcl)-1
529 CVU(I,L,nt) = ( UL(I,L,nt)-saveu(i,l,nt) )*sp(i)*tminv
530 ENDDO
531 ENDDO
532 DO L=nsubcl,lm
533 DO I=num,num+nindeces(nsubcl)-1
534 if( usubcl(i,nt).ne.0.0 ) then
535 cvu(i,L,nt) = ( ul(i,nsubcl,nt)-usubcl(i,nt) ) *
536 . ( saveu(i,L,nt)/usubcl(i,nt) )*sp(i)*tminv
537 else
538 cvu(i,L,nt) = 0.0
539 endif
540 ENDDO
541 ENDDO
542 enddo
543
544 c Compute Diagnostic PSUBCLD (Subcloud Layer Pressure)
545 c ----------------------------------------------------
546 do i=num,num+nindeces(nsubcl)-1
547 lras = .false.
548 do L=nltop,nsubcl
549 if( cvq(i,L).ne.0.0 ) lras = .true.
550 enddo
551 psubcld (i) = 0.0
552 psubcld_cnt(i) = 0.0
553 if( lras ) then
554 psubcld (i) = sp(i)+ptop-ple(i,nsubcl)
555 psubcld_cnt(i) = 1.0
556 endif
557 enddo
558
559
560 C End of subcloud layer depth loop (iloop)
561
562 num = num+nindeces(nsubcl)
563
564 ENDDO
565
566 C **********************************************************************
567 C **** TENDENCY UPDATES ****
568 C **** (Keep 'Gathered' tendencies in 'gather' arrays now) ****
569 C **********************************************************************
570
571 call paste( CVTH,deltgather,istrip,im*jm,lm,NN )
572 call paste( CVQ,delqgather,istrip,im*jm,lm,NN )
573 do nt = 1,ntracer-ptracer
574 call paste( CVU(1,1,nt),delugather(1,1,nt),istrip,im*jm,lm,NN )
575 enddo
576
577 C **********************************************************************
578 C And now paste some arrays for filling diagnostics
579 C (use pkegather to hold detrainment and tmpgather for cloud mass flux)
580 C **********************************************************************
581
582 if(icldmas .gt.0) call paste( cldmas,tmpgather,istrip,im*jm,lm,NN)
583 if(idtrain .gt.0) call paste(detrain,pkegather,istrip,im*jm,lm,NN)
584 if(ipsubcld.gt.0) then
585 call paste(psubcld ,psubcldg ,istrip,im*jm,1,NN)
586 call paste(psubcld_cnt,psubcldgc,istrip,im*jm,1,NN)
587 endif
588
589 C *********************************************************************
590 C **** RE-EVAPORATION OF PENETRATING CONVECTIVE RAIN ****
591 C *********************************************************************
592
593 CALL STRIP ( thgather,TH ,im*jm,ISTRIP,lm,NN)
594 CALL STRIP ( shgather,SHL,im*jm,ISTRIP,lm,NN)
595 DO L=1,lm
596 DO I=1,ISTRIP
597 TH(I,L) = TH(I,L) + CVTH(I,L)*tmstp/SP(I)
598 SHL(I,L) = SHL(I,L) + CVQ(I,L)*tmstp/SP(I)
599 TL(I,L) = TH(I,L)*PLK(I,L)
600 saveth(I,L) = th(I,L)
601 saveq (I,L) = SHL(I,L)
602 ENDDO
603 ENDDO
604
605 CALL RNEVP (NN,ISTRIP,lm,TL,SHL,PCPEN,PL,CLFRAC,SP,DP,PLKE,
606 . PLK,TH,TMP1,TMP2,TMP3,ITMP1,ITMP2,PCNET,PRECIP,
607 . CLSBTH,TMSTP,1.,cp,grav,alhl,gamfac,cldlz,rhcrit,offset,alpha)
608
609 C **********************************************************************
610 C **** TENDENCY UPDATES ****
611 C **********************************************************************
612
613 DO L=1,lm
614
615 DO I =1,ISTRIP
616 TMP1(I,L) = sp(i) * (SHL(I,L)-saveq(I,L)) * tminv
617 ENDDO
618 CALL PSTBMP(TMP1(1,L),delqgather(1,L),ISTRIP,im*jm,1,NN)
619
620 DO I =1,ISTRIP
621 TMP1(I,L) = sp(i) * ((TL(I,L)/PLK(I,L))-saveth(i,l)) * tminv
622 ENDDO
623 CALL PSTBMP(TMP1(1,L),deltgather(1,L),ISTRIP,im*jm,1,NN)
624
625 C Paste rain evap tendencies into arrays for diagnostic output
626 c ------------------------------------------------------------
627 if(idtls.gt.0)then
628 DO I =1,ISTRIP
629 TMP1(I,L) = ((TL(I,L)/PLK(I,L))-saveth(i,l))*plk(i,l)*sday*tminv
630 ENDDO
631 call paste(tmp1(1,L),deltrnev(1,L),istrip,im*jm,1,NN)
632 endif
633
634 if(idqls.gt.0)then
635 DO I =1,ISTRIP
636 TMP1(I,L) = (SHL(I,L)-saveq(I,L)) * 1000. * sday * tminv
637 ENDDO
638 call paste(tmp1(1,L),delqrnev(1,L),istrip,im*jm,1,NN)
639 endif
640
641 ENDDO
642
643 C *********************************************************************
644 C Add Non-Precipitating Clouds where the relative
645 C humidity is less than 100%
646 C Apply Cloud Top Entrainment Instability
647 C *********************************************************************
648
649 do L=1,lm
650 do i=1,istrip
651 srcld(i,L) = -clsbth(i,L)
652 enddo
653 enddo
654
655 call srclouds (saveth,saveq,plk,pl,plke,clsbth,cldlz,istrip,lm,
656 . rhcrit,offset,alpha)
657
658 do L=1,lm
659 do i=1,istrip
660 srcld(i,L) = srcld(i,L)+clsbth(i,L)
661 enddo
662 enddo
663
664 C *********************************************************************
665 C **** PASTE CLOUD AMOUNTS ****
666 C *********************************************************************
667
668 call paste ( srcld, cldsr,istrip,im*jm,lm,nn )
669 call paste ( cldlz,cldwater,istrip,im*jm,lm,nn )
670 call paste ( clsbth, cldls,istrip,im*jm,lm,nn )
671 call paste ( clboth, cpen ,istrip,im*jm,lm,nn )
672
673 c compute Total Accumulated Precip for Landsurface Model
674 c ------------------------------------------------------
675 do i = 1,istrip
676 C Initialize Rainlsp, Rainconv and Snowfall
677 tmp1(i,1) = 0.0
678 tmp1(i,2) = 0.0
679 tmp1(i,3) = 0.0
680 enddo
681
682 do i = 1,istrip
683 prep(i) = PRECIP(I) + PCNET(I)
684 tmp1(i,1) = PRECIP(I)
685 tmp1(i,2) = pcnet(i)
686 enddo
687 c
688 c check whether there is snow
689 c-------------------------------------------------------
690 c snow algorthm:
691 c if temperature profile from the surface level to 700 mb
692 c uniformaly c below zero, then precipitation (total) is
693 c snowfall. Else there is no snow.
694 c-------------------------------------------------------
695
696 do i = 1,istrip
697 snowcrit=0
698 do l=lup,lm
699 if (saveth(i,l)*plk(i,l).le. tice ) then
700 snowcrit=snowcrit+1
701 endif
702 enddo
703 if (snowcrit .eq. (lm-lup+1)) then
704 tmp1(i,3) = prep(i)
705 tmp1(i,1)=0.0
706 tmp1(i,2)=0.0
707 endif
708 enddo
709
710 CALL paste (tmp1(1,1), lsp_new,ISTRIP,im*jm,1,NN)
711 CALL paste (tmp1(1,2),conv_new,ISTRIP,im*jm,1,NN)
712 CALL paste (tmp1(1,3),snow_new,ISTRIP,im*jm,1,NN)
713
714 if(iprecon.gt.0) then
715 CALL paste (pcnet,raincgath,ISTRIP,im*jm,1,NN)
716 endif
717
718 C *********************************************************************
719 C **** End Major Stripped Region ****
720 C *********************************************************************
721
722 1000 CONTINUE
723
724 C Large Scale Rainfall, Conv rain, and snowfall
725 c ---------------------------------------------
726 call back2grd ( lsp_new,pblindex, lsp_new,im*jm)
727 call back2grd (conv_new,pblindex,conv_new,im*jm)
728 call back2grd (snow_new,pblindex,snow_new,im*jm)
729
730 if(iprecon.gt.0) then
731 call back2grd (raincgath,pblindex,raincgath,im*jm)
732 endif
733
734 c Subcloud Layer Pressure
735 c -----------------------
736 if(ipsubcld.gt.0) then
737 call back2grd (psubcldg ,pblindex,psubcldg ,im*jm)
738 call back2grd (psubcldgc,pblindex,psubcldgc,im*jm)
739 endif
740
741 do L = 1,lm
742 C Delta theta,q, convective, max and ls clouds
743 c --------------------------------------------
744 call back2grd (deltgather(1,L),pblindex, dtmoist(1,1,L) ,im*jm)
745 call back2grd (delqgather(1,L),pblindex, dqmoist(1,1,L,1),im*jm)
746 call back2grd ( cpen(1,1,L),pblindex, cpen(1,1,L) ,im*jm)
747 call back2grd ( cldls(1,1,L),pblindex, cldls(1,1,L) ,im*jm)
748 call back2grd (cldwater(1,1,L),pblindex,cldwater(1,1,L) ,im*jm)
749 call back2grd ( pkzgather(1,L),pblindex, pkzgather(1,L) ,im*jm)
750
751 C Diagnostics:
752 c ------------
753 if(icldmas.gt.0)call back2grd(tmpgather(1,L),pblindex,
754 . tmpgather(1,L),im*jm)
755 if(idtrain.gt.0)call back2grd(pkegather(1,L),pblindex,
756 . pkegather(1,L),im*jm)
757 if(idtls.gt.0)call back2grd(deltrnev(1,L),pblindex,
758 . deltrnev(1,L),im*jm)
759 if(idqls.gt.0)call back2grd(delqrnev(1,L),pblindex,
760 . delqrnev(1,L),im*jm)
761 if(icldnp.gt.0)call back2grd(cldsr(1,1,L),pblindex,
762 . cldsr(1,1,L),im*jm)
763 enddo
764
765 c Tracers
766 c -------
767 do nt = 1,ntracer-ptracer
768 do L = 1,lm
769 call back2grd (delugather(1,L,nt),pblindex,
770 . dqmoist(1,1,L,ptracer+nt),im*jm)
771 enddo
772 enddo
773
774
775 C **********************************************************************
776 C BUMP DIAGNOSTICS
777 C **********************************************************************
778
779 c Sub-Cloud Layer
780 c -------------------------
781 if( ipsubcld.ne.0 ) then
782 do j = 1,jm
783 do i = 1,im
784 qdiag(i,j,ipsubcld,bi,bj) = qdiag(i,j,ipsubcld,bi,bj) +
785 . psubcldg (i,j)
786 qdiag(i,j,ipsubcldc,bi,bj) = qdiag(i,j,ipsubcldc,bi,bj) +
787 . psubcldgc(i,j)
788 enddo
789 enddo
790 endif
791
792 c Non-Precipitating Cloud Fraction
793 c --------------------------------
794 if( icldnp.ne.0 ) then
795 do L = 1,lm
796 do j = 1,jm
797 do i = 1,im
798 qdiag(i,j,icldnp+L-1,bi,bj) = qdiag(i,j,icldnp+L-1,bi,bj) +
799 . cldsr(i,j,L)
800 enddo
801 enddo
802 enddo
803 ncldnp = ncldnp + 1
804 endif
805
806 c Moist Processes Heating Rate
807 c ----------------------------
808 if(imoistt.gt.0) then
809 do L = 1,lm
810 do i = 1,im*jm
811 qdiag(i,1,imoistt+L-1,bi,bj) = qdiag(i,1,imoistt+L-1,bi,bj) +
812 . (dtmoist(i,1,L)*sday*pkzgather(i,L)/pz(i,1))
813 enddo
814 enddo
815 endif
816
817 c Moist Processes Moistening Rate
818 c -------------------------------
819 if(imoistq.gt.0) then
820 do L = 1,lm
821 do j = 1,jm
822 do i = 1,im
823 qdiag(i,j,imoistq+L-1,bi,bj) = qdiag(i,j,imoistq+L-1,bi,bj) +
824 . (dqmoist(i,j,L,1)*sday*1000.0/pz(i,j))
825 enddo
826 enddo
827 enddo
828 endif
829
830 c Cloud Mass Flux
831 c ---------------
832 if(icldmas.gt.0) then
833 do L = 1,lm
834 do i = 1,im*jm
835 qdiag(i,1,icldmas+L-1,bi,bj) = qdiag(i,1,icldmas+L-1,bi,bj) +
836 . tmpgather(i,L)
837 enddo
838 enddo
839 endif
840
841 c Detrained Cloud Mass Flux
842 c -------------------------
843 if(idtrain.gt.0) then
844 do L = 1,lm
845 do i = 1,im*jm
846 qdiag(i,1,idtrain+L-1,bi,bj) = qdiag(i,1,idtrain+L-1,bi,bj) +
847 . pkegather(i,L)
848 enddo
849 enddo
850 endif
851
852 c Grid-Scale Condensational Heating Rate
853 c --------------------------------------
854 if(idtls.gt.0) then
855 do L = 1,lm
856 do i = 1,im*jm
857 qdiag(i,1,idtls+L-1,bi,bj) = qdiag(i,1,idtls+L-1,bi,bj) +
858 . deltrnev(i,L)
859 enddo
860 enddo
861 endif
862
863 c Grid-Scale Condensational Moistening Rate
864 c -----------------------------------------
865 if(idqls.gt.0) then
866 do L = 1,lm
867 do i = 1,im*jm
868 qdiag(i,1,idqls+L-1,bi,bj) = qdiag(i,1,idqls+L-1,bi,bj) +
869 . delqrnev(i,L)
870 enddo
871 enddo
872 endif
873
874 c Total Precipitation
875 c -------------------
876 if(ipreacc.gt.0) then
877 do j = 1,jm
878 do i = 1,im
879 qdiag(i,j,ipreacc,bi,bj) = qdiag(i,j,ipreacc,bi,bj)
880 . + ( lsp_new(I,j)
881 . + snow_new(I,j)
882 . + conv_new(i,j) ) *sday*tminv
883 enddo
884 enddo
885 endif
886
887 c Convective Precipitation
888 c ------------------------
889 if(iprecon.gt.0) then
890 do i = 1,im*jm
891 qdiag(i,1,iprecon,bi,bj) = qdiag(i,1,iprecon,bi,bj) +
892 . raincgath(i)*sday*tminv
893 enddo
894 endif
895
896 C **********************************************************************
897 C **** Fill Rainfall and Snowfall Arrays for Land Surface Model ****
898 C **** Note: Precip Rates work when DT(turb)<DT(moist) ****
899 C **********************************************************************
900
901 do j = 1,jm
902 do i = 1,im
903 rainlsp (i,j) = rainlsp (i,j) + lsp_new(i,j)*tminv
904 rainconv(i,j) = rainconv(i,j) + conv_new(i,j)*tminv
905 snowfall(i,j) = snowfall(i,j) + snow_new(i,j)*tminv
906 enddo
907 enddo
908
909 C **********************************************************************
910 C *** Compute Time-averaged Quantities for Radiation ***
911 C *** CPEN => Cloud Fraction from RAS ***
912 C *** CLDLS => Cloud Fraction from RNEVP ***
913 C **********************************************************************
914
915 do j = 1,jm
916 do i = 1,im
917 cldhi (i,j) = 0.
918 cldmid(i,j) = 0.
919 cldlow(i,j) = 0.
920 cldtmp(i,j) = 0.
921 cldprs(i,j) = 0.
922 tmpimjm(i,j) = 0.
923 enddo
924 enddo
925
926 c Set Moist-Process Memory Coefficient
927 c ------------------------------------
928 cldras_mem = 1.0-tmstp/ 3600.0
929 cldlsp_mem = 1.0-tmstp/(3600.0*3)
930
931 do L = 1,lm
932 do i = 1,im*jm
933 plev = pl(i,L)
934
935 c Compute Time-averaged Cloud and Water Amounts for Longwave Radiation
936 c --------------------------------------------------------------------
937 watnow = cldwater(i,1,L)
938 if( plev.le.500.0 ) then
939 cldras = min( max( cldras_lw(i,1,L)*cldras_mem,cpen(i,1,L)),1.0)
940 else
941 cldras = 0.0
942 endif
943 cldlsp = min( max( cldlsp_lw(i,1,L)*cldlsp_mem,cldls(i,1,L)),1.0)
944
945 if( cldras.lt.cldmin ) cldras = 0.0
946 if( cldlsp.lt.cldmin ) cldlsp = 0.0
947
948 cldnow = max( cldlsp,cldras )
949
950 lwlz(i,1,L) = ( nlwlz*lwlz(i,1,L) + watnow)/(nlwlz +1)
951 cldtot_lw(i,1,L) = (nlwcld*cldtot_lw(i,1,L) + cldnow)/(nlwcld+1)
952 cldlsp_lw(i,1,L) = (nlwcld*cldlsp_lw(i,1,L) + cldlsp)/(nlwcld+1)
953 cldras_lw(i,1,L) = (nlwcld*cldras_lw(i,1,L) + cldras)/(nlwcld+1)
954
955
956 c Compute Time-averaged Cloud and Water Amounts for Shortwave Radiation
957 c ---------------------------------------------------------------------
958 watnow = cldwater(i,1,L)
959 if( plev.le.500.0 ) then
960 cldras = min( max(cldras_sw(i,1,L)*cldras_mem, cpen(i,1,L)),1.0)
961 else
962 cldras = 0.0
963 endif
964 cldlsp = min( max(cldlsp_sw(i,1,L)*cldlsp_mem,cldls(i,1,L)),1.0)
965
966 if( cldras.lt.cldmin ) cldras = 0.0
967 if( cldlsp.lt.cldmin ) cldlsp = 0.0
968
969 cldnow = max( cldlsp,cldras )
970
971 swlz(i,1,L) = ( nswlz*swlz(i,1,L) + watnow)/(nswlz +1)
972 cldtot_sw(i,1,L) = (nswcld*cldtot_sw(i,1,L) + cldnow)/(nswcld+1)
973 cldlsp_sw(i,1,L) = (nswcld*cldlsp_sw(i,1,L) + cldlsp)/(nswcld+1)
974 cldras_sw(i,1,L) = (nswcld*cldras_sw(i,1,L) + cldras)/(nswcld+1)
975
976
977 c Compute Instantaneous Low-Mid-High Maximum Overlap Cloud Fractions
978 c ----------------------------------------------------------------------
979
980 if( L.lt.midlevel ) cldhi (i,1) = max( cldnow,cldhi (i,1) )
981 if( L.ge.midlevel .and.
982 . L.lt.lowlevel ) cldmid(i,1) = max( cldnow,cldmid(i,1) )
983 if( L.ge.lowlevel ) cldlow(i,1) = max( cldnow,cldlow(i,1) )
984
985 c Compute Cloud-Top Temperature and Pressure
986 c ------------------------------------------
987 cldtmp(i,1) = cldtmp(i,1) + cldnow*pkzgather(i,L)
988 . * ( tz(i,1,L) + dtmoist(i,1,L)*tmstp/pz(i,1) )
989 cldprs(i,1) = cldprs(i,1) + cldnow*plev
990 tmpimjm(i,1) = tmpimjm(i,1) + cldnow
991
992 enddo
993 enddo
994
995 c Compute Instantanious Total 2-D Cloud Fraction
996 c ----------------------------------------------
997 do j = 1,jm
998 do i = 1,im
999 totcld(i,j) = 1.0 - (1.-cldhi (i,j))
1000 . * (1.-cldmid(i,j))
1001 . * (1.-cldlow(i,j))
1002 enddo
1003 enddo
1004
1005
1006 C **********************************************************************
1007 C *** Fill Cloud Top Pressure and Temperature Diagnostic ***
1008 C **********************************************************************
1009
1010 if(icldtmp.gt.0) then
1011 do j = 1,jm
1012 do i = 1,im
1013 if( cldtmp(i,j).gt.0.0 ) then
1014 qdiag(i,j,icldtmp,bi,bj) = qdiag(i,j,icldtmp,bi,bj) +
1015 . cldtmp(i,j)*totcld(i,j)/tmpimjm(i,j)
1016 qdiag(i,j,icttcnt,bi,bj) = qdiag(i,j,icttcnt,bi,bj) +
1017 . totcld(i,j)
1018 endif
1019 enddo
1020 enddo
1021 endif
1022
1023 if(icldprs.gt.0) then
1024 do j = 1,jm
1025 do i = 1,im
1026 if( cldprs(i,j).gt.0.0 ) then
1027 qdiag(i,j,icldprs,bi,bj) = qdiag(i,j,icldprs,bi,bj) +
1028 . cldprs(i,j)*totcld(i,j)/tmpimjm(i,j)
1029 qdiag(i,j,ictpcnt,bi,bj) = qdiag(i,j,ictpcnt,bi,bj) +
1030 . totcld(i,j)
1031 endif
1032 enddo
1033 enddo
1034 endif
1035
1036 C **********************************************************************
1037 C **** INCREMENT COUNTERS ****
1038 C **********************************************************************
1039
1040 nlwlz = nlwlz + 1
1041 nswlz = nswlz + 1
1042
1043 nlwcld = nlwcld + 1
1044 nswcld = nswcld + 1
1045
1046 nmoistt = nmoistt + 1
1047 nmoistq = nmoistq + 1
1048 npreacc = npreacc + 1
1049 nprecon = nprecon + 1
1050
1051 ncldmas = ncldmas + 1
1052 ndtrain = ndtrain + 1
1053
1054 ndtls = ndtls + 1
1055 ndqls = ndqls + 1
1056
1057 RETURN
1058 END
1059 SUBROUTINE RAS( NN, LEN, LENC, K, NLTOP, nlayr, DT
1060 *, UOI, ntracer, POI, QOI, PRS, PRJ, rnd, ncrnd
1061 *, RAINS, CLN, CLF, cldmas, detrain
1062 *, cp,grav,rkappa,alhl,rhfrac,rasmax )
1063 C
1064 C*********************************************************************
1065 C*********************** ARIES MODEL *******************************
1066 C********************* SUBROUTINE RAS *****************************
1067 C********************** 16 MARCH 1988 ******************************
1068 C*********************************************************************
1069 C
1070 implicit none
1071
1072 integer ntracer
1073 integer nltop,nlayr
1074 real UOI(len,nlayr,ntracer), POI(len,K)
1075 real QOI(len,K), PRS(len,K+1), PRJ(len,K+1)
1076 real rnd(ncrnd)
1077 C
1078 real RAINS(len,K), CLN(len,K), CLF(len,K)
1079 real cldmas(len,K), detrain(len,K)
1080 real TCU(len,K), QCU(len,K)
1081 real ucu(len,K,ntracer)
1082 real ALF(len,K), BET(len,K), GAM(len,K)
1083 *, ETA(len,K), HOI(len,K)
1084 *, PRH(len,K), PRI(len,K)
1085 real HST(len,K), QOL(len,K), GMH(len,K)
1086
1087 real TX1(len), TX2(len), TX3(len), TX4(len), TX5(len)
1088 *, TX6(len), TX7(len), TX8(len), TX9(len)
1089 *, TX11(len), TX12(len), TX13(len), TX14(len,ntracer)
1090 *, TX15(len), TX16(len)
1091 *, WFN(len), IA1(len), IA2(len), IA3(len)
1092 real cloudn(len), pcu(len)
1093
1094 real rhfrac(len),rasmax
1095
1096 integer IC(ICM), IRND(icm)
1097 real cmass(len,K)
1098 LOGICAL SETRAS
1099
1100 integer krmin,icm
1101 real rknob, cmb2pa
1102 PARAMETER (KRMIN=01)
1103 PARAMETER (ICM=1000)
1104 PARAMETER (CMB2PA=100.0)
1105 PARAMETER (rknob = 10.)
1106 C
1107 integer i,L,nc
1108 integer km1,kp1,kprv,kcr,kfx,ncmx
1109 real p00, crtmsf, frac, rasblf
1110
1111 do L = 1,k
1112 do I = 1,LENC
1113 rains(i,l) = 0.
1114 enddo
1115 enddo
1116
1117 p00 = 1000.
1118 crtmsf = 0.
1119
1120 C The numerator here is the fraction of the subcloud layer mass flux
1121 C allowed to entrain into the cloud
1122
1123 CCC FRAC = 1./dt
1124 FRAC = 0.5/dt
1125
1126 KM1 = K - 1
1127 KP1 = K + 1
1128 C we want the ras adjustment time scale to be one hour (indep of dt)
1129 RASBLF = 1./3600.
1130 C
1131 KPRV = KM1
1132 C Removed KRMAX parameter
1133 KCR = MIN(KM1,nlayr-2)
1134 KFX = KM1 - KCR
1135 NCMX = KFX + NCRND
1136 C
1137 IF (KFX .GT. 0) THEN
1138 DO NC=1,KFX
1139 IC(NC) = K - NC
1140 ENDDO
1141 ENDIF
1142 C
1143 IF (NCRND .GT. 0) THEN
1144 DO I=1,ncrnd
1145 IRND(I) = (RND(I)-0.0005)*(KCR-KRMIN+1)
1146 IRND(I) = IRND(I) + KRMIN
1147 ENDDO
1148 C
1149 DO NC=1,NCRND
1150 IC(KFX+NC) = IRND(NC)
1151 ENDDO
1152 ENDIF
1153 C
1154 DO 100 NC=1,NCMX
1155 C
1156 IF (NC .EQ. 1 ) THEN
1157 SETRAS = .TRUE.
1158 ELSE
1159 SETRAS = .FALSE.
1160 ENDIF
1161 IB = IC(NC)
1162
1163 c Initialize Cloud Fraction Array
1164 c -------------------------------
1165 do i = 1,lenc
1166 cloudn(i) = 0.0
1167 enddo
1168
1169 CALL CLOUD(nn,LEN, LENC, K, NLTOP, nlayr, IB, RASBLF,SETRAS,FRAC
1170 *, CP, ALHL, RKAPPA, GRAV, P00, CRTMSF
1171 *, POI, QOI, UOI, Ntracer, PRS, PRJ
1172 *, PCU, CLOUDN, TCU, QCU, UCU, CMASS
1173 *, ALF, BET, GAM, PRH, PRI, HOI, ETA
1174 *, HST, QOL, GMH
1175 *, TX1, TX2, TX3, TX4, TX5, TX6, TX7, TX8, TX9
1176 *, WFN, TX11, TX12, TX13, TX14, TX15
1177 *, IA1,IA2,IA3,rhfrac)
1178
1179 C Compute fraction of grid box into which rain re-evap occurs (clf)
1180 c -----------------------------------------------------------------
1181 do i = 1,lenc
1182
1183 c mass in detrainment layer
1184 c -------------------------
1185 tx1(i) = cmb2pa * (prs(i,ib+1) - prs(i,ib))/(grav*dt)
1186
1187 c ratio of detraining cloud mass to mass in detrainment layer
1188 c -----------------------------------------------------------
1189 tx1(i) = rhfrac(i)*rknob * cmass(i,ib) / tx1(i)
1190 if(cmass(i,K).gt.0.) clf(i,ib) = clf(i,ib) + tx1(i)
1191 if( clf(i,ib).gt.1.) clf(i,ib) = 1.
1192 enddo
1193
1194 c Compute Total Cloud Mass Flux
1195 c *****************************
1196 do L=ib,k
1197 do i=1,lenc
1198 cmass(i,L) = rhfrac(i)*cmass(i,L) * dt
1199 enddo
1200 enddo
1201
1202 do L=ib,k
1203 do i=1,lenc
1204 cldmas(i,L) = cldmas(i,L) + cmass(i,L)
1205 enddo
1206 enddo
1207
1208 do i=1,lenc
1209 detrain(i,ib) = detrain(i,ib) + cmass(i,ib)
1210 enddo
1211
1212 DO L=IB,K
1213 DO I=1,LENC
1214 POI(I,L) = POI(I,L) + TCU(I,L) * DT * rhfrac(i)
1215 QOI(I,L) = QOI(I,L) + QCU(I,L) * DT * rhfrac(i)
1216 ENDDO
1217 ENDDO
1218 DO NT=1,Ntracer
1219 DO L=IB,K
1220 DO I=1,LENC
1221 UOI(I,L+nltop-1,NT)=UOI(I,L+nltop-1,NT)+UCU(I,L,NT)*DT*rhfrac(i)
1222 ENDDO
1223 ENDDO
1224 ENDDO
1225 DO I=1,LENC
1226 rains(I,ib) = rains(I,ib) + PCU(I)*dt * rhfrac(i)
1227 ENDDO
1228
1229 100 CONTINUE
1230
1231 c Fill Convective Cloud Fractions based on 3-D Rain Amounts
1232 c ---------------------------------------------------------
1233 do L=k-1,1,-1
1234 do i=1,lenc
1235 tx1(i) = 100*(prs(i,L+1)-prs(i,L))/grav
1236 cln(i,L) = min(1600*rains(i,L)/tx1(i),rasmax )
1237 enddo
1238 enddo
1239
1240 RETURN
1241 END
1242
1243 subroutine rndcloud (iras,nrnd,rnd,myid)
1244 implicit none
1245 integer n,iras,nrnd,myid
1246 real random_numbx
1247 real rnd(nrnd)
1248 integer irm
1249 parameter (irm = 1000)
1250 real random(irm)
1251 integer i,mcheck,numrand,iseed,index
1252 logical first
1253 data first /.true./
1254 integer iras0
1255 data iras0 /0/
1256 save random, iras0
1257
1258 if(nrnd.eq.0.)then
1259 do i = 1,nrnd
1260 rnd(i) = 0
1261 enddo
1262 if(first .and. myid.eq.0) print *,' NO RANDOM CLOUDS IN RAS '
1263 go to 100
1264 endif
1265
1266 mcheck = mod(iras-1,irm/nrnd)
1267
1268 c First Time In From a Continuing RESTART (IRAS.GT.1) or Reading a New RESTART
1269 c ----------------------------------------------------------------------------
1270 if( first.and.(iras.gt.1) .or. iras.ne.iras0+1 )then
1271 if( myid.eq.0 ) print *, 'Recreating Rand Numb Array in RNDCLOUD'
1272 if( myid.eq.0 ) print *, 'IRAS: ',iras,' IRAS0: ',iras0
1273 numrand = mod(iras,irm/nrnd) * nrnd
1274 iseed = iras * nrnd - numrand
1275 call random_seedx(iseed)
1276 do i = 1,irm
1277 random(i) = random_numbx()
1278 enddo
1279 index = (iras-1)*nrnd
1280
1281 c Multiple Time In But have Used Up all 1000 numbers (MCHECK.EQ.0)
1282 c ----------------------------------------------------------------
1283 else if (mcheck.eq.0) then
1284 iseed = (iras-1)*nrnd
1285 call random_seedx(iseed)
1286 do i = 1,irm
1287 random(i) = random_numbx()
1288 enddo
1289 index = iseed
1290
1291 c Multiple Time In But have NOT Used Up all 1000 numbers (MCHECK.NE.0)
1292 c --------------------------------------------------------------------
1293 else
1294 index = (iras-1)*nrnd
1295 endif
1296
1297 index = mod(index,irm)
1298 if( index+nrnd.gt.1000 ) index=1000-nrnd
1299
1300 do n = 1,nrnd
1301 rnd(n) = random(index+n)
1302 enddo
1303
1304 100 continue
1305 first = .false.
1306 iras0 = iras
1307 return
1308 end
1309
1310 real function random_numbx()
1311 implicit none
1312 #if CRAY
1313 real ranf
1314 random_numbx = ranf()
1315 #endif
1316 #if SGI
1317 real rand
1318 random_numbx = rand()
1319 #endif
1320 return
1321 end
1322 subroutine random_seedx (iseed)
1323 implicit none
1324 integer iseed
1325 #if CRAY
1326 call ranset (iseed)
1327 #endif
1328 #if SGI
1329 integer*4 seed
1330 seed = iseed
1331 call srand (seed)
1332 #endif
1333 return
1334 end
1335
1336 SUBROUTINE CLOUD(nn,LEN, LENC, K, NLTOP, nlayr, IC, RASALF,
1337 *, SETRAS, FRAC
1338 *, CP, ALHL, RKAP, GRAV, P00, CRTMSF
1339 *, POI, QOI, UOI, Ntracer, PRS, PRJ
1340 *, PCU, CLN, TCU, QCU, UCU, CMASS
1341 *, ALF, BET, GAM, PRH, PRI, HOL, ETA
1342 *, HST, QOL, GMH
1343 *, TX1, TX2, TX3, TX4, TX5, TX6, TX7, TX8, ALM
1344 *, WFN, AKM, QS1, CLF, UHT, WLQ
1345 *, IA, I1, I2,rhfrac)
1346 C
1347 C*********************************************************************
1348 C******************** Relaxed Arakawa-Schubert ***********************
1349 C********************* Plug Compatible Version **********************
1350 C************************ SUBROUTINE CLOUD ***************************
1351 C************************* 23 JULY 1992 ***************************
1352 C*********************************************************************
1353 C*********************************************************************
1354 C*********************************************************************
1355 C************************** Developed By *****************************
1356 C************************** *****************************
1357 C************************ Shrinivas Moorthi **************************
1358 C************************ and **************************
1359 C************************ Max J. Suarez *****************************
1360 C************************ *****************************
1361 C******************** Laboratory for Atmospheres *********************
1362 C****************** NASA/GSFC, Greenbelt, MD 20771 *******************
1363 C*********************************************************************
1364 C*********************************************************************
1365 C
1366 C The calculations of Moorthi and Suarez (1992, MWR) are
1367 C contained in the CLOUD routine.
1368 C It is probably advisable, at least initially, to treat CLOUD
1369 C as a black box that computes the single cloud adjustments. RAS,
1370 C on the other hand, can be tailored to each GCMs configuration
1371 C (ie, number and placement of levels, nature of boundary layer,
1372 C time step and frequency with which RAS is called).
1373 C
1374 C
1375 C Input:
1376 C ------
1377 C
1378 C LEN : The inner dimension of update and input arrays.
1379 C
1380 C LENC : The run: the number of soundings processes in a single call.
1381 C RAS works on the first LENC of the LEN soundings
1382 C passed. This allows working on pieces of the world
1383 C say for multitasking, without declaring temporary arrays
1384 C and copying the data to and from them. This is an f77
1385 C version. An F90 version would have to allow more
1386 C flexibility in the argument declarations. Obviously
1387 C (LENC<=LEN).
1388 C
1389 C K : Number of vertical layers (increasing downwards).
1390 C Need not be the same as the number of layers in the
1391 C GCM, since it is the outer dimension. The bottom layer
1392 C (K) is the subcloud layer.
1393 C
1394 C IC : Detrainment level to check for presence of convection
1395 C
1396 C RASALF : Relaxation parameter (< 1.) for present cloud-type
1397 C
1398 C SETRAS : Logical parameter to control re-calculation of
1399 C saturation specific humidity and mid level P**kappa
1400 C
1401 C FRAC : Fraction of the PBL (layer K) mass allowed to be used
1402 C by a cloud-type in time DT
1403 C
1404 C CP : Specific heat at constant pressure
1405 C
1406 C ALHL : Latent Heat of condensation
1407 C
1408 C RKAP : R/Cp, where R is the gas constant
1409 C
1410 C GRAV : Acceleration due to gravity
1411 C
1412 C P00 : A reference pressure in hPa, useually 1000 hPa
1413 C
1414 C CRTMSF : Critical value of mass flux above which cloudiness at
1415 C the detrainment layer of that cloud-type is assumed.
1416 C Affects only cloudiness calculation.
1417 C
1418 C POI : 2D array of dimension (LEN,K) containing potential
1419 C temperature. Updated but not initialized by RAS.
1420 C
1421 C QOI : 2D array of dimension (LEN,K) containing specific
1422 C humidity. Updated but not initialized by RAS.
1423 C
1424 C UOI : 3D array of dimension (LEN,K,NTRACER) containing tracers
1425 C Updated but not initialized by RAS.
1426 C
1427 C PRS : 2D array of dimension (LEN,K+1) containing pressure
1428 C in hPa at the interfaces of K-layers from top of the
1429 C atmosphere to the bottom. Not modified.
1430 C
1431 C PRJ : 2D array of dimension (LEN,K+1) containing (PRS/P00) **
1432 C RKAP. i.e. Exner function at layer edges. Not modified.
1433 C
1434 C rhfrac : 1D array of dimension (LEN) containing a rel.hum. scaling
1435 C fraction. Not modified.
1436 C
1437 C Output:
1438 C -------
1439 C
1440 C PCU : 1D array of length LEN containing accumulated
1441 C precipitation in mm/sec.
1442 C
1443 C CLN : 2D array of dimension (LEN,K) containing cloudiness
1444 C Note: CLN is bumped but NOT initialized
1445 C
1446 C TCU : 2D array of dimension (LEN,K) containing accumulated
1447 C convective heating (K/sec).
1448 C
1449 C QCU : 2D array of dimension (LEN,K) containing accumulated
1450 C convective drying (kg/kg/sec).
1451 C
1452 C CMASS : 2D array of dimension (LEN,K) containing the
1453 C cloud mass flux (kg/sec). Filled from cloud top
1454 C to base.
1455 C
1456 C Temporaries:
1457 C
1458 C ALF, BET, GAM, ETA, PRH, PRI, HOI, HST, QOL, GMH are temporary
1459 C 2D real arrays of dimension of at least (LENC,K) where LENC is
1460 C the horizontal dimension over which convection is invoked.
1461 C
1462 C
1463 C TX1, TX2, TX3, TX4, TX5, TX6, TX7, TX8, TX9, AKM, QS1, CLF, UHT
1464 C VHT, WLQ WFN are temporary real arrays of length at least LENC
1465 C
1466 C IA, I1, and I2 are temporary integer arrays of length LENC
1467 C
1468 C
1469 C************************************************************************
1470 C
1471 C
1472
1473 PARAMETER (DAYLEN=86400.0, HALF=0.5, ONE=1.0, ZERO=0.0)
1474 PARAMETER (CMB2PA=100.0)
1475 PARAMETER (RHMAX=0.9999)
1476 C
1477 integer nltop,ntracer,nlayr
1478 DIMENSION POI(LEN,K), QOI(LEN,K), PRS(LEN,K+1)
1479 *, PRJ(LEN,K+1)
1480 *, TCU(LEN,K), QCU(LEN,K), CMASS(LEN,K), CLN(LEN)
1481 real uoi(len,nlayr,ntracer)
1482 DIMENSION ALF(LEN,K), BET(LEN,K), GAM(LEN,K)
1483 *, PRH(LEN,K), PRI(LEN,K)
1484 DIMENSION AKM(LENC), WFN(LENC)
1485 DIMENSION HOL(LENC,K), QOL(LENC,K), ETA(LENC,K), HST(LENC,K)
1486 *, GMH(LENC,K), ALM(LENC), WLQ(LENC), QS1(LENC)
1487 *, TX1(LENC), TX2(LENC), TX3(LENC), TX4(LENC)
1488 *, TX5(LENC), TX6(LENC), TX7(LENC), TX8(LENC)
1489 *, CLF(LENC), PCU(LENC)
1490 DIMENSION IA(LENC), I1(LENC), I2(LENC)
1491 real rhfrac(len)
1492 real ucu(len,k,ntracer),uht(len,ntracer)
1493 LOGICAL SETRAS
1494
1495 integer nt
1496
1497 c Explicit Inline Directives
1498 c --------------------------
1499 #if CRAY
1500 #if f77
1501 cfpp$ expand (qsat)
1502 #endif
1503 #endif
1504
1505 RKAPP1 = 1.0 + RKAP
1506 ONEBCP = 1.0 / CP
1507 ALBCP = ALHL * ONEBCP
1508 ONEBG = 1.0 / GRAV
1509 CPBG = CP * ONEBG
1510 TWOBAL = 2.0 / ALHL
1511 C
1512 KM1 = K - 1
1513 IC1 = IC + 1
1514 C
1515 C SETTIING ALF, BET, GAM, PRH, AND PRI : DONE ONLY WHEN SETRAS=.T.
1516 C
1517
1518 IF (SETRAS) THEN
1519
1520 DO 2050 L=1,K
1521 DO 2030 I=1,LENC
1522 PRH(I,L) = (PRJ(I,L+1)*PRS(I,L+1) - PRJ(I,L)*PRS(I,L))
1523 * / ((PRS(I,L+1)-PRS(I,L)) * RKAPP1)
1524 2030 CONTINUE
1525 2050 CONTINUE
1526
1527 DO 2070 L=1,K
1528 DO 2060 I=1,LENC
1529 TX5(I) = POI(I,L) * PRH(I,L)
1530 TX1(I) = (PRS(I,L) + PRS(I,L+1)) * 0.5
1531 TX3(I) = TX5(I)
1532 CALL QSAT(TX3(I), TX1(I), TX2(I), TX4(I), .TRUE.)
1533 ALF(I,L) = TX2(I) - TX4(I) * TX5(I)
1534 BET(I,L) = TX4(I) * PRH(I,L)
1535 GAM(I,L) = 1.0 / ((1.0 + TX4(I)*ALBCP) * PRH(I,L))
1536 PRI(I,L) = (CP/CMB2PA) / (PRS(I,L+1) - PRS(I,L))
1537 2060 CONTINUE
1538 2070 CONTINUE
1539
1540 ENDIF
1541 C
1542 C
1543 DO 10 L=1,K
1544 DO 10 I=1,LEN
1545 TCU(I,L) = 0.0
1546 QCU(I,L) = 0.0
1547 CMASS(I,L) = 0.0
1548 10 CONTINUE
1549
1550 do nt = 1,ntracer
1551 do L=1,K
1552 do I=1,LENC
1553 ucu(I,L,nt) = 0.0
1554 enddo
1555 enddo
1556 enddo
1557 C
1558 DO 30 I=1,LENC
1559 TX1(I) = PRJ(I,K+1) * POI(I,K)
1560 QS1(I) = ALF(I,K) + BET(I,K)*POI(I,K)
1561 QOL(I,K) = MIN(QS1(I)*RHMAX,QOI(I,K))
1562
1563 HOL(I,K) = TX1(I)*CP + QOL(I,K)*ALHL
1564 ETA(I,K) = ZERO
1565 TX2(I) = (PRJ(I,K+1) - PRJ(I,K)) * POI(I,K) * CP
1566 30 CONTINUE
1567 C
1568 IF (IC .LT. KM1) THEN
1569 DO 3703 L=KM1,IC1,-1
1570 DO 50 I=1,LENC
1571 QS1(I) = ALF(I,L) + BET(I,L)*POI(I,L)
1572 QOL(I,L) = MIN(QS1(I)*RHMAX,QOI(I,L))
1573 C
1574 TEM1 = TX2(I) + PRJ(I,L+1) * POI(I,L) * CP
1575 HOL(I,L) = TEM1 + QOL(I,L )* ALHL
1576 HST(I,L) = TEM1 + QS1(I) * ALHL
1577
1578 TX1(I) = (PRJ(I,L+1) - PRJ(I,L)) * POI(I,L)
1579 ETA(I,L) = ETA(I,L+1) + TX1(I)*CPBG
1580 TX2(I) = TX2(I) + TX1(I)*CP
1581 50 CONTINUE
1582 C
1583 3703 CONTINUE
1584 ENDIF
1585
1586
1587 DO 70 I=1,LENC
1588 HOL(I,IC) = TX2(I)
1589 QS1(I) = ALF(I,IC) + BET(I,IC)*POI(I,IC)
1590 QOL(I,IC) = MIN(QS1(I)*RHMAX,QOI(I,IC))
1591 c
1592 TEM1 = TX2(I) + PRJ(I,IC1) * POI(I,IC) * CP
1593 HOL(I,IC) = TEM1 + QOL(I,IC) * ALHL
1594 HST(I,IC) = TEM1 + QS1(I) * ALHL
1595 C
1596 TX3(I ) = (PRJ(I,IC1) - PRH(I,IC)) * POI(I,IC)
1597 ETA(I,IC) = ETA(I,IC1) + CPBG * TX3(I)
1598 70 CONTINUE
1599 C
1600 DO 130 I=1,LENC
1601 TX2(I) = HOL(I,K) - HST(I,IC)
1602 TX1(I) = ZERO
1603
1604 130 CONTINUE
1605 C
1606 C ENTRAINMENT PARAMETER
1607 C
1608 DO 160 L=IC,KM1
1609 DO 160 I=1,LENC
1610 TX1(I) = TX1(I) + (HST(I,IC) - HOL(I,L)) * (ETA(I,L) - ETA(I,L+1))
1611 160 CONTINUE
1612 C
1613 LEN1 = 0
1614 LEN2 = 0
1615 ISAV = 0
1616 DO 195 I=1,LENC
1617 IF (TX1(I) .GT. ZERO .AND. TX2(I) .GT. ZERO
1618 . .AND. rhfrac(i).ne.0.0 ) THEN
1619 LEN1 = LEN1 + 1
1620 IA(LEN1) = I
1621 ALM(LEN1) = TX2(I) / TX1(I)
1622 ENDIF
1623 195 CONTINUE
1624 C
1625 LEN2 = LEN1
1626 if (IC1 .lt. K) then
1627 DO 196 I=1,LENC
1628 IF (TX2(I) .LE. 0.0 .AND. (HOL(I,K) .GT. HST(I,IC1))
1629 . .AND. rhfrac(i).ne.0.0 ) THEN
1630 LEN2 = LEN2 + 1
1631 IA(LEN2) = I
1632 ALM(LEN2) = 0.0
1633 ENDIF
1634 196 CONTINUE
1635 endif
1636 C
1637 IF (LEN2 .EQ. 0) THEN
1638 DO 5010 I=1,LENC*K
1639 HST(I,1) = 0.0
1640 QOL(I,1) = 0.0
1641 5010 CONTINUE
1642 DO 5020 I=1,LENC
1643 PCU(I) = 0.0
1644 5020 CONTINUE
1645 RETURN
1646 ENDIF
1647 LEN11 = LEN1 + 1
1648 C
1649 C NORMALIZED MASSFLUX
1650 C
1651 DO 250 I=1,LEN2
1652 ETA(I,K) = 1.0
1653 II = IA(I)
1654 TX2(I) = 0.5 * (PRS(II,IC) + PRS(II,IC1))
1655 TX4(I) = PRS(II,K)
1656 250 CONTINUE
1657 C
1658 DO 252 I=LEN11,LEN2
1659 WFN(I) = 0.0
1660 II = IA(I)
1661 IF (HST(II,IC1) .LT. HST(II,IC)) THEN
1662 TX6(I) = (HST(II,IC1)-HOL(II,K))/(HST(II,IC1)-HST(II,IC))
1663 ELSE
1664 TX6(I) = 0.0
1665 ENDIF
1666 TX2(I) = 0.5 * (PRS(II,IC1)+PRS(II,IC1+1)) * (1.0-TX6(I))
1667 * + TX2(I) * TX6(I)
1668 252 CONTINUE
1669 C
1670 CALL ACRITN(LEN2, TX2, TX4, TX3)
1671 C
1672 DO 260 L=KM1,IC,-1
1673 DO 255 I=1,LEN2
1674 TX1(I) = ETA(IA(I),L)
1675 255 CONTINUE
1676 DO 260 I=1,LEN2
1677 ETA(I,L) = 1.0 + ALM(I) * TX1(I)
1678 260 CONTINUE
1679 C
1680 C CLOUD WORKFUNCTION
1681 C
1682 IF (LEN1 .GT. 0) THEN
1683 DO 270 I=1,LEN1
1684 II = IA(I)
1685 WFN(I) = - GAM(II,IC) * (PRJ(II,IC1) - PRH(II,IC))
1686 * * HST(II,IC) * ETA(I,IC1)
1687 270 CONTINUE
1688 ENDIF
1689 C
1690 DO 290 I=1,LEN2
1691 II = IA(I)
1692 TX1(I) = HOL(II,K)
1693 290 CONTINUE
1694 C
1695 IF (IC1 .LE. KM1) THEN
1696
1697 DO 380 L=KM1,IC1,-1
1698 DO 380 I=1,LEN2
1699 II = IA(I)
1700 TEM = TX1(I) + (ETA(I,L) - ETA(I,L+1)) * HOL(II,L)
1701 C
1702 PCU(I) = PRJ(II,L+1) - PRH(II,L)
1703 TEM1 = ETA(I,L+1) * PCU(I)
1704 TX1(I) = TX1(I)*PCU(I)
1705 C
1706 PCU(I) = PRH(II,L) - PRJ(II,L)
1707 TEM1 = (TEM1 + ETA(I,L) * PCU(I)) * HST(II,L)
1708 TX1(I) = TX1(I) + TEM*PCU(I)
1709 C
1710 WFN(I) = WFN(I) + (TX1(I) - TEM1) * GAM(II,L)
1711 TX1(I) = TEM
1712 380 CONTINUE
1713 ENDIF
1714 C
1715 LENA = 0
1716 IF (LEN1 .GT. 0) THEN
1717 DO 512 I=1,LEN1
1718 II = IA(I)
1719 WFN(I) = WFN(I) + TX1(I) * GAM(II,IC)*(PRJ(II,IC1)-PRH(II,IC))
1720 * - TX3(I)
1721 IF (WFN(I) .GT. 0.0) THEN
1722 LENA = LENA + 1
1723 I1(LENA) = IA(I)
1724 I2(LENA) = I
1725 TX1(LENA) = WFN(I)
1726 TX2(LENA) = QS1(IA(I))
1727 TX6(LENA) = 1.0
1728 ENDIF
1729 512 CONTINUE
1730 ENDIF
1731 LENB = LENA
1732 DO 515 I=LEN11,LEN2
1733 WFN(I) = WFN(I) - TX3(I)
1734 IF (WFN(I) .GT. 0.0 .AND. TX6(I) .GT. 0.0) THEN
1735 LENB = LENB + 1
1736 I1(LENB) = IA(I)
1737 I2(LENB) = I
1738 TX1(LENB) = WFN(I)
1739 TX2(LENB) = QS1(IA(I))
1740 TX4(LENB) = TX6(I)
1741 ENDIF
1742 515 CONTINUE
1743 C
1744 IF (LENB .LE. 0) THEN
1745 DO 5030 I=1,LENC*K
1746 HST(I,1) = 0.0
1747 QOL(I,1) = 0.0
1748 5030 CONTINUE
1749 DO 5040 I=1,LENC
1750 PCU(I) = 0.0
1751 5040 CONTINUE
1752 RETURN
1753 ENDIF
1754
1755 C
1756 DO 516 I=1,LENB
1757 WFN(I) = TX1(I)
1758 QS1(I) = TX2(I)
1759 516 CONTINUE
1760 C
1761 DO 520 L=IC,K
1762 DO 517 I=1,LENB
1763 TX1(I) = ETA(I2(I),L)
1764 517 CONTINUE
1765 DO 520 I=1,LENB
1766 ETA(I,L) = TX1(I)
1767 520 CONTINUE
1768 C
1769 LENA1 = LENA + 1
1770 C
1771 DO 510 I=1,LENA
1772 II = I1(I)
1773 TX8(I) = HST(II,IC) - HOL(II,IC)
1774 510 CONTINUE
1775 DO 530 I=LENA1,LENB
1776 II = I1(I)
1777 TX6(I) = TX4(I)
1778 TEM = TX6(I) * (HOL(II,IC)-HOL(II,IC1)) + HOL(II,IC1)
1779 TX8(I) = HOL(II,K) - TEM
1780
1781 TEM1 = TX6(I) * (QOL(II,IC)-QOL(II,IC1)) + QOL(II,IC1)
1782 TX5(I) = TEM - TEM1 * ALHL
1783 QS1(I) = TEM1 + TX8(I)*(ONE/ALHL)
1784 TX3(I) = HOL(II,IC)
1785 530 CONTINUE
1786 C
1787 C
1788 DO 620 I=1,LENB
1789 II = I1(I)
1790 WLQ(I) = QOL(II,K) - QS1(I) * ETA(I,IC)
1791 TX7(I) = HOL(II,K)
1792 620 CONTINUE
1793 DO NT=1,Ntracer
1794 DO 621 I=1,LENB
1795 II = I1(I)
1796 UHT(I,NT) = UOI(II,K+nltop-1,NT)-UOI(II,IC+nltop-1,NT) * ETA(I,IC)
1797 621 CONTINUE
1798 ENDDO
1799 C
1800 DO 635 L=KM1,IC,-1
1801 DO 630 I=1,LENB
1802 II = I1(I)
1803 TEM = ETA(I,L) - ETA(I,L+1)
1804 WLQ(I) = WLQ(I) + TEM * QOL(II,L)
1805 630 CONTINUE
1806 635 CONTINUE
1807 DO NT=1,Ntracer
1808 DO L=KM1,IC,-1
1809 DO I=1,LENB
1810 II = I1(I)
1811 TEM = ETA(I,L) - ETA(I,L+1)
1812 UHT(I,NT) = UHT(I,NT) + TEM * UOI(II,L+nltop-1,NT)
1813 ENDDO
1814 ENDDO
1815 ENDDO
1816 C
1817 C CALCULATE GS AND PART OF AKM (THAT REQUIRES ETA)
1818 C
1819 DO 690 I=1,LENB
1820 II = I1(I)
1821 c TX7(I) = HOL(II,K)
1822 TEM = (POI(II,KM1) - POI(II,K)) / (PRH(II,K) - PRH(II,KM1))
1823 HOL(I,K) = TEM * (PRJ(II,K)-PRH(II,KM1))*PRH(II,K)*PRI(II,K)
1824 HOL(I,KM1) = TEM * (PRH(II,K)-PRJ(II,K))*PRH(II,KM1)*PRI(II,KM1)
1825 AKM(I) = ZERO
1826 TX2(I) = 0.5 * (PRS(II,IC) + PRS(II,IC1))
1827 690 CONTINUE
1828
1829 IF (IC1 .LE. KM1) THEN
1830 DO 750 L=KM1,IC1,-1
1831 DO 750 I=1,LENB
1832 II = I1(I)
1833 TEM = (POI(II,L-1) - POI(II,L)) * ETA(I,L)
1834 * / (PRH(II,L) - PRH(II,L-1))
1835 C
1836 HOL(I,L) = TEM * (PRJ(II,L)-PRH(II,L-1)) * PRH(II,L)
1837 * * PRI(II,L) + HOL(I,L)
1838 HOL(I,L-1) = TEM * (PRH(II,L)-PRJ(II,L)) * PRH(II,L-1)
1839 * * PRI(II,L-1)
1840 C
1841 AKM(I) = AKM(I) - HOL(I,L)
1842 * * (ETA(I,L) * (PRH(II,L)-PRJ(II,L)) +
1843 * ETA(I,L+1) * (PRJ(II,L+1)-PRH(II,L))) / PRH(II,L)
1844 750 CONTINUE
1845 ENDIF
1846 C
1847 C
1848 CALL RNCL(LENB, TX2, TX1, CLF)
1849
1850 DO 770 I=1,LENB
1851 TX2(I) = (ONE - TX1(I)) * WLQ(I)
1852 WLQ(I) = TX1(I) * WLQ(I)
1853 C
1854 TX1(I) = HOL(I,IC)
1855 770 CONTINUE
1856 DO 790 I=LENA1, LENB
1857 II = I1(I)
1858 TX1(I) = TX1(I) + (TX5(I)-TX3(I)+QOL(II,IC)*ALHL)*(PRI(II,IC)/CP)
1859 790 CONTINUE
1860
1861 DO 800 I=1,LENB
1862 HOL(I,IC) = TX1(I) - TX2(I) * ALBCP * PRI(I1(I),IC)
1863 800 CONTINUE
1864
1865 IF (LENA .GT. 0) THEN
1866 DO 810 I=1,LENA
1867 II = I1(I)
1868 AKM(I) = AKM(I) - ETA(I,IC1) * (PRJ(II,IC1) - PRH(II,IC))
1869 * * TX1(I) / PRH(II,IC)
1870 810 CONTINUE
1871 ENDIF
1872 c
1873 C CALCULATE GH
1874 C
1875 DO 830 I=1,LENB
1876 II = I1(I)
1877 TX3(I) = QOL(II,KM1) - QOL(II,K)
1878 GMH(I,K) = HOL(I,K) + TX3(I) * PRI(II,K) * (ALBCP)
1879
1880 AKM(I) = AKM(I) + GAM(II,KM1)*(PRJ(II,K)-PRH(II,KM1))
1881 * * GMH(I,K)
1882 TX3(I) = zero
1883 830 CONTINUE
1884 C
1885 IF (IC1 .LE. KM1) THEN
1886 DO 840 L=KM1,IC1,-1
1887 DO 840 I=1,LENB
1888 II = I1(I)
1889 TX2(I) = TX3(I)
1890 TX3(I) = (QOL(II,L-1) - QOL(II,L)) * ETA(I,L)
1891 TX2(I) = TX2(I) + TX3(I)
1892 C
1893 GMH(I,L) = HOL(I,L) + TX2(I) * PRI(II,L) * (ALBCP*HALF)
1894 840 CONTINUE
1895 C
1896 C
1897 ENDIF
1898 DO 850 I=LENA1,LENB
1899 TX3(I) = TX3(I) + TWOBAL
1900 * * (TX7(I) - TX8(I) - TX5(I) - QOL(I1(I),IC)*ALHL)
1901 850 CONTINUE
1902 DO 860 I=1,LENB
1903 GMH(I,IC) = TX1(I) + PRI(I1(I),IC) * ONEBCP
1904 * * (TX3(I)*(ALHL*HALF) + ETA(I,IC) * TX8(I))
1905 860 CONTINUE
1906 C
1907 C CALCULATE HC PART OF AKM
1908 C
1909 IF (IC1 .LE. KM1) THEN
1910 DO 870 I=1,LENB
1911 TX1(I) = GMH(I,K)
1912 870 CONTINUE
1913 DO 3725 L=KM1,IC1,-1
1914 DO 880 I=1,LENB
1915 II = I1(I)
1916 TX1(I) = TX1(I) + (ETA(I,L) - ETA(I,L+1)) * GMH(I,L)
1917 TX2(I) = GAM(II,L-1) * (PRJ(II,L) - PRH(II,L-1))
1918 880 CONTINUE
1919 C
1920 IF (L .EQ. IC1) THEN
1921 DO 890 I=LENA1,LENB
1922 TX2(I) = ZERO
1923 890 CONTINUE
1924 ENDIF
1925 DO 900 I=1,LENB
1926 II = I1(I)
1927 AKM(I) = AKM(I) + TX1(I) *
1928 * (TX2(I) + GAM(II,L)*(PRH(II,L)-PRJ(II,L)))
1929 900 CONTINUE
1930 3725 CONTINUE
1931 ENDIF
1932 C
1933 DO 920 I=LENA1,LENB
1934 II = I1(I)
1935 TX2(I) = 0.5 * (PRS(II,IC) + PRS(II,IC1))
1936 * + 0.5*(PRS(II,IC+2) - PRS(II,IC)) * (ONE-TX6(I))
1937 c
1938 TX1(I) = PRS(II,IC1)
1939 TX5(I) = 0.5 * (PRS(II,IC1) + PRS(II,IC+2))
1940 C
1941 IF ((TX2(I) .GE. TX1(I)) .AND. (TX2(I) .LT. TX5(I))) THEN
1942 TX6(I) = ONE - (TX2(I) - TX1(I)) / (TX5(I) - TX1(I))
1943 C
1944 TEM = PRI(II,IC1) / PRI(II,IC)
1945 HOL(I,IC1) = HOL(I,IC1) + HOL(I,IC) * TEM
1946 HOL(I,IC) = ZERO
1947 C
1948 GMH(I,IC1) = GMH(I,IC1) + GMH(I,IC) * TEM
1949 GMH(I,IC) = ZERO
1950 ELSEIF (TX2(I) .LT. TX1(I)) THEN
1951 TX6(I) = 1.0
1952 ELSE
1953 TX6(I) = 0.0
1954 ENDIF
1955 920 CONTINUE
1956 C
1957 C
1958 DO I=1,LENC
1959 PCU(I) = 0.0
1960 ENDDO
1961
1962 DO 970 I=1,LENB
1963 II = I1(I)
1964 IF (AKM(I) .LT. ZERO .AND. WLQ(I) .GE. 0.0) THEN
1965 WFN(I) = - TX6(I) * WFN(I) * RASALF / AKM(I)
1966 ELSE
1967 WFN(I) = ZERO
1968 ENDIF
1969 TEM = (PRS(II,K+1)-PRS(II,K))*(CMB2PA*FRAC)
1970 WFN(I) = MIN(WFN(I), TEM)
1971 C
1972 C compute cloud amount
1973 C
1974 CC TX1(I) = CLN(II)
1975 CC IF (WFN(I) .GT. CRTMSF) TX1(I) = TX1(I) + CLF(I)
1976 CC IF (TX1(I) .GT. ONE) TX1(I) = ONE
1977 C
1978 C PRECIPITATION
1979 C
1980 PCU(II) = WLQ(I) * WFN(I) * ONEBG
1981 C
1982 C CUMULUS FRICTION AT THE BOTTOM LAYER
1983 C
1984 TX4(I) = WFN(I) * (1.0/ALHL)
1985 TX5(I) = WFN(I) * ONEBCP
1986 970 CONTINUE
1987 C
1988 C compute cloud mass flux for diagnostic output
1989 C
1990 DO L = IC,K
1991 DO I=1,LENB
1992 II = I1(I)
1993 if(L.lt.K)then
1994 CMASS(II,L) = ETA(I,L+1) * WFN(I) * ONEBG
1995 else
1996 CMASS(II,L) = WFN(I) * ONEBG
1997 endif
1998 ENDDO
1999 ENDDO
2000 C
2001 CC DO 975 I=1,LENB
2002 CC II = I1(I)
2003 CC CLN(II) = TX1(I)
2004 CC975 CONTINUE
2005 C
2006 C THETA AND Q CHANGE DUE TO CLOUD TYPE IC
2007 C
2008
2009 c TEMA = 0.0
2010 c TEMB = 0.0
2011 DO 990 L=IC,K
2012 DO 980 I=1,LENB
2013 II = I1(I)
2014 TEM = (GMH(I,L) - HOL(I,L)) * TX4(I)
2015 TEM1 = HOL(I,L) * TX5(I)
2016 C
2017 TCU(II,L) = TEM1 / PRH(II,L)
2018 QCU(II,L) = TEM
2019 980 CONTINUE
2020
2021 c I = I1(IP1)
2022 c
2023 c TEM = (PRS(I,L+1)-PRS(I,L)) * (ONEBG*100.0)
2024 c TEMA = TEMA + TCU(I,L) * PRH(I,L) * TEM * (CP/ALHL)
2025 c TEMB = TEMB + QCU(I,L) * TEM
2026 C
2027 990 CONTINUE
2028 C
2029 c Compute Tracer Tendencies
2030 c -------------------------
2031 do nt = 1,ntracer
2032
2033 c Tracer Tendency at the Bottom Layer
2034 c -----------------------------------
2035 DO 995 I=1,LENB
2036 II = I1(I)
2037 TEM = half*TX5(I) * PRI(II,K)
2038 TX1(I) = (UOI(II,KM1+nltop-1,nt) - UOI(II,K+nltop-1,nt))
2039 ucu(II,K,nt) = TEM * TX1(I)
2040 995 CONTINUE
2041
2042 c Tracer Tendency at all other Levels
2043 c -----------------------------------
2044 DO 1020 L=KM1,IC1,-1
2045 DO 1010 I=1,LENB
2046 II = I1(I)
2047 TEM = half*TX5(I) * PRI(II,L)
2048 TEM1 = TX1(I)
2049 TX1(I) = (UOI(II,L-1+nltop-1,nt)-UOI(II,L+nltop-1,nt)) * ETA(I,L)
2050 TX3(I) = (TX1(I) + TEM1) * TEM
2051 1010 CONTINUE
2052 DO 1020 I=1,LENB
2053 II = I1(I)
2054 ucu(II,L,nt) = TX3(I)
2055 1020 CONTINUE
2056
2057 DO 1030 I=1,LENB
2058 II = I1(I)
2059 IF (TX6(I) .GE. 1.0) THEN
2060 TEM = half*TX5(I) * PRI(II,IC)
2061 ELSE
2062 TEM = 0.0
2063 ENDIF
2064 TX1(I) = (TX1(I) + UHT(I,nt) + UHT(I,nt)) * TEM
2065 1030 CONTINUE
2066 DO 1040 I=1,LENB
2067 II = I1(I)
2068 ucu(II,IC,nt) = TX1(I)
2069 1040 CONTINUE
2070
2071 enddo
2072 C
2073 C PENETRATIVE CONVECTION CALCULATION OVER
2074 C
2075
2076 RETURN
2077 END
2078 SUBROUTINE RNCL(LEN, PL, RNO, CLF)
2079 C
2080 C
2081 C*********************************************************************
2082 C********************** Relaxed Arakawa-Schubert *********************
2083 C************************ SUBROUTINE RNCL ************************
2084 C**************************** 23 July 1992 ***************************
2085 C*********************************************************************
2086
2087 PARAMETER (P5=500.0, P8=800.0, PT8=0.8, PT2=0.2)
2088 PARAMETER (PFAC=PT2/(P8-P5))
2089 C
2090 PARAMETER (P4=400.0, P6=401.0)
2091 PARAMETER (P7=700.0, P9=900.0)
2092 PARAMETER (CUCLD=0.5,CFAC=CUCLD/(P6-P4))
2093 C
2094 DIMENSION PL(LEN), RNO(LEN), CLF(LEN)
2095
2096 DO 10 I=1,LEN
2097 rno(i) = 1.0
2098 ccc if( pl(i).le.400.0 ) rno(i) = max( 0.75, 1.0-0.0025*(400.0-pl(i)) )
2099
2100 ccc IF ( PL(I).GE.P7 .AND. PL(I).LE.P9 ) THEN
2101 ccc RNO(I) = ((P9-PL(I))/(P9-P7)) **2
2102 ccc ELSE IF (PL(I).GT.P9) THEN
2103 ccc RNO(I) = 0.
2104 ccc ENDIF
2105
2106 CLF(I) = CUCLD
2107 C
2108 CARIESIF (PL(I) .GE. P5 .AND. PL(I) .LE. P8) THEN
2109 CARIES RNO(I) = (P8-PL(I))*PFAC + PT8
2110 CARIESELSEIF (PL(I) .GT. P8 ) THEN
2111 CARIES RNO(I) = PT8
2112 CARIESENDIF
2113 CARIES
2114 IF (PL(I) .GE. P4 .AND. PL(I) .LE. P6) THEN
2115 CLF(I) = (P6-PL(I))*CFAC
2116 ELSEIF (PL(I) .GT. P6 ) THEN
2117 CLF(I) = 0.0
2118 ENDIF
2119 10 CONTINUE
2120 C
2121 RETURN
2122 END
2123 SUBROUTINE ACRITN ( LEN,PL,PLB,ACR )
2124
2125 C*********************************************************************
2126 C********************** Relaxed Arakawa-Schubert *********************
2127 C************************** SUBROUTINE ACRIT *********************
2128 C****************** modified August 28, 1996 L.Takacs ************
2129 C**** *****
2130 C**** Note: Data obtained from January Mean After-Analysis *****
2131 C**** from 4x5 46-layer GEOS Assimilation *****
2132 C**** *****
2133 C*********************************************************************
2134
2135 real PL(LEN), PLB(LEN), ACR(LEN)
2136
2137 parameter (lma=18)
2138 real p(lma)
2139 real a(lma)
2140
2141 data p / 93.81, 111.65, 133.46, 157.80, 186.51,
2142 . 219.88, 257.40, 301.21, 352.49, 409.76,
2143 . 471.59, 535.04, 603.33, 672.79, 741.12,
2144 . 812.52, 875.31, 930.20/
2145
2146 data a / 3.35848, 3.13645, 2.48072, 2.08277, 1.53364,
2147 . 1.01971, .65846, .45867, .38687, .31002,
2148 . .25574, .20347, .17254, .15260, .16756,
2149 . .09916, .10360, .05880/
2150
2151
2152 do L=1,lma-1
2153 do i=1,len
2154 if( pl(i).ge.p(L) .and.
2155 . pl(i).le.p(L+1)) then
2156 temp = ( pl(i)-p(L) )/( p(L+1)-p(L) )
2157 acr(i) = a(L+1)*temp + a(L)*(1-temp)
2158 endif
2159 enddo
2160 enddo
2161
2162 do i=1,len
2163 if( pl(i).lt.p(1) ) acr(i) = a(1)
2164 if( pl(i).gt.p(lma) ) acr(i) = a(lma)
2165 enddo
2166
2167 do i=1,len
2168 acr(i) = acr(i) * (plb(i)-pl(i))
2169 enddo
2170
2171 RETURN
2172 END
2173 SUBROUTINE RNEVP(NN,IRUN,NLAY,TL,QL,RAIN,PL,CLFRAC,SP,DP,PLKE,
2174 1 PLK,TH,TEMP1,TEMP2,TEMP3,ITMP1,ITMP2,RCON,RLAR,CLSBTH,tmscl,
2175 2 tmfrc,cp,gravity,alhl,gamfac,cldlz,RHCRIT,offset,alpha)
2176
2177 PARAMETER (ZM1P04 = -1.04E-4 )
2178 PARAMETER (ZERO = 0.)
2179 PARAMETER (TWO89= 2.89E-5)
2180 PARAMETER ( ZP44= 0.44)
2181 PARAMETER ( ZP01= 0.01)
2182 PARAMETER ( ZP1 = 0.1 )
2183 PARAMETER ( ZP001= 0.001)
2184 PARAMETER ( HALF= 0.5)
2185 PARAMETER ( ZP578 = 0.578 )
2186 PARAMETER ( ONE = 1.)
2187 PARAMETER ( THOUSAND = 1000.)
2188 PARAMETER ( Z3600 = 3600.)
2189 C
2190 DIMENSION TL(IRUN,NLAY),QL(IRUN,NLAY),RAIN(IRUN,NLAY),
2191 $ PL(IRUN,NLAY),CLFRAC(IRUN,NLAY),SP(IRUN),TEMP1(IRUN,NLAY),
2192 $ TEMP2(IRUN,NLAY),PLKE(IRUN,NLAY+1),
2193 $ RCON(IRUN),RLAR(IRUN),DP(IRUN,NLAY),PLK(IRUN,NLAY),TH(IRUN,NLAY),
2194 $ TEMP3(IRUN,NLAY),ITMP1(IRUN,NLAY),
2195 $ ITMP2(IRUN,NLAY),CLSBTH(IRUN,NLAY)
2196 C
2197 DIMENSION EVP9(IRUN,NLAY)
2198 real water(irun),crystal(irun)
2199 real watevap(irun),iceevap(irun)
2200 real fracwat,fracice, tice,rh,fact,dum
2201
2202 real cldlz(irun,nlay)
2203 real rhcrit(irun,nlay), rainmax(irun)
2204 real offset, alpha
2205
2206 c Explicit Inline Directives
2207 c --------------------------
2208 #if CRAY
2209 #if f77
2210 cfpp$ expand (qsat)
2211 #endif
2212 #endif
2213
2214 tice = getcon('FREEZING-POINT')
2215
2216 fracwat = 0.70
2217 fracice = 0.01
2218
2219 NLAYM1 = NLAY - 1
2220 IRNLAY = IRUN*NLAY
2221 IRNLM1 = IRUN*(NLAY-1)
2222
2223 RPHF = Z3600/tmscl
2224
2225 ELOCP = alhl/cp
2226 CPOG = cp/gravity
2227
2228 DO I = 1,IRUN
2229 RLAR(I) = 0.
2230 water(i) = 0.
2231 crystal(i) = 0.
2232 ENDDO
2233
2234 do L = 1,nlay
2235 do i = 1,irun
2236 EVP9(i,L) = 0.
2237 TEMP1(i,L) = 0.
2238 TEMP2(i,L) = 0.
2239 TEMP3(i,L) = 0.
2240 CLSBTH(i,L) = 0.
2241 cldlz(i,L) = 0.
2242 enddo
2243 enddo
2244
2245 C RHO(ZERO) / RHO FOR TERMINAL VELOCITY APPROX.
2246 c ---------------------------------------------
2247 DO L = 1,NLAY
2248 DO I = 1,IRUN
2249 TEMP2(I,L) = PL(I,L)*ZP001
2250 TEMP2(I,L) = SQRT(TEMP2(I,L))
2251 ENDDO
2252 ENDDO
2253
2254 C INVERSE OF MASS IN EACH LAYER
2255 c -----------------------------
2256 DO L = 1,NLAY
2257 DO I = 1,IRUN
2258 TEMP3(I,L) = GRAVITY*ZP01 / DP(I,L)
2259 ENDDO
2260 ENDDO
2261
2262 C DO LOOP FOR MOISTURE EVAPORATION ABILITY AND CONVEC EVAPORATION.
2263 c ----------------------------------------------------------------
2264 DO 100 L=1,NLAY
2265
2266 DO I = 1,IRUN
2267 TEMP1(I,3) = TL(I,L)
2268 TEMP1(I,4) = QL(I,L)
2269 ENDDO
2270
2271 DO 50 N=1,2
2272 IF(N.EQ.1)RELAX=HALF
2273 IF(N.GT.1)RELAX=ONE
2274
2275 DO I = 1,IRUN
2276 call qsat ( temp1(i,3),pl(i,L),temp1(i,2),temp1(i,6),.true. )
2277 TEMP1(I,5)=TEMP1(I,2)-TEMP1(I,4)
2278 TEMP1(I,6)=TEMP1(I,6)*ELOCP
2279 TEMP1(I,5)=TEMP1(I,5)/(ONE+TEMP1(I,6))
2280 TEMP1(I,4)=TEMP1(I,4)+TEMP1(I,5)*RELAX
2281 TEMP1(I,3)=TEMP1(I,3)-TEMP1(I,5)*ELOCP*RELAX
2282 ENDDO
2283 50 CONTINUE
2284
2285 DO I = 1,IRUN
2286 EVP9(I,L) = (TEMP1(I,4) - QL(I,L))/TEMP3(I,L)
2287 C convective detrained water
2288 cldlz(i,L) = rain(i,L)*temp3(i,L)
2289 if( tl(i,L).gt.tice-20.) then
2290 water(i) = water(i) + rain(i,L)
2291 else
2292 crystal(i) = crystal(i) + rain(i,L)
2293 endif
2294 ENDDO
2295
2296 C**********************************************************************
2297 C FOR CONVECTIVE PRECIP, FIND THE "EVAPORATION EFFICIENCY" USING *
2298 C KESSLERS PARAMETERIZATION *
2299 C**********************************************************************
2300
2301 DO 20 I=1,IRUN
2302
2303 iceevap(i) = 0.
2304 watevap(i) = 0.
2305
2306 if( (evp9(i,L).gt.0.) .and. (crystal(i).gt.0.) ) then
2307 iceevap(I) = EVP9(I,L)*fracice
2308 IF(iceevap(i).GE.crystal(i)) iceevap(i) = crystal(i)
2309 EVP9(I,L)=EVP9(I,L)-iceevap(I)
2310 crystal(i) = crystal(i) - iceevap(i)
2311 endif
2312
2313 C and now warm precipitate
2314 if( (evp9(i,L).gt.0.) .and. (water(i).gt.0.) ) then
2315 exparg = ZM1P04*tmscl*((water(i)*RPHF*TEMP2(I,L))**ZP578)
2316 AREARAT = ONE-(EXP(EXPARG))
2317 watevap(I) = EVP9(I,L)*AREARAT*fracwat
2318 IF(watevap(I).GE.water(i)) watevap(I) = water(i)
2319 EVP9(I,L)=EVP9(I,L)-watevap(I)
2320 water(i) = water(i) - watevap(i)
2321 endif
2322
2323 QL(I,L) = QL(I,L)+(iceevap(i)+watevap(i))*TEMP3(I,L)
2324 TL(I,L) = TL(I,L)-(iceevap(i)+watevap(i))*TEMP3(I,L)*ELOCP
2325
2326 20 CONTINUE
2327
2328 100 CONTINUE
2329
2330 do i = 1,irun
2331 rcon(i) = water(i) + crystal(i)
2332 enddo
2333
2334 C**********************************************************************
2335 C Large Scale Precip
2336 C**********************************************************************
2337
2338 DO 200 L=1,NLAY
2339 DO I = 1,IRUN
2340 rainmax(i) = rhcrit(i,L)*evp9(i,L) +
2341 . ql(i,L)*(rhcrit(i,L)-1.)/temp3(i,L)
2342
2343 if (rainmax(i).LE.0.0) then
2344 call qsat( tl(i,L),pl(i,L),rh,dum,.false.)
2345 rh = ql(i,L)/rh
2346
2347 if( rhcrit(i,L).eq.1.0 ) then
2348 fact = 1.0
2349 else
2350 fact = min( 1.0, alpha + (1.0-alpha)*( rh-rhcrit(i,L)) /
2351 1 (1.0-rhcrit(i,L)) )
2352 endif
2353
2354 C Do not allow clouds above 10 mb
2355 if( pl(i,L).ge.10.0 ) CLSBTH(I,L) = fact
2356 RLAR(I) = RLAR(I)-rainmax(I)
2357 QL(I,L) = QL(I,L)+rainmax(I)*TEMP3(I,L)
2358 TL(I,L) = TL(I,L)-rainmax(I)*TEMP3(I,L)*ELOCP
2359 C Large-scale water
2360 cldlz(i,L) = cldlz(i,L) - rainmax(i)*temp3(i,L)
2361 ENDIF
2362 ENDDO
2363
2364 DO I=1,IRUN
2365 IF((RLAR(I).GT.0.0).AND.(rainmax(I).GT.0.0))THEN
2366 RPOW=(RLAR(I)*RPHF*TEMP2(I,L))**ZP578
2367 EXPARG = ZM1P04*tmscl*RPOW
2368 AREARAT = ONE-(EXP(EXPARG))
2369 TEMP1(I,7) = rainmax(I)*AREARAT
2370 IF(TEMP1(I,7).GE.RLAR(I)) TEMP1(I,7) = RLAR(I)
2371 RLAR(I) = RLAR(I)-TEMP1(I,7)
2372 QL(I,L) = QL(I,L)+TEMP1(I,7)*TEMP3(I,L)
2373 TL(I,L) = TL(I,L)-TEMP1(I,7)*TEMP3(I,L)*ELOCP
2374 ENDIF
2375 ENDDO
2376
2377 200 CONTINUE
2378
2379 RETURN
2380 END
2381 subroutine srclouds (th,q,plk,pl,plke,cloud,cldwat,irun,irise,
2382 1 rhc,offset,alpha)
2383 C***********************************************************************
2384 C
2385 C PURPOSE:
2386 C ========
2387 C Compute non-precipitating cloud fractions
2388 C based on Slingo and Ritter (1985).
2389 C Remove cloudiness where conditionally unstable.
2390 C
2391 C INPUT:
2392 C ======
2393 C th ......... Potential Temperature (irun,irise)
2394 C q .......... Specific Humidity (irun,irise)
2395 C plk ........ P**Kappa at mid-layer (irun,irise)
2396 C pl ......... Pressure at mid-layer (irun,irise)
2397 C plke ....... P**Kappa at edge (irun,irise+1)
2398 C irun ....... Horizontal dimension
2399 C irise ...... Vertical dimension
2400 C
2401 C OUTPUT:
2402 C =======
2403 C cloud ...... Cloud Fraction (irun,irise)
2404 C
2405 C***********************************************************************
2406 C* GODDARD LABORATORY FOR ATMOSPHERES *
2407 C***********************************************************************
2408
2409 implicit none
2410 integer irun,irise
2411
2412 real th(irun,irise)
2413 real q(irun,irise)
2414 real plk(irun,irise)
2415 real pl(irun,irise)
2416 real plke(irun,irise+1)
2417
2418 real tempth(irun)
2419 real tempqs(irun)
2420 real dhstar(irun)
2421 real cloud(irun,irise)
2422 real cldwat(irun,irise)
2423 real qs(irun,irise)
2424
2425 real cp, alhl, getcon, akap, pcheck
2426 real ratio, temp, pke, elocp
2427 real rhcrit,rh,dum,pbar,tbar
2428 integer i,L,ntradesu,ntradesl
2429
2430 real factor
2431 real rhc(irun,irise)
2432 real offset,alpha
2433
2434 c Explicit Inline Directives
2435 c --------------------------
2436 #if CRAY
2437 #if f77
2438 cfpp$ expand (qsat)
2439 #endif
2440 #endif
2441
2442 cp = getcon('CP')
2443 alhl = getcon('LATENT HEAT COND')
2444 elocp = alhl/cp
2445 akap = getcon('KAPPA')
2446
2447 do L = 1,irise
2448 do i = 1,irun
2449 temp = th(i,L)*plk(i,L)
2450 call qsat ( temp,pl(i,L),qs(i,L),dum,.false. )
2451 enddo
2452 enddo
2453
2454 do L = 2,irise
2455 do i = 1,irun
2456 rh = q(i,L)/qs(i,L)
2457
2458 rhcrit = rhc(i,L) - offset
2459 ratio = alpha*(rh-rhcrit)/offset
2460
2461 if(cloud(i,L).eq. 0.0 .and. ratio.gt.0.0 ) then
2462 cloud(i,L) = min( ratio,1.0 )
2463 endif
2464
2465 enddo
2466 enddo
2467
2468 c Reduce clouds from conditionally unstable layer
2469 c -----------------------------------------------
2470 call ctei ( th,q,cloud,cldwat,pl,plk,plke,irun,irise )
2471
2472 return
2473 end
2474
2475 subroutine ctei ( th,q,cldfrc,cldwat,pl,plk,plke,im,lm )
2476 implicit none
2477 integer im,lm
2478 real th(im,lm),q(im,lm),plke(im,lm+1),cldwat(im,lm)
2479 real plk(im,lm),pl(im,lm),cldfrc(im,lm)
2480 integer i,L
2481 real getcon,cp,alhl,elocp,cpoel,t,p,s,qs,dqsdt,dq
2482 real k,krd,kmm,f
2483
2484 cp = getcon('CP')
2485 alhl = getcon('LATENT HEAT COND')
2486 cpoel = cp/alhl
2487 elocp = alhl/cp
2488
2489 do L=lm,2,-1
2490 do i=1,im
2491 dq = q(i,L)+cldwat(i,L)-q(i,L-1)-cldwat(i,L-1)
2492 if( dq.eq.0.0 ) dq = 1.0e-20
2493 k = 1.0 + cpoel*plke(i,L)*( th(i,L)-th(i,L-1) ) / dq
2494
2495 t = th(i,L)*plk(i,L)
2496 p = pl(i,L)
2497 call qsat ( t,p,qs,dqsdt,.true. )
2498
2499 krd = ( cpoel*t*(1+elocp*dqsdt) )/( 1 + 1.608*dqsdt*t )
2500
2501 kmm = ( 1+elocp*dqsdt )*( 1 + 0.392*cpoel*t )
2502 . / ( 2+(1+1.608*cpoel*t)*elocp*dqsdt )
2503
2504 s = ( (k-krd)/(kmm-krd) )
2505 f = 1.0 - min( 1.0, max(0.0,1.0-exp(-s)) )
2506
2507 cldfrc(i,L) = cldfrc(i,L)*f
2508 cldwat(i,L) = cldwat(i,L)*f
2509
2510 enddo
2511 enddo
2512
2513 return
2514 end
2515
2516 subroutine back2grd(gathered,indeces,scattered,irun)
2517 implicit none
2518 integer i,irun,indeces(irun)
2519 real gathered(irun),scattered(irun)
2520 real temp(irun)
2521 do i = 1,irun
2522 temp(indeces(i)) = gathered(i)
2523 enddo
2524 do i = 1,irun
2525 scattered(i) = temp(i)
2526 enddo
2527 return
2528 end

  ViewVC Help
Powered by ViewVC 1.1.22