/[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.5 - (show annotations) (download)
Wed Jul 7 20:07:47 2004 UTC (19 years, 11 months ago) by molod
Branch: MAIN
Changes since 1.4: +23 -41 lines
Still removing sigma references....

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

  ViewVC Help
Powered by ViewVC 1.1.22