/[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.17 - (show annotations) (download)
Thu Aug 12 15:21:22 2004 UTC (19 years, 10 months ago) by molod
Branch: MAIN
Changes since 1.16: +1 -3 lines
Debugging - Code for work-around for input vegetation dataset

1 C $Header: /u/gcmpack/MITgcm/pkg/fizhi/fizhi_moist.F,v 1.16 2004/08/10 15:13:47 molod Exp $
2 C $Name: $
3
4 #include "FIZHI_OPTIONS.h"
5 subroutine moistio (ndmoist,istrip,npcs,
6 . lowlevel,midlevel,nltop,nsubmin,nsubmax,Lup,
7 . pz,plz,plze,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 implicit none
16
17 #ifdef ALLOW_DIAGNOSTICS
18 #include "SIZE.h"
19 #include "diagnostics_SIZE.h"
20 #include "diagnostics.h"
21 #endif
22
23 c Input Variables
24 c ---------------
25 integer im,jm,lm
26 integer ndmoist,istrip,npcs
27 integer bi,bj,ntracer,ptracer
28 integer lowlevel,midlevel,nltop,nsubmin,nsubmax,Lup
29 _RL pz(im,jm),plz(im,jm,lm),plze(im,jm,lm+1),dpres(im,jm,lm)
30 _RL pkht(im,jm,lm+1),pkl(im,jm,lm)
31 _RL tz(im,jm,lm),qz(im,jm,lm,ntracer)
32 _RL qqz(im,jm,lm)
33 _RL dumoist(im,jm,lm),dvmoist(im,jm,lm)
34 _RL dtmoist(im,jm,lm),dqmoist(im,jm,lm,ntracer)
35 _RL ptop
36 integer iras
37 _RL rainlsp(im,jm),rainconv(im,jm),snowfall(im,jm)
38 integer nswcld,nswlz
39 _RL cldlsp_sw(im,jm,lm),cldras_sw(im,jm,lm)
40 _RL cldtot_sw(im,jm,lm),swlz(im,jm,lm)
41 integer nlwcld,nlwlz
42 _RL cldlsp_lw(im,jm,lm),cldras_lw(im,jm,lm)
43 _RL cldtot_lw(im,jm,lm),lwlz(im,jm,lm)
44 logical lpnt
45 integer myid
46
47 c Local Variables
48 c ---------------
49 integer ncrnd,nsecf
50
51 _RL fracqq, dum
52 integer snowcrit
53 parameter (fracqq = 0.1)
54 _RL one
55 parameter (one=1.)
56
57 _RL cldsr(im,jm,lm)
58 _RL srcld(istrip,lm)
59
60 _RL plev
61 _RL cldnow,cldlsp_mem,cldlsp,cldras_mem,cldras
62 _RL watnow,watmin,cldmin
63 _RL cldprs(im,jm),cldtmp(im,jm)
64 _RL cldhi (im,jm),cldlow(im,jm)
65 _RL cldmid(im,jm),totcld(im,jm)
66
67 _RL CLDLS(im,jm,lm) , CPEN(im,jm,lm)
68 _RL tmpimjm(im,jm)
69 _RL lsp_new(im,jm)
70 _RL conv_new(im,jm)
71 _RL snow_new(im,jm)
72
73 _RL qqcolmin(im,jm)
74 _RL qqcolmax(im,jm)
75 integer levpbl(im,jm)
76
77 c Gathered Arrays for Variable Cloud Base
78 c ---------------------------------------
79 _RL raincgath(im*jm)
80 _RL pigather(im*jm)
81 _RL thgather(im*jm,lm)
82 _RL shgather(im*jm,lm)
83 _RL pkzgather(im*jm,lm)
84 _RL pkegather(im*jm,lm+1)
85 _RL plzgather(im*jm,lm)
86 _RL plegather(im*jm,lm+1)
87 _RL dpgather(im*jm,lm)
88 _RL tmpgather(im*jm,lm)
89 _RL deltgather(im*jm,lm)
90 _RL delqgather(im*jm,lm)
91 _RL ugather(im*jm,lm,ntracer)
92 _RL delugather(im*jm,lm,ntracer)
93 _RL deltrnev(im*jm,lm)
94 _RL delqrnev(im*jm,lm)
95
96 integer nindeces(lm)
97 integer pblindex(im*jm)
98 integer levgather(im*jm)
99
100 c Stripped Arrays
101 c ---------------
102 _RL saveth (istrip,lm)
103 _RL saveq (istrip,lm)
104 _RL saveu (istrip,lm,ntracer)
105 _RL usubcl (istrip, ntracer)
106
107 _RL ple(istrip,lm+1)
108 _RL dp(istrip,lm)
109 _RL TL(ISTRIP,lm) , SHL(ISTRIP,lm)
110 _RL PL(ISTRIP,lm) , PLK(ISTRIP,lm)
111 _RL PLKE(ISTRIP,lm+1)
112 _RL TH(ISTRIP,lm) ,CVTH(ISTRIP,lm)
113 _RL CVQ(ISTRIP,lm)
114 _RL UL(ISTRIP,lm,ntracer)
115 _RL cvu(istrip,lm,ntracer)
116 _RL CLMAXO(ISTRIP,lm),CLBOTH(ISTRIP,lm)
117 _RL CLSBTH(ISTRIP,lm)
118 _RL TMP1(ISTRIP,lm), TMP2(ISTRIP,lm)
119 _RL TMP3(ISTRIP,lm), TMP4(ISTRIP,lm+1)
120 _RL TMP5(ISTRIP,lm+1)
121 integer ITMP1(ISTRIP,lm), ITMP2(ISTRIP,lm)
122
123 _RL PRECIP(ISTRIP), PCNET(ISTRIP)
124 _RL SP(ISTRIP), PREP(ISTRIP)
125 _RL PCPEN (ISTRIP,lm)
126 integer pbl(istrip),depths(lm)
127
128 _RL cldlz(istrip,lm), cldwater(im,jm,lm)
129 _RL rhfrac(istrip), rhmin, pup, ppbl, rhcrit(istrip,lm)
130 _RL offset, alpha, rasmax
131
132 logical first
133 logical lras
134 _RL clfrac (istrip,lm)
135 _RL cldmas (istrip,lm)
136 _RL detrain(istrip,lm)
137 _RL psubcld (istrip), psubcldg (im,jm)
138 _RL psubcld_cnt(istrip), psubcldgc(im,jm)
139 _RL rnd(lm/2)
140 DATA FIRST /.TRUE./
141
142 integer imstp,nsubcl,nlras
143 integer i,j,iloop,indx,indgath,l,nn,num,numdeps,nt
144 _RL tmstp,tminv,sday,grav,alhl,cp,elocp,gamfac
145 _RL rkappa,p0kappa,p0kinv,ptopkap,pcheck
146 _RL 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.1) 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 indx = 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 indx = indx + 1
239 pblindex(indx) = (j-1)*im + i
240 endif
241 enddo
242 enddo
243 enddo
244
245 do indx = 1,im*jm
246 levgather(indx) = levpbl(pblindex(indx),1)
247 pigather(indx) = pz(pblindex(indx),1)
248 pkegather(indx,lm+1) = pkht(pblindex(indx),1,lm+1)
249 plegather(indx,lm+1) = plze(pblindex(indx),1,lm+1)
250 enddo
251
252 do L = 1,lm
253 do indx = 1,im*jm
254 thgather(indx,L) = tz(pblindex(indx),1,L)
255 shgather(indx,L) = qz(pblindex(indx),1,L,1)
256 pkegather(indx,L) = pkht(pblindex(indx),1,L)
257 pkzgather(indx,L) = pkl(pblindex(indx),1,L)
258 plegather(indx,L) = plze(pblindex(indx),1,L)
259 plzgather(indx,L) = plz(pblindex(indx),1,L)
260 dpgather(indx,L) = dpres(pblindex(indx),1,L)
261 enddo
262 enddo
263 do nt = 1,ntracer-ptracer
264 do L = 1,lm
265 do indx = 1,im*jm
266 ugather(indx,L,nt) = qz(pblindex(indx),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,one,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
780 c Sub-Cloud Layer
781 c -------------------------
782 if( ipsubcld.ne.0 ) then
783 do j = 1,jm
784 do i = 1,im
785 qdiag(i,j,ipsubcld,bi,bj) = qdiag(i,j,ipsubcld,bi,bj) +
786 . psubcldg (i,j)
787 qdiag(i,j,ipsubcldc,bi,bj) = qdiag(i,j,ipsubcldc,bi,bj) +
788 . psubcldgc(i,j)
789 enddo
790 enddo
791 endif
792
793 c Non-Precipitating Cloud Fraction
794 c --------------------------------
795 if( icldnp.ne.0 ) then
796 do L = 1,lm
797 do j = 1,jm
798 do i = 1,im
799 qdiag(i,j,icldnp+L-1,bi,bj) = qdiag(i,j,icldnp+L-1,bi,bj) +
800 . cldsr(i,j,L)
801 enddo
802 enddo
803 enddo
804 ncldnp = ncldnp + 1
805 endif
806
807 c Moist Processes Heating Rate
808 c ----------------------------
809 if(imoistt.gt.0) then
810 do L = 1,lm
811 do j = 1,jm
812 do i = 1,im
813 indgath = (j-1)*im + i
814 qdiag(i,j,imoistt+L-1,bi,bj) = qdiag(i,j,imoistt+L-1,bi,bj) +
815 . (dtmoist(i,j,L)*sday*pkzgather(indgath,L)/pz(i,j))
816 enddo
817 enddo
818 enddo
819 endif
820
821 c Moist Processes Moistening Rate
822 c -------------------------------
823 if(imoistq.gt.0) then
824 do L = 1,lm
825 do j = 1,jm
826 do i = 1,im
827 qdiag(i,j,imoistq+L-1,bi,bj) = qdiag(i,j,imoistq+L-1,bi,bj) +
828 . (dqmoist(i,j,L,1)*sday*1000.0/pz(i,j))
829 enddo
830 enddo
831 enddo
832 endif
833
834 c Cloud Mass Flux
835 c ---------------
836 if(icldmas.gt.0) then
837 do L = 1,lm
838 do j = 1,jm
839 do i = 1,im
840 indgath = (j-1)*im + i
841 qdiag(i,j,icldmas+L-1,bi,bj) = qdiag(i,j,icldmas+L-1,bi,bj) +
842 . tmpgather(indgath,L)
843 enddo
844 enddo
845 enddo
846 endif
847
848 c Detrained Cloud Mass Flux
849 c -------------------------
850 if(idtrain.gt.0) then
851 do L = 1,lm
852 do j = 1,jm
853 do i = 1,im
854 indgath = (j-1)*im + i
855 qdiag(i,j,idtrain+L-1,bi,bj) = qdiag(i,j,idtrain+L-1,bi,bj) +
856 . pkegather(indgath,L)
857 enddo
858 enddo
859 enddo
860 endif
861
862 c Grid-Scale Condensational Heating Rate
863 c --------------------------------------
864 if(idtls.gt.0) then
865 do L = 1,lm
866 do j = 1,jm
867 do i = 1,im
868 indgath = (j-1)*im + i
869 qdiag(i,j,idtls+L-1,bi,bj) = qdiag(i,j,idtls+L-1,bi,bj) +
870 . deltrnev(indgath,L)
871 enddo
872 enddo
873 enddo
874 endif
875
876 c Grid-Scale Condensational Moistening Rate
877 c -----------------------------------------
878 if(idqls.gt.0) then
879 do L = 1,lm
880 do j = 1,jm
881 do i = 1,im
882 indgath = (j-1)*im + i
883 qdiag(i,j,idqls+L-1,bi,bj) = qdiag(i,j,idqls+L-1,bi,bj) +
884 . delqrnev(indgath,L)
885 enddo
886 enddo
887 enddo
888 endif
889
890 c Total Precipitation
891 c -------------------
892 if(ipreacc.gt.0) then
893 do j = 1,jm
894 do i = 1,im
895 qdiag(i,j,ipreacc,bi,bj) = qdiag(i,j,ipreacc,bi,bj)
896 . + ( lsp_new(I,j)
897 . + snow_new(I,j)
898 . + conv_new(i,j) ) *sday*tminv
899 enddo
900 enddo
901 endif
902
903 c Convective Precipitation
904 c ------------------------
905 if(iprecon.gt.0) then
906 do j = 1,jm
907 do i = 1,im
908 indgath = (j-1)*im + i
909 qdiag(i,j,iprecon,bi,bj) = qdiag(i,j,iprecon,bi,bj) +
910 . raincgath(indgath)*sday*tminv
911 enddo
912 enddo
913 endif
914
915 C **********************************************************************
916 C **** Fill Rainfall and Snowfall Arrays for Land Surface Model ****
917 C **** Note: Precip Rates work when DT(turb)<DT(moist) ****
918 C **********************************************************************
919
920 do j = 1,jm
921 do i = 1,im
922 rainlsp (i,j) = rainlsp (i,j) + lsp_new(i,j)*tminv
923 rainconv(i,j) = rainconv(i,j) + conv_new(i,j)*tminv
924 snowfall(i,j) = snowfall(i,j) + snow_new(i,j)*tminv
925 enddo
926 enddo
927
928 C **********************************************************************
929 C *** Compute Time-averaged Quantities for Radiation ***
930 C *** CPEN => Cloud Fraction from RAS ***
931 C *** CLDLS => Cloud Fraction from RNEVP ***
932 C **********************************************************************
933
934 do j = 1,jm
935 do i = 1,im
936 cldhi (i,j) = 0.
937 cldmid(i,j) = 0.
938 cldlow(i,j) = 0.
939 cldtmp(i,j) = 0.
940 cldprs(i,j) = 0.
941 tmpimjm(i,j) = 0.
942 enddo
943 enddo
944
945 c Set Moist-Process Memory Coefficient
946 c ------------------------------------
947 cldras_mem = 1.0-tmstp/ 3600.0
948 cldlsp_mem = 1.0-tmstp/(3600.0*3)
949
950 do L = 1,lm
951 do i = 1,im*jm
952 plev = pl(i,L)
953
954 c Compute Time-averaged Cloud and Water Amounts for Longwave Radiation
955 c --------------------------------------------------------------------
956 watnow = cldwater(i,1,L)
957 if( plev.le.500.0 ) then
958 cldras = min( max( cldras_lw(i,1,L)*cldras_mem,cpen(i,1,L)),1.0)
959 else
960 cldras = 0.0
961 endif
962 cldlsp = min( max( cldlsp_lw(i,1,L)*cldlsp_mem,cldls(i,1,L)),1.0)
963
964 if( cldras.lt.cldmin ) cldras = 0.0
965 if( cldlsp.lt.cldmin ) cldlsp = 0.0
966
967 cldnow = max( cldlsp,cldras )
968
969 lwlz(i,1,L) = ( nlwlz*lwlz(i,1,L) + watnow)/(nlwlz +1)
970 cldtot_lw(i,1,L) = (nlwcld*cldtot_lw(i,1,L) + cldnow)/(nlwcld+1)
971 cldlsp_lw(i,1,L) = (nlwcld*cldlsp_lw(i,1,L) + cldlsp)/(nlwcld+1)
972 cldras_lw(i,1,L) = (nlwcld*cldras_lw(i,1,L) + cldras)/(nlwcld+1)
973
974
975 c Compute Time-averaged Cloud and Water Amounts for Shortwave Radiation
976 c ---------------------------------------------------------------------
977 watnow = cldwater(i,1,L)
978 if( plev.le.500.0 ) then
979 cldras = min( max(cldras_sw(i,1,L)*cldras_mem, cpen(i,1,L)),1.0)
980 else
981 cldras = 0.0
982 endif
983 cldlsp = min( max(cldlsp_sw(i,1,L)*cldlsp_mem,cldls(i,1,L)),1.0)
984
985 if( cldras.lt.cldmin ) cldras = 0.0
986 if( cldlsp.lt.cldmin ) cldlsp = 0.0
987
988 cldnow = max( cldlsp,cldras )
989
990 swlz(i,1,L) = ( nswlz*swlz(i,1,L) + watnow)/(nswlz +1)
991 cldtot_sw(i,1,L) = (nswcld*cldtot_sw(i,1,L) + cldnow)/(nswcld+1)
992 cldlsp_sw(i,1,L) = (nswcld*cldlsp_sw(i,1,L) + cldlsp)/(nswcld+1)
993 cldras_sw(i,1,L) = (nswcld*cldras_sw(i,1,L) + cldras)/(nswcld+1)
994
995
996 c Compute Instantaneous Low-Mid-High Maximum Overlap Cloud Fractions
997 c ----------------------------------------------------------------------
998
999 if( L.lt.midlevel ) cldhi (i,1) = max( cldnow,cldhi (i,1) )
1000 if( L.ge.midlevel .and.
1001 . L.lt.lowlevel ) cldmid(i,1) = max( cldnow,cldmid(i,1) )
1002 if( L.ge.lowlevel ) cldlow(i,1) = max( cldnow,cldlow(i,1) )
1003
1004 c Compute Cloud-Top Temperature and Pressure
1005 c ------------------------------------------
1006 cldtmp(i,1) = cldtmp(i,1) + cldnow*pkzgather(i,L)
1007 . * ( tz(i,1,L) + dtmoist(i,1,L)*tmstp/pz(i,1) )
1008 cldprs(i,1) = cldprs(i,1) + cldnow*plev
1009 tmpimjm(i,1) = tmpimjm(i,1) + cldnow
1010
1011 enddo
1012 enddo
1013
1014 c Compute Instantanious Total 2-D Cloud Fraction
1015 c ----------------------------------------------
1016 do j = 1,jm
1017 do i = 1,im
1018 totcld(i,j) = 1.0 - (1.-cldhi (i,j))
1019 . * (1.-cldmid(i,j))
1020 . * (1.-cldlow(i,j))
1021 enddo
1022 enddo
1023
1024
1025 C **********************************************************************
1026 C *** Fill Cloud Top Pressure and Temperature Diagnostic ***
1027 C **********************************************************************
1028
1029 if(icldtmp.gt.0) then
1030 do j = 1,jm
1031 do i = 1,im
1032 if( cldtmp(i,j).gt.0.0 ) then
1033 qdiag(i,j,icldtmp,bi,bj) = qdiag(i,j,icldtmp,bi,bj) +
1034 . cldtmp(i,j)*totcld(i,j)/tmpimjm(i,j)
1035 qdiag(i,j,icttcnt,bi,bj) = qdiag(i,j,icttcnt,bi,bj) +
1036 . totcld(i,j)
1037 endif
1038 enddo
1039 enddo
1040 endif
1041
1042 if(icldprs.gt.0) then
1043 do j = 1,jm
1044 do i = 1,im
1045 if( cldprs(i,j).gt.0.0 ) then
1046 qdiag(i,j,icldprs,bi,bj) = qdiag(i,j,icldprs,bi,bj) +
1047 . cldprs(i,j)*totcld(i,j)/tmpimjm(i,j)
1048 qdiag(i,j,ictpcnt,bi,bj) = qdiag(i,j,ictpcnt,bi,bj) +
1049 . totcld(i,j)
1050 endif
1051 enddo
1052 enddo
1053 endif
1054
1055 C **********************************************************************
1056 C **** INCREMENT COUNTERS ****
1057 C **********************************************************************
1058
1059 nlwlz = nlwlz + 1
1060 nswlz = nswlz + 1
1061
1062 nlwcld = nlwcld + 1
1063 nswcld = nswcld + 1
1064
1065 #ifdef ALLOW_DIAGNOSTICS
1066 if( (bi.eq.1) .and. (bj.eq.1) ) then
1067 nmoistt = nmoistt + 1
1068 nmoistq = nmoistq + 1
1069 npreacc = npreacc + 1
1070 nprecon = nprecon + 1
1071
1072 ncldmas = ncldmas + 1
1073 ndtrain = ndtrain + 1
1074
1075 ndtls = ndtls + 1
1076 ndqls = ndqls + 1
1077 endif
1078 #endif
1079
1080 RETURN
1081 END
1082 SUBROUTINE RAS( NN, LNG, LENC, K, NLTOP, nlayr, DT
1083 *, UOI, ntracer, POI, QOI, PRS, PRJ, rnd, ncrnd
1084 *, RAINS, CLN, CLF, cldmas, detrain
1085 *, cp,grav,rkappa,alhl,rhfrac,rasmax )
1086 C
1087 C*********************************************************************
1088 C********************* SUBROUTINE RAS *****************************
1089 C********************** 16 MARCH 1988 ******************************
1090 C*********************************************************************
1091 C
1092 implicit none
1093
1094 C Argument List
1095 integer nn,lng,lenc,k,nltop,nlayr
1096 integer ntracer
1097 integer ncrnd
1098 _RL dt
1099 _RL UOI(lng,nlayr,ntracer), POI(lng,K)
1100 _RL QOI(lng,K), PRS(lng,K+1), PRJ(lng,K+1)
1101 _RL rnd(ncrnd)
1102 _RL RAINS(lng,K), CLN(lng,K), CLF(lng,K)
1103 _RL cldmas(lng,K), detrain(lng,K)
1104 _RL cp,grav,rkappa,alhl,rhfrac(lng),rasmax
1105
1106 C Local Variables
1107 _RL TCU(lng,K), QCU(lng,K)
1108 _RL ucu(lng,K,ntracer)
1109 _RL ALF(lng,K), BET(lng,K), GAM(lng,K)
1110 *, ETA(lng,K), HOI(lng,K)
1111 *, PRH(lng,K), PRI(lng,K)
1112 _RL HST(lng,K), QOL(lng,K), GMH(lng,K)
1113
1114 _RL TX1(lng), TX2(lng), TX3(lng), TX4(lng), TX5(lng)
1115 *, TX6(lng), TX7(lng), TX8(lng), TX9(lng)
1116 *, TX11(lng), TX12(lng), TX13(lng), TX14(lng,ntracer)
1117 *, TX15(lng)
1118 *, WFN(lng)
1119 integer IA1(lng), IA2(lng), IA3(lng)
1120 _RL cloudn(lng), pcu(lng)
1121
1122 integer krmin,icm
1123 _RL rknob, cmb2pa
1124 PARAMETER (KRMIN=01)
1125 PARAMETER (ICM=1000)
1126 PARAMETER (CMB2PA=100.0)
1127 PARAMETER (rknob = 10.)
1128
1129 integer IC(ICM), IRND(icm)
1130 _RL cmass(lng,K)
1131 LOGICAL SETRAS
1132
1133 integer i,L,nc,ib,nt
1134 integer km1,kp1,kprv,kcr,kfx,ncmx
1135 _RL p00, crtmsf, frac, rasblf
1136
1137 do L = 1,k
1138 do I = 1,LENC
1139 rains(i,l) = 0.
1140 enddo
1141 enddo
1142
1143 p00 = 1000.
1144 crtmsf = 0.
1145
1146 C The numerator here is the fraction of the subcloud layer mass flux
1147 C allowed to entrain into the cloud
1148
1149 CCC FRAC = 1./dt
1150 FRAC = 0.5/dt
1151
1152 KM1 = K - 1
1153 KP1 = K + 1
1154 C we want the ras adjustment time scale to be one hour (indep of dt)
1155 RASBLF = 1./3600.
1156 C
1157 KPRV = KM1
1158 C Removed KRMAX parameter
1159 KCR = MIN(KM1,nlayr-2)
1160 KFX = KM1 - KCR
1161 NCMX = KFX + NCRND
1162 C
1163 IF (KFX .GT. 0) THEN
1164 DO NC=1,KFX
1165 IC(NC) = K - NC
1166 ENDDO
1167 ENDIF
1168 C
1169 IF (NCRND .GT. 0) THEN
1170 DO I=1,ncrnd
1171 IRND(I) = (RND(I)-0.0005)*(KCR-KRMIN+1)
1172 IRND(I) = IRND(I) + KRMIN
1173 ENDDO
1174 C
1175 DO NC=1,NCRND
1176 IC(KFX+NC) = IRND(NC)
1177 ENDDO
1178 ENDIF
1179 C
1180 DO 100 NC=1,NCMX
1181 C
1182 IF (NC .EQ. 1 ) THEN
1183 SETRAS = .TRUE.
1184 ELSE
1185 SETRAS = .FALSE.
1186 ENDIF
1187 IB = IC(NC)
1188
1189 c Initialize Cloud Fraction Array
1190 c -------------------------------
1191 do i = 1,lenc
1192 cloudn(i) = 0.0
1193 enddo
1194
1195 CALL CLOUD(nn,lng, LENC, K, NLTOP, nlayr, IB, RASBLF,SETRAS,FRAC
1196 *, CP, ALHL, RKAPPA, GRAV, P00, CRTMSF
1197 *, POI, QOI, UOI, Ntracer, PRS, PRJ
1198 *, PCU, CLOUDN, TCU, QCU, UCU, CMASS
1199 *, ALF, BET, GAM, PRH, PRI, HOI, ETA
1200 *, HST, QOL, GMH
1201 *, TX1, TX2, TX3, TX4, TX5, TX6, TX7, TX8, TX9
1202 *, WFN, TX11, TX12, TX13, TX14, TX15
1203 *, IA1,IA2,IA3,rhfrac)
1204
1205 C Compute fraction of grid box into which rain re-evap occurs (clf)
1206 c -----------------------------------------------------------------
1207 do i = 1,lenc
1208
1209 c mass in detrainment layer
1210 c -------------------------
1211 tx1(i) = cmb2pa * (prs(i,ib+1) - prs(i,ib))/(grav*dt)
1212
1213 c ratio of detraining cloud mass to mass in detrainment layer
1214 c -----------------------------------------------------------
1215 tx1(i) = rhfrac(i)*rknob * cmass(i,ib) / tx1(i)
1216 if(cmass(i,K).gt.0.) clf(i,ib) = clf(i,ib) + tx1(i)
1217 if( clf(i,ib).gt.1.) clf(i,ib) = 1.
1218 enddo
1219
1220 c Compute Total Cloud Mass Flux
1221 c *****************************
1222 do L=ib,k
1223 do i=1,lenc
1224 cmass(i,L) = rhfrac(i)*cmass(i,L) * dt
1225 enddo
1226 enddo
1227
1228 do L=ib,k
1229 do i=1,lenc
1230 cldmas(i,L) = cldmas(i,L) + cmass(i,L)
1231 enddo
1232 enddo
1233
1234 do i=1,lenc
1235 detrain(i,ib) = detrain(i,ib) + cmass(i,ib)
1236 enddo
1237
1238 DO L=IB,K
1239 DO I=1,LENC
1240 POI(I,L) = POI(I,L) + TCU(I,L) * DT * rhfrac(i)
1241 QOI(I,L) = QOI(I,L) + QCU(I,L) * DT * rhfrac(i)
1242 ENDDO
1243 ENDDO
1244 DO NT=1,Ntracer
1245 DO L=IB,K
1246 DO I=1,LENC
1247 UOI(I,L+nltop-1,NT)=UOI(I,L+nltop-1,NT)+UCU(I,L,NT)*DT*rhfrac(i)
1248 ENDDO
1249 ENDDO
1250 ENDDO
1251 DO I=1,LENC
1252 rains(I,ib) = rains(I,ib) + PCU(I)*dt * rhfrac(i)
1253 ENDDO
1254
1255 100 CONTINUE
1256
1257 c Fill Convective Cloud Fractions based on 3-D Rain Amounts
1258 c ---------------------------------------------------------
1259 do L=k-1,1,-1
1260 do i=1,lenc
1261 tx1(i) = 100*(prs(i,L+1)-prs(i,L))/grav
1262 cln(i,L) = min(1600*rains(i,L)/tx1(i),rasmax )
1263 enddo
1264 enddo
1265
1266 RETURN
1267 END
1268 subroutine rndcloud (iras,nrnd,rnd,myid)
1269 implicit none
1270 integer n,iras,nrnd,myid
1271 _RL random_numbx
1272 _RL rnd(nrnd)
1273 integer irm
1274 parameter (irm = 1000)
1275 _RL random(irm)
1276 integer i,mcheck,numrand,iseed,indx
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.1) 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.1 ) print *, 'Recreating Rand Numb Array in RNDCLOUD'
1297 if( myid.eq.1 ) 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(iseed)
1303 enddo
1304 indx = (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(iseed)
1313 enddo
1314 indx = iseed
1315
1316 c Multiple Time In But have NOT Used Up all 1000 numbers (MCHECK.NE.0)
1317 c --------------------------------------------------------------------
1318 else
1319 indx = (iras-1)*nrnd
1320 endif
1321
1322 indx = mod(indx,irm)
1323 if( indx+nrnd.gt.1000 ) indx=1000-nrnd
1324
1325 do n = 1,nrnd
1326 rnd(n) = random(indx+n)
1327 enddo
1328
1329 100 continue
1330 first = .false.
1331 iras0 = iras
1332 return
1333 end
1334 function random_numbx(iseed)
1335 implicit none
1336 integer iseed
1337 real *8 seed,port_rand
1338 _RL random_numbx
1339 random_numbx = 0
1340 #ifdef CRAY
1341 _RL ranf
1342 random_numbx = ranf()
1343 #else
1344 #ifdef SGI
1345 _RL rand
1346 random_numbx = rand()
1347 #endif
1348 random_numbx = port_rand(seed)
1349 #endif
1350 return
1351 end
1352 subroutine random_seedx (iseed)
1353 implicit none
1354 integer iseed
1355 #ifdef CRAY
1356 call ranset (iseed)
1357 #endif
1358 #ifdef SGI
1359 integer*4 seed
1360 seed = iseed
1361 call srand (seed)
1362 #endif
1363 return
1364 end
1365 SUBROUTINE CLOUD(nn,lng, LENC, K, NLTOP, nlayr, IC, RASALF
1366 *, SETRAS, FRAC
1367 *, CP, ALHL, RKAP, GRAV, P00, CRTMSF
1368 *, POI, QOI, UOI, Ntracer, PRS, PRJ
1369 *, PCU, CLN, TCU, QCU, UCU, CMASS
1370 *, ALF, BET, GAM, PRH, PRI, HOL, ETA
1371 *, HST, QOL, GMH
1372 *, TX1, TX2, TX3, TX4, TX5, TX6, TX7, TX8, ALM
1373 *, WFN, AKM, QS1, CLF, UHT, WLQ
1374 *, IA, I1, I2,rhfrac)
1375 C
1376 C*********************************************************************
1377 C******************** Relaxed Arakawa-Schubert ***********************
1378 C********************* Plug Compatible Version **********************
1379 C************************ SUBROUTINE CLOUD ***************************
1380 C************************* 23 JULY 1992 ***************************
1381 C*********************************************************************
1382 C*********************************************************************
1383 C*********************************************************************
1384 C************************** Developed By *****************************
1385 C************************** *****************************
1386 C************************ Shrinivas Moorthi **************************
1387 C************************ and **************************
1388 C************************ Max J. Suarez *****************************
1389 C************************ *****************************
1390 C******************** Laboratory for Atmospheres *********************
1391 C****************** NASA/GSFC, Greenbelt, MD 20771 *******************
1392 C*********************************************************************
1393 C*********************************************************************
1394 C
1395 C The calculations of Moorthi and Suarez (1992, MWR) are
1396 C contained in the CLOUD routine.
1397 C It is probably advisable, at least initially, to treat CLOUD
1398 C as a black box that computes the single cloud adjustments. RAS,
1399 C on the other hand, can be tailored to each GCMs configuration
1400 C (ie, number and placement of levels, nature of boundary layer,
1401 C time step and frequency with which RAS is called).
1402 C
1403 C
1404 C Input:
1405 C ------
1406 C
1407 C lng : The inner dimension of update and input arrays.
1408 C
1409 C LENC : The run: the number of soundings processes in a single call.
1410 C RAS works on the first LENC of the lng soundings
1411 C passed. This allows working on pieces of the world
1412 C say for multitasking, without declaring temporary arrays
1413 C and copying the data to and from them. This is an f77
1414 C version. An F90 version would have to allow more
1415 C flexibility in the argument declarations. Obviously
1416 C (LENC<=lng).
1417 C
1418 C K : Number of vertical layers (increasing downwards).
1419 C Need not be the same as the number of layers in the
1420 C GCM, since it is the outer dimension. The bottom layer
1421 C (K) is the subcloud layer.
1422 C
1423 C IC : Detrainment level to check for presence of convection
1424 C
1425 C RASALF : Relaxation parameter (< 1.) for present cloud-type
1426 C
1427 C SETRAS : Logical parameter to control re-calculation of
1428 C saturation specific humidity and mid level P**kappa
1429 C
1430 C FRAC : Fraction of the PBL (layer K) mass allowed to be used
1431 C by a cloud-type in time DT
1432 C
1433 C CP : Specific heat at constant pressure
1434 C
1435 C ALHL : Latent Heat of condensation
1436 C
1437 C RKAP : R/Cp, where R is the gas constant
1438 C
1439 C GRAV : Acceleration due to gravity
1440 C
1441 C P00 : A reference pressure in hPa, useually 1000 hPa
1442 C
1443 C CRTMSF : Critical value of mass flux above which cloudiness at
1444 C the detrainment layer of that cloud-type is assumed.
1445 C Affects only cloudiness calculation.
1446 C
1447 C POI : 2D array of dimension (lng,K) containing potential
1448 C temperature. Updated but not initialized by RAS.
1449 C
1450 C QOI : 2D array of dimension (lng,K) containing specific
1451 C humidity. Updated but not initialized by RAS.
1452 C
1453 C UOI : 3D array of dimension (lng,K,NTRACER) containing tracers
1454 C Updated but not initialized by RAS.
1455 C
1456 C PRS : 2D array of dimension (lng,K+1) containing pressure
1457 C in hPa at the interfaces of K-layers from top of the
1458 C atmosphere to the bottom. Not modified.
1459 C
1460 C PRJ : 2D array of dimension (lng,K+1) containing (PRS/P00) **
1461 C RKAP. i.e. Exner function at layer edges. Not modified.
1462 C
1463 C rhfrac : 1D array of dimension (lng) containing a rel.hum. scaling
1464 C fraction. Not modified.
1465 C
1466 C Output:
1467 C -------
1468 C
1469 C PCU : 1D array of length lng containing accumulated
1470 C precipitation in mm/sec.
1471 C
1472 C CLN : 2D array of dimension (lng,K) containing cloudiness
1473 C Note: CLN is bumped but NOT initialized
1474 C
1475 C TCU : 2D array of dimension (lng,K) containing accumulated
1476 C convective heating (K/sec).
1477 C
1478 C QCU : 2D array of dimension (lng,K) containing accumulated
1479 C convective drying (kg/kg/sec).
1480 C
1481 C CMASS : 2D array of dimension (lng,K) containing the
1482 C cloud mass flux (kg/sec). Filled from cloud top
1483 C to base.
1484 C
1485 C Temporaries:
1486 C
1487 C ALF, BET, GAM, ETA, PRH, PRI, HOI, HST, QOL, GMH are temporary
1488 C 2D real arrays of dimension of at least (LENC,K) where LENC is
1489 C the horizontal dimension over which convection is invoked.
1490 C
1491 C
1492 C TX1, TX2, TX3, TX4, TX5, TX6, TX7, TX8, TX9, AKM, QS1, CLF, UHT
1493 C VHT, WLQ WFN are temporary real arrays of length at least LENC
1494 C
1495 C IA, I1, and I2 are temporary integer arrays of length LENC
1496 C
1497 C
1498 C************************************************************************
1499 implicit none
1500 C Argument List declarations
1501 integer nn,lng,LENC,K,NLTOP,nlayr,ic,ntracer
1502 _RL rasalf
1503 LOGICAL SETRAS
1504 _RL frac, cp, alhl, rkap, grav, p00, crtmsf
1505 _RL POI(lng,K),QOI(lng,K),PRS(lng,K+1),PRJ(lng,K+1)
1506 _RL uoi(lng,nlayr,ntracer)
1507 _RL PCU(LENC), CLN(lng)
1508 _RL TCU(lng,K), QCU(lng,K), ucu(lng,k,ntracer), CMASS(lng,K)
1509 _RL ALF(lng,K), BET(lng,K), GAM(lng,K), PRH(lng,K), PRI(lng,K)
1510 _RL HOL(LENC,K), ETA(LENC,K), HST(LENC,K), QOL(LENC,K)
1511 _RL GMH(LENC,K)
1512 _RL TX1(LENC), TX2(LENC), TX3(LENC), TX4(LENC)
1513 _RL TX5(LENC), TX6(LENC), TX7(LENC), TX8(LENC)
1514 _RL ALM(LENC), WFN(LENC), AKM(LENC), QS1(LENC)
1515 _RL WLQ(LENC), CLF(LENC)
1516 _RL uht(lng,ntracer)
1517 integer IA(LENC), I1(LENC),I2(LENC)
1518 _RL rhfrac(lng)
1519
1520 C Local Variables
1521 _RL daylen,half,one,zero,cmb2pa,rhmax
1522 PARAMETER (DAYLEN=86400.0, HALF=0.5, ONE=1.0, ZERO=0.0)
1523 PARAMETER (CMB2PA=100.0)
1524 PARAMETER (RHMAX=0.9999)
1525 _RL rkapp1,onebcp,albcp,onebg,cpbg,twobal
1526 C
1527 integer nt,km1,ic1,i,L,len1,len2,isav,len11,ii
1528 integer lena,lena1,lenb
1529 _RL tem,tem1
1530
1531 c Explicit Inline Directives
1532 c --------------------------
1533 #ifdef CRAY
1534 #ifdef f77
1535 cfpp$ expand (qsat)
1536 #endif
1537 #endif
1538
1539 RKAPP1 = 1.0 + RKAP
1540 ONEBCP = 1.0 / CP
1541 ALBCP = ALHL * ONEBCP
1542 ONEBG = 1.0 / GRAV
1543 CPBG = CP * ONEBG
1544 TWOBAL = 2.0 / ALHL
1545 C
1546 KM1 = K - 1
1547 IC1 = IC + 1
1548 C
1549 C SETTING ALF, BET, GAM, PRH, AND PRI : DONE ONLY WHEN SETRAS=.T.
1550 C
1551
1552 IF (SETRAS) THEN
1553
1554 DO 2050 L=1,K
1555 DO 2030 I=1,LENC
1556 PRH(I,L) = (PRJ(I,L+1)*PRS(I,L+1) - PRJ(I,L)*PRS(I,L))
1557 * / ((PRS(I,L+1)-PRS(I,L)) * RKAPP1)
1558 2030 CONTINUE
1559 2050 CONTINUE
1560
1561 DO 2070 L=1,K
1562 DO 2060 I=1,LENC
1563 TX5(I) = POI(I,L) * PRH(I,L)
1564 TX1(I) = (PRS(I,L) + PRS(I,L+1)) * 0.5
1565 TX3(I) = TX5(I)
1566 CALL QSAT(TX3(I), TX1(I), TX2(I), TX4(I), .TRUE.)
1567 ALF(I,L) = TX2(I) - TX4(I) * TX5(I)
1568 BET(I,L) = TX4(I) * PRH(I,L)
1569 GAM(I,L) = 1.0 / ((1.0 + TX4(I)*ALBCP) * PRH(I,L))
1570 PRI(I,L) = (CP/CMB2PA) / (PRS(I,L+1) - PRS(I,L))
1571 2060 CONTINUE
1572 2070 CONTINUE
1573
1574 ENDIF
1575 C
1576 C
1577 DO 10 L=1,K
1578 DO 10 I=1,lng
1579 TCU(I,L) = 0.0
1580 QCU(I,L) = 0.0
1581 CMASS(I,L) = 0.0
1582 10 CONTINUE
1583
1584 do nt = 1,ntracer
1585 do L=1,K
1586 do I=1,LENC
1587 ucu(I,L,nt) = 0.0
1588 enddo
1589 enddo
1590 enddo
1591 C
1592 DO 30 I=1,LENC
1593 TX1(I) = PRJ(I,K+1) * POI(I,K)
1594 QS1(I) = ALF(I,K) + BET(I,K)*POI(I,K)
1595 QOL(I,K) = MIN(QS1(I)*RHMAX,QOI(I,K))
1596
1597 HOL(I,K) = TX1(I)*CP + QOL(I,K)*ALHL
1598 ETA(I,K) = ZERO
1599 TX2(I) = (PRJ(I,K+1) - PRJ(I,K)) * POI(I,K) * CP
1600 30 CONTINUE
1601 C
1602 IF (IC .LT. KM1) THEN
1603 DO 3703 L=KM1,IC1,-1
1604 DO 50 I=1,LENC
1605 QS1(I) = ALF(I,L) + BET(I,L)*POI(I,L)
1606 QOL(I,L) = MIN(QS1(I)*RHMAX,QOI(I,L))
1607 C
1608 TEM1 = TX2(I) + PRJ(I,L+1) * POI(I,L) * CP
1609 HOL(I,L) = TEM1 + QOL(I,L )* ALHL
1610 HST(I,L) = TEM1 + QS1(I) * ALHL
1611
1612 TX1(I) = (PRJ(I,L+1) - PRJ(I,L)) * POI(I,L)
1613 ETA(I,L) = ETA(I,L+1) + TX1(I)*CPBG
1614 TX2(I) = TX2(I) + TX1(I)*CP
1615 50 CONTINUE
1616 C
1617 3703 CONTINUE
1618 ENDIF
1619
1620
1621 DO 70 I=1,LENC
1622 HOL(I,IC) = TX2(I)
1623 QS1(I) = ALF(I,IC) + BET(I,IC)*POI(I,IC)
1624 QOL(I,IC) = MIN(QS1(I)*RHMAX,QOI(I,IC))
1625 c
1626 TEM1 = TX2(I) + PRJ(I,IC1) * POI(I,IC) * CP
1627 HOL(I,IC) = TEM1 + QOL(I,IC) * ALHL
1628 HST(I,IC) = TEM1 + QS1(I) * ALHL
1629 C
1630 TX3(I ) = (PRJ(I,IC1) - PRH(I,IC)) * POI(I,IC)
1631 ETA(I,IC) = ETA(I,IC1) + CPBG * TX3(I)
1632 70 CONTINUE
1633 C
1634 DO 130 I=1,LENC
1635 TX2(I) = HOL(I,K) - HST(I,IC)
1636 TX1(I) = ZERO
1637
1638 130 CONTINUE
1639 C
1640 C ENTRAINMENT PARAMETER
1641 C
1642 DO 160 L=IC,KM1
1643 DO 160 I=1,LENC
1644 TX1(I) = TX1(I) + (HST(I,IC) - HOL(I,L)) * (ETA(I,L) - ETA(I,L+1))
1645 160 CONTINUE
1646 C
1647 LEN1 = 0
1648 LEN2 = 0
1649 ISAV = 0
1650 DO 195 I=1,LENC
1651 IF (TX1(I) .GT. ZERO .AND. TX2(I) .GT. ZERO
1652 . .AND. rhfrac(i).ne.0.0 ) THEN
1653 LEN1 = LEN1 + 1
1654 IA(LEN1) = I
1655 ALM(LEN1) = TX2(I) / TX1(I)
1656 ENDIF
1657 195 CONTINUE
1658 C
1659 LEN2 = LEN1
1660 if (IC1 .lt. K) then
1661 DO 196 I=1,LENC
1662 IF (TX2(I) .LE. 0.0 .AND. (HOL(I,K) .GT. HST(I,IC1))
1663 . .AND. rhfrac(i).ne.0.0 ) THEN
1664 LEN2 = LEN2 + 1
1665 IA(LEN2) = I
1666 ALM(LEN2) = 0.0
1667 ENDIF
1668 196 CONTINUE
1669 endif
1670 C
1671 IF (LEN2 .EQ. 0) THEN
1672 DO 5010 I=1,LENC*K
1673 HST(I,1) = 0.0
1674 QOL(I,1) = 0.0
1675 5010 CONTINUE
1676 DO 5020 I=1,LENC
1677 PCU(I) = 0.0
1678 5020 CONTINUE
1679 RETURN
1680 ENDIF
1681 LEN11 = LEN1 + 1
1682 C
1683 C NORMALIZED MASSFLUX
1684 C
1685 DO 250 I=1,LEN2
1686 ETA(I,K) = 1.0
1687 II = IA(I)
1688 TX2(I) = 0.5 * (PRS(II,IC) + PRS(II,IC1))
1689 TX4(I) = PRS(II,K)
1690 250 CONTINUE
1691 C
1692 DO 252 I=LEN11,LEN2
1693 WFN(I) = 0.0
1694 II = IA(I)
1695 IF (HST(II,IC1) .LT. HST(II,IC)) THEN
1696 TX6(I) = (HST(II,IC1)-HOL(II,K))/(HST(II,IC1)-HST(II,IC))
1697 ELSE
1698 TX6(I) = 0.0
1699 ENDIF
1700 TX2(I) = 0.5 * (PRS(II,IC1)+PRS(II,IC1+1)) * (1.0-TX6(I))
1701 * + TX2(I) * TX6(I)
1702 252 CONTINUE
1703 C
1704 CALL ACRITN(LEN2, TX2, TX4, TX3)
1705 C
1706 DO 260 L=KM1,IC,-1
1707 DO 255 I=1,LEN2
1708 TX1(I) = ETA(IA(I),L)
1709 255 CONTINUE
1710 DO 260 I=1,LEN2
1711 ETA(I,L) = 1.0 + ALM(I) * TX1(I)
1712 260 CONTINUE
1713 C
1714 C CLOUD WORKFUNCTION
1715 C
1716 IF (LEN1 .GT. 0) THEN
1717 DO 270 I=1,LEN1
1718 II = IA(I)
1719 WFN(I) = - GAM(II,IC) * (PRJ(II,IC1) - PRH(II,IC))
1720 * * HST(II,IC) * ETA(I,IC1)
1721 270 CONTINUE
1722 ENDIF
1723 C
1724 DO 290 I=1,LEN2
1725 II = IA(I)
1726 TX1(I) = HOL(II,K)
1727 290 CONTINUE
1728 C
1729 IF (IC1 .LE. KM1) THEN
1730
1731 DO 380 L=KM1,IC1,-1
1732 DO 380 I=1,LEN2
1733 II = IA(I)
1734 TEM = TX1(I) + (ETA(I,L) - ETA(I,L+1)) * HOL(II,L)
1735 C
1736 PCU(I) = PRJ(II,L+1) - PRH(II,L)
1737 TEM1 = ETA(I,L+1) * PCU(I)
1738 TX1(I) = TX1(I)*PCU(I)
1739 C
1740 PCU(I) = PRH(II,L) - PRJ(II,L)
1741 TEM1 = (TEM1 + ETA(I,L) * PCU(I)) * HST(II,L)
1742 TX1(I) = TX1(I) + TEM*PCU(I)
1743 C
1744 WFN(I) = WFN(I) + (TX1(I) - TEM1) * GAM(II,L)
1745 TX1(I) = TEM
1746 380 CONTINUE
1747 ENDIF
1748 C
1749 LENA = 0
1750 IF (LEN1 .GT. 0) THEN
1751 DO 512 I=1,LEN1
1752 II = IA(I)
1753 WFN(I) = WFN(I) + TX1(I) * GAM(II,IC)*(PRJ(II,IC1)-PRH(II,IC))
1754 * - TX3(I)
1755 IF (WFN(I) .GT. 0.0) THEN
1756 LENA = LENA + 1
1757 I1(LENA) = IA(I)
1758 I2(LENA) = I
1759 TX1(LENA) = WFN(I)
1760 TX2(LENA) = QS1(IA(I))
1761 TX6(LENA) = 1.0
1762 ENDIF
1763 512 CONTINUE
1764 ENDIF
1765 LENB = LENA
1766 DO 515 I=LEN11,LEN2
1767 WFN(I) = WFN(I) - TX3(I)
1768 IF (WFN(I) .GT. 0.0 .AND. TX6(I) .GT. 0.0) THEN
1769 LENB = LENB + 1
1770 I1(LENB) = IA(I)
1771 I2(LENB) = I
1772 TX1(LENB) = WFN(I)
1773 TX2(LENB) = QS1(IA(I))
1774 TX4(LENB) = TX6(I)
1775 ENDIF
1776 515 CONTINUE
1777 C
1778 IF (LENB .LE. 0) THEN
1779 DO 5030 I=1,LENC*K
1780 HST(I,1) = 0.0
1781 QOL(I,1) = 0.0
1782 5030 CONTINUE
1783 DO 5040 I=1,LENC
1784 PCU(I) = 0.0
1785 5040 CONTINUE
1786 RETURN
1787 ENDIF
1788
1789 C
1790 DO 516 I=1,LENB
1791 WFN(I) = TX1(I)
1792 QS1(I) = TX2(I)
1793 516 CONTINUE
1794 C
1795 DO 520 L=IC,K
1796 DO 517 I=1,LENB
1797 TX1(I) = ETA(I2(I),L)
1798 517 CONTINUE
1799 DO 520 I=1,LENB
1800 ETA(I,L) = TX1(I)
1801 520 CONTINUE
1802 C
1803 LENA1 = LENA + 1
1804 C
1805 DO 510 I=1,LENA
1806 II = I1(I)
1807 TX8(I) = HST(II,IC) - HOL(II,IC)
1808 510 CONTINUE
1809 DO 530 I=LENA1,LENB
1810 II = I1(I)
1811 TX6(I) = TX4(I)
1812 TEM = TX6(I) * (HOL(II,IC)-HOL(II,IC1)) + HOL(II,IC1)
1813 TX8(I) = HOL(II,K) - TEM
1814
1815 TEM1 = TX6(I) * (QOL(II,IC)-QOL(II,IC1)) + QOL(II,IC1)
1816 TX5(I) = TEM - TEM1 * ALHL
1817 QS1(I) = TEM1 + TX8(I)*(ONE/ALHL)
1818 TX3(I) = HOL(II,IC)
1819 530 CONTINUE
1820 C
1821 C
1822 DO 620 I=1,LENB
1823 II = I1(I)
1824 WLQ(I) = QOL(II,K) - QS1(I) * ETA(I,IC)
1825 TX7(I) = HOL(II,K)
1826 620 CONTINUE
1827 DO NT=1,Ntracer
1828 DO 621 I=1,LENB
1829 II = I1(I)
1830 UHT(I,NT) = UOI(II,K+nltop-1,NT)-UOI(II,IC+nltop-1,NT) * ETA(I,IC)
1831 621 CONTINUE
1832 ENDDO
1833 C
1834 DO 635 L=KM1,IC,-1
1835 DO 630 I=1,LENB
1836 II = I1(I)
1837 TEM = ETA(I,L) - ETA(I,L+1)
1838 WLQ(I) = WLQ(I) + TEM * QOL(II,L)
1839 630 CONTINUE
1840 635 CONTINUE
1841 DO NT=1,Ntracer
1842 DO L=KM1,IC,-1
1843 DO I=1,LENB
1844 II = I1(I)
1845 TEM = ETA(I,L) - ETA(I,L+1)
1846 UHT(I,NT) = UHT(I,NT) + TEM * UOI(II,L+nltop-1,NT)
1847 ENDDO
1848 ENDDO
1849 ENDDO
1850 C
1851 C CALCULATE GS AND PART OF AKM (THAT REQUIRES ETA)
1852 C
1853 DO 690 I=1,LENB
1854 II = I1(I)
1855 c TX7(I) = HOL(II,K)
1856 TEM = (POI(II,KM1) - POI(II,K)) / (PRH(II,K) - PRH(II,KM1))
1857 HOL(I,K) = TEM * (PRJ(II,K)-PRH(II,KM1))*PRH(II,K)*PRI(II,K)
1858 HOL(I,KM1) = TEM * (PRH(II,K)-PRJ(II,K))*PRH(II,KM1)*PRI(II,KM1)
1859 AKM(I) = ZERO
1860 TX2(I) = 0.5 * (PRS(II,IC) + PRS(II,IC1))
1861 690 CONTINUE
1862
1863 IF (IC1 .LE. KM1) THEN
1864 DO 750 L=KM1,IC1,-1
1865 DO 750 I=1,LENB
1866 II = I1(I)
1867 TEM = (POI(II,L-1) - POI(II,L)) * ETA(I,L)
1868 * / (PRH(II,L) - PRH(II,L-1))
1869 C
1870 HOL(I,L) = TEM * (PRJ(II,L)-PRH(II,L-1)) * PRH(II,L)
1871 * * PRI(II,L) + HOL(I,L)
1872 HOL(I,L-1) = TEM * (PRH(II,L)-PRJ(II,L)) * PRH(II,L-1)
1873 * * PRI(II,L-1)
1874 C
1875 AKM(I) = AKM(I) - HOL(I,L)
1876 * * (ETA(I,L) * (PRH(II,L)-PRJ(II,L)) +
1877 * ETA(I,L+1) * (PRJ(II,L+1)-PRH(II,L))) / PRH(II,L)
1878 750 CONTINUE
1879 ENDIF
1880 C
1881 C
1882 CALL RNCL(LENB, TX2, TX1, CLF)
1883
1884 DO 770 I=1,LENB
1885 TX2(I) = (ONE - TX1(I)) * WLQ(I)
1886 WLQ(I) = TX1(I) * WLQ(I)
1887 C
1888 TX1(I) = HOL(I,IC)
1889 770 CONTINUE
1890 DO 790 I=LENA1, LENB
1891 II = I1(I)
1892 TX1(I) = TX1(I) + (TX5(I)-TX3(I)+QOL(II,IC)*ALHL)*(PRI(II,IC)/CP)
1893 790 CONTINUE
1894
1895 DO 800 I=1,LENB
1896 HOL(I,IC) = TX1(I) - TX2(I) * ALBCP * PRI(I1(I),IC)
1897 800 CONTINUE
1898
1899 IF (LENA .GT. 0) THEN
1900 DO 810 I=1,LENA
1901 II = I1(I)
1902 AKM(I) = AKM(I) - ETA(I,IC1) * (PRJ(II,IC1) - PRH(II,IC))
1903 * * TX1(I) / PRH(II,IC)
1904 810 CONTINUE
1905 ENDIF
1906 c
1907 C CALCULATE GH
1908 C
1909 DO 830 I=1,LENB
1910 II = I1(I)
1911 TX3(I) = QOL(II,KM1) - QOL(II,K)
1912 GMH(I,K) = HOL(I,K) + TX3(I) * PRI(II,K) * (ALBCP)
1913
1914 AKM(I) = AKM(I) + GAM(II,KM1)*(PRJ(II,K)-PRH(II,KM1))
1915 * * GMH(I,K)
1916 TX3(I) = zero
1917 830 CONTINUE
1918 C
1919 IF (IC1 .LE. KM1) THEN
1920 DO 840 L=KM1,IC1,-1
1921 DO 840 I=1,LENB
1922 II = I1(I)
1923 TX2(I) = TX3(I)
1924 TX3(I) = (QOL(II,L-1) - QOL(II,L)) * ETA(I,L)
1925 TX2(I) = TX2(I) + TX3(I)
1926 C
1927 GMH(I,L) = HOL(I,L) + TX2(I) * PRI(II,L) * (ALBCP*HALF)
1928 840 CONTINUE
1929 C
1930 C
1931 ENDIF
1932 DO 850 I=LENA1,LENB
1933 TX3(I) = TX3(I) + TWOBAL
1934 * * (TX7(I) - TX8(I) - TX5(I) - QOL(I1(I),IC)*ALHL)
1935 850 CONTINUE
1936 DO 860 I=1,LENB
1937 GMH(I,IC) = TX1(I) + PRI(I1(I),IC) * ONEBCP
1938 * * (TX3(I)*(ALHL*HALF) + ETA(I,IC) * TX8(I))
1939 860 CONTINUE
1940 C
1941 C CALCULATE HC PART OF AKM
1942 C
1943 IF (IC1 .LE. KM1) THEN
1944 DO 870 I=1,LENB
1945 TX1(I) = GMH(I,K)
1946 870 CONTINUE
1947 DO 3725 L=KM1,IC1,-1
1948 DO 880 I=1,LENB
1949 II = I1(I)
1950 TX1(I) = TX1(I) + (ETA(I,L) - ETA(I,L+1)) * GMH(I,L)
1951 TX2(I) = GAM(II,L-1) * (PRJ(II,L) - PRH(II,L-1))
1952 880 CONTINUE
1953 C
1954 IF (L .EQ. IC1) THEN
1955 DO 890 I=LENA1,LENB
1956 TX2(I) = ZERO
1957 890 CONTINUE
1958 ENDIF
1959 DO 900 I=1,LENB
1960 II = I1(I)
1961 AKM(I) = AKM(I) + TX1(I) *
1962 * (TX2(I) + GAM(II,L)*(PRH(II,L)-PRJ(II,L)))
1963 900 CONTINUE
1964 3725 CONTINUE
1965 ENDIF
1966 C
1967 DO 920 I=LENA1,LENB
1968 II = I1(I)
1969 TX2(I) = 0.5 * (PRS(II,IC) + PRS(II,IC1))
1970 * + 0.5*(PRS(II,IC+2) - PRS(II,IC)) * (ONE-TX6(I))
1971 c
1972 TX1(I) = PRS(II,IC1)
1973 TX5(I) = 0.5 * (PRS(II,IC1) + PRS(II,IC+2))
1974 C
1975 IF ((TX2(I) .GE. TX1(I)) .AND. (TX2(I) .LT. TX5(I))) THEN
1976 TX6(I) = ONE - (TX2(I) - TX1(I)) / (TX5(I) - TX1(I))
1977 C
1978 TEM = PRI(II,IC1) / PRI(II,IC)
1979 HOL(I,IC1) = HOL(I,IC1) + HOL(I,IC) * TEM
1980 HOL(I,IC) = ZERO
1981 C
1982 GMH(I,IC1) = GMH(I,IC1) + GMH(I,IC) * TEM
1983 GMH(I,IC) = ZERO
1984 ELSEIF (TX2(I) .LT. TX1(I)) THEN
1985 TX6(I) = 1.0
1986 ELSE
1987 TX6(I) = 0.0
1988 ENDIF
1989 920 CONTINUE
1990 C
1991 C
1992 DO I=1,LENC
1993 PCU(I) = 0.0
1994 ENDDO
1995
1996 DO 970 I=1,LENB
1997 II = I1(I)
1998 IF (AKM(I) .LT. ZERO .AND. WLQ(I) .GE. 0.0) THEN
1999 WFN(I) = - TX6(I) * WFN(I) * RASALF / AKM(I)
2000 ELSE
2001 WFN(I) = ZERO
2002 ENDIF
2003 TEM = (PRS(II,K+1)-PRS(II,K))*(CMB2PA*FRAC)
2004 WFN(I) = MIN(WFN(I), TEM)
2005 C
2006 C compute cloud amount
2007 C
2008 CC TX1(I) = CLN(II)
2009 CC IF (WFN(I) .GT. CRTMSF) TX1(I) = TX1(I) + CLF(I)
2010 CC IF (TX1(I) .GT. ONE) TX1(I) = ONE
2011 C
2012 C PRECIPITATION
2013 C
2014 PCU(II) = WLQ(I) * WFN(I) * ONEBG
2015 C
2016 C CUMULUS FRICTION AT THE BOTTOM LAYER
2017 C
2018 TX4(I) = WFN(I) * (1.0/ALHL)
2019 TX5(I) = WFN(I) * ONEBCP
2020 970 CONTINUE
2021 C
2022 C compute cloud mass flux for diagnostic output
2023 C
2024 DO L = IC,K
2025 DO I=1,LENB
2026 II = I1(I)
2027 if(L.lt.K)then
2028 CMASS(II,L) = ETA(I,L+1) * WFN(I) * ONEBG
2029 else
2030 CMASS(II,L) = WFN(I) * ONEBG
2031 endif
2032 ENDDO
2033 ENDDO
2034 C
2035 CC DO 975 I=1,LENB
2036 CC II = I1(I)
2037 CC CLN(II) = TX1(I)
2038 CC975 CONTINUE
2039 C
2040 C THETA AND Q CHANGE DUE TO CLOUD TYPE IC
2041 C
2042
2043 c TEMA = 0.0
2044 c TEMB = 0.0
2045 DO 990 L=IC,K
2046 DO 980 I=1,LENB
2047 II = I1(I)
2048 TEM = (GMH(I,L) - HOL(I,L)) * TX4(I)
2049 TEM1 = HOL(I,L) * TX5(I)
2050 C
2051 TCU(II,L) = TEM1 / PRH(II,L)
2052 QCU(II,L) = TEM
2053 980 CONTINUE
2054
2055 c I = I1(IP1)
2056 c
2057 c TEM = (PRS(I,L+1)-PRS(I,L)) * (ONEBG*100.0)
2058 c TEMA = TEMA + TCU(I,L) * PRH(I,L) * TEM * (CP/ALHL)
2059 c TEMB = TEMB + QCU(I,L) * TEM
2060 C
2061 990 CONTINUE
2062 C
2063 c Compute Tracer Tendencies
2064 c -------------------------
2065 do nt = 1,ntracer
2066
2067 c Tracer Tendency at the Bottom Layer
2068 c -----------------------------------
2069 DO 995 I=1,LENB
2070 II = I1(I)
2071 TEM = half*TX5(I) * PRI(II,K)
2072 TX1(I) = (UOI(II,KM1+nltop-1,nt) - UOI(II,K+nltop-1,nt))
2073 ucu(II,K,nt) = TEM * TX1(I)
2074 995 CONTINUE
2075
2076 c Tracer Tendency at all other Levels
2077 c -----------------------------------
2078 DO 1020 L=KM1,IC1,-1
2079 DO 1010 I=1,LENB
2080 II = I1(I)
2081 TEM = half*TX5(I) * PRI(II,L)
2082 TEM1 = TX1(I)
2083 TX1(I) = (UOI(II,L-1+nltop-1,nt)-UOI(II,L+nltop-1,nt)) * ETA(I,L)
2084 TX3(I) = (TX1(I) + TEM1) * TEM
2085 1010 CONTINUE
2086 DO 1020 I=1,LENB
2087 II = I1(I)
2088 ucu(II,L,nt) = TX3(I)
2089 1020 CONTINUE
2090
2091 DO 1030 I=1,LENB
2092 II = I1(I)
2093 IF (TX6(I) .GE. 1.0) THEN
2094 TEM = half*TX5(I) * PRI(II,IC)
2095 ELSE
2096 TEM = 0.0
2097 ENDIF
2098 TX1(I) = (TX1(I) + UHT(I,nt) + UHT(I,nt)) * TEM
2099 1030 CONTINUE
2100 DO 1040 I=1,LENB
2101 II = I1(I)
2102 ucu(II,IC,nt) = TX1(I)
2103 1040 CONTINUE
2104
2105 enddo
2106 C
2107 C PENETRATIVE CONVECTION CALCULATION OVER
2108 C
2109
2110 RETURN
2111 END
2112 SUBROUTINE RNCL(lng, PL, RNO, CLF)
2113 C
2114 C*********************************************************************
2115 C********************** Relaxed Arakawa-Schubert *********************
2116 C************************ SUBROUTINE RNCL ************************
2117 C**************************** 23 July 1992 ***************************
2118 C*********************************************************************
2119 implicit none
2120 C Argument List declarations
2121 integer lng
2122 _RL PL(lng), RNO(lng), CLF(lng)
2123
2124 C Local Variables
2125 _RL p5,p8,pt8,pt2,pfac,p4,p6,p7,p9,cucld,cfac
2126 PARAMETER (P5=500.0, P8=800.0, PT8=0.8, PT2=0.2)
2127 PARAMETER (PFAC=PT2/(P8-P5))
2128 PARAMETER (P4=400.0, P6=401.0)
2129 PARAMETER (P7=700.0, P9=900.0)
2130 PARAMETER (CUCLD=0.5,CFAC=CUCLD/(P6-P4))
2131
2132 integer i
2133 C
2134 DO 10 I=1,lng
2135 rno(i) = 1.0
2136 ccc if( pl(i).le.400.0 ) rno(i) = max( 0.75, 1.0-0.0025*(400.0-pl(i)) )
2137
2138 ccc IF ( PL(I).GE.P7 .AND. PL(I).LE.P9 ) THEN
2139 ccc RNO(I) = ((P9-PL(I))/(P9-P7)) **2
2140 ccc ELSE IF (PL(I).GT.P9) THEN
2141 ccc RNO(I) = 0.
2142 ccc ENDIF
2143
2144 CLF(I) = CUCLD
2145 C
2146 CARIESIF (PL(I) .GE. P5 .AND. PL(I) .LE. P8) THEN
2147 CARIES RNO(I) = (P8-PL(I))*PFAC + PT8
2148 CARIESELSEIF (PL(I) .GT. P8 ) THEN
2149 CARIES RNO(I) = PT8
2150 CARIESENDIF
2151 CARIES
2152 IF (PL(I) .GE. P4 .AND. PL(I) .LE. P6) THEN
2153 CLF(I) = (P6-PL(I))*CFAC
2154 ELSEIF (PL(I) .GT. P6 ) THEN
2155 CLF(I) = 0.0
2156 ENDIF
2157 10 CONTINUE
2158 C
2159 RETURN
2160 END
2161 SUBROUTINE ACRITN ( lng,PL,PLB,ACR )
2162
2163 C*********************************************************************
2164 C********************** Relaxed Arakawa-Schubert *********************
2165 C************************** SUBROUTINE ACRIT *********************
2166 C****************** modified August 28, 1996 L.Takacs ************
2167 C**** *****
2168 C**** Note: Data obtained from January Mean After-Analysis *****
2169 C**** from 4x5 46-layer GEOS Assimilation *****
2170 C**** *****
2171 C*********************************************************************
2172 implicit none
2173 C Argument List declarations
2174 integer lng
2175 _RL PL(lng), PLB(lng), ACR(lng)
2176
2177 C Local variables
2178 integer lma
2179 parameter (lma=18)
2180 _RL p(lma)
2181 _RL a(lma)
2182 integer i,L
2183 _RL temp
2184
2185 data p / 93.81, 111.65, 133.46, 157.80, 186.51,
2186 . 219.88, 257.40, 301.21, 352.49, 409.76,
2187 . 471.59, 535.04, 603.33, 672.79, 741.12,
2188 . 812.52, 875.31, 930.20/
2189
2190 data a / 3.35848, 3.13645, 2.48072, 2.08277, 1.53364,
2191 . 1.01971, .65846, .45867, .38687, .31002,
2192 . .25574, .20347, .17254, .15260, .16756,
2193 . .09916, .10360, .05880/
2194
2195
2196 do L=1,lma-1
2197 do i=1,lng
2198 if( pl(i).ge.p(L) .and.
2199 . pl(i).le.p(L+1)) then
2200 temp = ( pl(i)-p(L) )/( p(L+1)-p(L) )
2201 acr(i) = a(L+1)*temp + a(L)*(1-temp)
2202 endif
2203 enddo
2204 enddo
2205
2206 do i=1,lng
2207 if( pl(i).lt.p(1) ) acr(i) = a(1)
2208 if( pl(i).gt.p(lma) ) acr(i) = a(lma)
2209 enddo
2210
2211 do i=1,lng
2212 acr(i) = acr(i) * (plb(i)-pl(i))
2213 enddo
2214
2215 RETURN
2216 END
2217 SUBROUTINE RNEVP(NN,IRUN,NLAY,TL,QL,RAIN,PL,CLFRAC,SP,DP,PLKE,
2218 1 PLK,TH,TEMP1,TEMP2,TEMP3,ITMP1,ITMP2,RCON,RLAR,CLSBTH,tmscl,
2219 2 tmfrc,cp,gravity,alhl,gamfac,cldlz,RHCRIT,offset,alpha)
2220
2221 implicit none
2222 C Argument List declarations
2223 integer nn,irun,nlay
2224 _RL TL(IRUN,NLAY),QL(IRUN,NLAY),RAIN(IRUN,NLAY),
2225 . PL(IRUN,NLAY),CLFRAC(IRUN,NLAY),SP(IRUN),TEMP1(IRUN,NLAY),
2226 . TEMP2(IRUN,NLAY),PLKE(IRUN,NLAY+1),
2227 . RCON(IRUN),RLAR(IRUN),DP(IRUN,NLAY),PLK(IRUN,NLAY),TH(IRUN,NLAY),
2228 . TEMP3(IRUN,NLAY)
2229 integer ITMP1(IRUN,NLAY),ITMP2(IRUN,NLAY)
2230 _RL CLSBTH(IRUN,NLAY)
2231 _RL tmscl,tmfrc,cp,gravity,alhl,gamfac,offset,alpha
2232 _RL cldlz(irun,nlay)
2233 _RL rhcrit(irun,nlay)
2234 C
2235 C Local Variables
2236 _RL zm1p04,zero,two89,zp44,zp01,half,zp578,one,thousand,z3600
2237 _RL zp1,zp001
2238 PARAMETER (ZM1P04 = -1.04E-4 )
2239 PARAMETER (ZERO = 0.)
2240 PARAMETER (TWO89= 2.89E-5)
2241 PARAMETER ( ZP44= 0.44)
2242 PARAMETER ( ZP01= 0.01)
2243 PARAMETER ( ZP1 = 0.1 )
2244 PARAMETER ( ZP001= 0.001)
2245 PARAMETER ( HALF= 0.5)
2246 PARAMETER ( ZP578 = 0.578 )
2247 PARAMETER ( ONE = 1.)
2248 PARAMETER ( THOUSAND = 1000.)
2249 PARAMETER ( Z3600 = 3600.)
2250 C
2251 _RL EVP9(IRUN,NLAY)
2252 _RL water(irun),crystal(irun)
2253 _RL watevap(irun),iceevap(irun)
2254 _RL fracwat,fracice, tice,rh,fact,dum
2255 _RL rainmax(irun)
2256 _RL getcon,rphf,elocp,cpog,relax
2257 _RL exparg,arearat,rpow
2258
2259 integer i,L,n,nlaym1,irnlay,irnlm1
2260
2261 c Explicit Inline Directives
2262 c --------------------------
2263 #ifdef CRAY
2264 #ifdef f77
2265 cfpp$ expand (qsat)
2266 #endif
2267 #endif
2268
2269 tice = getcon('FREEZING-POINT')
2270
2271 fracwat = 0.70
2272 fracice = 0.01
2273
2274 NLAYM1 = NLAY - 1
2275 IRNLAY = IRUN*NLAY
2276 IRNLM1 = IRUN*(NLAY-1)
2277
2278 RPHF = Z3600/tmscl
2279
2280 ELOCP = alhl/cp
2281 CPOG = cp/gravity
2282
2283 DO I = 1,IRUN
2284 RLAR(I) = 0.
2285 water(i) = 0.
2286 crystal(i) = 0.
2287 ENDDO
2288
2289 do L = 1,nlay
2290 do i = 1,irun
2291 EVP9(i,L) = 0.
2292 TEMP1(i,L) = 0.
2293 TEMP2(i,L) = 0.
2294 TEMP3(i,L) = 0.
2295 CLSBTH(i,L) = 0.
2296 cldlz(i,L) = 0.
2297 enddo
2298 enddo
2299
2300 C RHO(ZERO) / RHO FOR TERMINAL VELOCITY APPROX.
2301 c ---------------------------------------------
2302 DO L = 1,NLAY
2303 DO I = 1,IRUN
2304 TEMP2(I,L) = PL(I,L)*ZP001
2305 TEMP2(I,L) = SQRT(TEMP2(I,L))
2306 ENDDO
2307 ENDDO
2308
2309 C INVERSE OF MASS IN EACH LAYER
2310 c -----------------------------
2311 DO L = 1,NLAY
2312 DO I = 1,IRUN
2313 TEMP3(I,L) = GRAVITY*ZP01 / DP(I,L)
2314 ENDDO
2315 ENDDO
2316
2317 C DO LOOP FOR MOISTURE EVAPORATION ABILITY AND CONVEC EVAPORATION.
2318 c ----------------------------------------------------------------
2319 DO 100 L=1,NLAY
2320
2321 DO I = 1,IRUN
2322 TEMP1(I,3) = TL(I,L)
2323 TEMP1(I,4) = QL(I,L)
2324 ENDDO
2325
2326 DO 50 N=1,2
2327 IF(N.EQ.1)RELAX=HALF
2328 IF(N.GT.1)RELAX=ONE
2329
2330 DO I = 1,IRUN
2331 call qsat ( temp1(i,3),pl(i,L),temp1(i,2),temp1(i,6),.true. )
2332 TEMP1(I,5)=TEMP1(I,2)-TEMP1(I,4)
2333 TEMP1(I,6)=TEMP1(I,6)*ELOCP
2334 TEMP1(I,5)=TEMP1(I,5)/(ONE+TEMP1(I,6))
2335 TEMP1(I,4)=TEMP1(I,4)+TEMP1(I,5)*RELAX
2336 TEMP1(I,3)=TEMP1(I,3)-TEMP1(I,5)*ELOCP*RELAX
2337 ENDDO
2338 50 CONTINUE
2339
2340 DO I = 1,IRUN
2341 EVP9(I,L) = (TEMP1(I,4) - QL(I,L))/TEMP3(I,L)
2342 C convective detrained water
2343 cldlz(i,L) = rain(i,L)*temp3(i,L)
2344 if( tl(i,L).gt.tice-20.) then
2345 water(i) = water(i) + rain(i,L)
2346 else
2347 crystal(i) = crystal(i) + rain(i,L)
2348 endif
2349 ENDDO
2350
2351 C**********************************************************************
2352 C FOR CONVECTIVE PRECIP, FIND THE "EVAPORATION EFFICIENCY" USING *
2353 C KESSLERS PARAMETERIZATION *
2354 C**********************************************************************
2355
2356 DO 20 I=1,IRUN
2357
2358 iceevap(i) = 0.
2359 watevap(i) = 0.
2360
2361 if( (evp9(i,L).gt.0.) .and. (crystal(i).gt.0.) ) then
2362 iceevap(I) = EVP9(I,L)*fracice
2363 IF(iceevap(i).GE.crystal(i)) iceevap(i) = crystal(i)
2364 EVP9(I,L)=EVP9(I,L)-iceevap(I)
2365 crystal(i) = crystal(i) - iceevap(i)
2366 endif
2367
2368 C and now warm precipitate
2369 if( (evp9(i,L).gt.0.) .and. (water(i).gt.0.) ) then
2370 exparg = ZM1P04*tmscl*((water(i)*RPHF*TEMP2(I,L))**ZP578)
2371 AREARAT = ONE-(EXP(EXPARG))
2372 watevap(I) = EVP9(I,L)*AREARAT*fracwat
2373 IF(watevap(I).GE.water(i)) watevap(I) = water(i)
2374 EVP9(I,L)=EVP9(I,L)-watevap(I)
2375 water(i) = water(i) - watevap(i)
2376 endif
2377
2378 QL(I,L) = QL(I,L)+(iceevap(i)+watevap(i))*TEMP3(I,L)
2379 TL(I,L) = TL(I,L)-(iceevap(i)+watevap(i))*TEMP3(I,L)*ELOCP
2380
2381 20 CONTINUE
2382
2383 100 CONTINUE
2384
2385 do i = 1,irun
2386 rcon(i) = water(i) + crystal(i)
2387 enddo
2388
2389 C**********************************************************************
2390 C Large Scale Precip
2391 C**********************************************************************
2392
2393 DO 200 L=1,NLAY
2394 DO I = 1,IRUN
2395 rainmax(i) = rhcrit(i,L)*evp9(i,L) +
2396 . ql(i,L)*(rhcrit(i,L)-1.)/temp3(i,L)
2397
2398 if (rainmax(i).LE.0.0) then
2399 call qsat( tl(i,L),pl(i,L),rh,dum,.false.)
2400 rh = ql(i,L)/rh
2401
2402 if( rhcrit(i,L).eq.1.0 ) then
2403 fact = 1.0
2404 else
2405 fact = min( 1.0, alpha + (1.0-alpha)*( rh-rhcrit(i,L)) /
2406 1 (1.0-rhcrit(i,L)) )
2407 endif
2408
2409 C Do not allow clouds above 10 mb
2410 if( pl(i,L).ge.10.0 ) CLSBTH(I,L) = fact
2411 RLAR(I) = RLAR(I)-rainmax(I)
2412 QL(I,L) = QL(I,L)+rainmax(I)*TEMP3(I,L)
2413 TL(I,L) = TL(I,L)-rainmax(I)*TEMP3(I,L)*ELOCP
2414 C Large-scale water
2415 cldlz(i,L) = cldlz(i,L) - rainmax(i)*temp3(i,L)
2416 ENDIF
2417 ENDDO
2418
2419 DO I=1,IRUN
2420 IF((RLAR(I).GT.0.0).AND.(rainmax(I).GT.0.0))THEN
2421 RPOW=(RLAR(I)*RPHF*TEMP2(I,L))**ZP578
2422 EXPARG = ZM1P04*tmscl*RPOW
2423 AREARAT = ONE-(EXP(EXPARG))
2424 TEMP1(I,7) = rainmax(I)*AREARAT
2425 IF(TEMP1(I,7).GE.RLAR(I)) TEMP1(I,7) = RLAR(I)
2426 RLAR(I) = RLAR(I)-TEMP1(I,7)
2427 QL(I,L) = QL(I,L)+TEMP1(I,7)*TEMP3(I,L)
2428 TL(I,L) = TL(I,L)-TEMP1(I,7)*TEMP3(I,L)*ELOCP
2429 ENDIF
2430 ENDDO
2431
2432 200 CONTINUE
2433
2434 RETURN
2435 END
2436 subroutine srclouds (th,q,plk,pl,plke,cloud,cldwat,irun,irise,
2437 1 rhc,offset,alpha)
2438 C***********************************************************************
2439 C
2440 C PURPOSE:
2441 C ========
2442 C Compute non-precipitating cloud fractions
2443 C based on Slingo and Ritter (1985).
2444 C Remove cloudiness where conditionally unstable.
2445 C
2446 C INPUT:
2447 C ======
2448 C th ......... Potential Temperature (irun,irise)
2449 C q .......... Specific Humidity (irun,irise)
2450 C plk ........ P**Kappa at mid-layer (irun,irise)
2451 C pl ......... Pressure at mid-layer (irun,irise)
2452 C plke ....... P**Kappa at edge (irun,irise+1)
2453 C irun ....... Horizontal dimension
2454 C irise ...... Vertical dimension
2455 C
2456 C OUTPUT:
2457 C =======
2458 C cloud ...... Cloud Fraction (irun,irise)
2459 C
2460 C***********************************************************************
2461
2462 implicit none
2463 integer irun,irise
2464
2465 _RL th(irun,irise)
2466 _RL q(irun,irise)
2467 _RL plk(irun,irise)
2468 _RL pl(irun,irise)
2469 _RL plke(irun,irise+1)
2470
2471 _RL cloud(irun,irise)
2472 _RL cldwat(irun,irise)
2473 _RL qs(irun,irise)
2474
2475 _RL cp, alhl, getcon, akap
2476 _RL ratio, temp, elocp
2477 _RL rhcrit,rh,dum
2478 integer i,L
2479
2480 _RL rhc(irun,irise)
2481 _RL offset,alpha
2482
2483 c Explicit Inline Directives
2484 c --------------------------
2485 #ifdef CRAY
2486 #ifdef f77
2487 cfpp$ expand (qsat)
2488 #endif
2489 #endif
2490
2491 cp = getcon('CP')
2492 alhl = getcon('LATENT HEAT COND')
2493 elocp = alhl/cp
2494 akap = getcon('KAPPA')
2495
2496 do L = 1,irise
2497 do i = 1,irun
2498 temp = th(i,L)*plk(i,L)
2499 call qsat ( temp,pl(i,L),qs(i,L),dum,.false. )
2500 enddo
2501 enddo
2502
2503 do L = 2,irise
2504 do i = 1,irun
2505 rh = q(i,L)/qs(i,L)
2506
2507 rhcrit = rhc(i,L) - offset
2508 ratio = alpha*(rh-rhcrit)/offset
2509
2510 if(cloud(i,L).eq. 0.0 .and. ratio.gt.0.0 ) then
2511 cloud(i,L) = min( ratio,1.0 )
2512 endif
2513
2514 enddo
2515 enddo
2516
2517 c Reduce clouds from conditionally unstable layer
2518 c -----------------------------------------------
2519 call ctei ( th,q,cloud,cldwat,pl,plk,plke,irun,irise )
2520
2521 return
2522 end
2523
2524 subroutine ctei ( th,q,cldfrc,cldwat,pl,plk,plke,im,lm )
2525 implicit none
2526 integer im,lm
2527 _RL th(im,lm),q(im,lm),plke(im,lm+1),cldwat(im,lm)
2528 _RL plk(im,lm),pl(im,lm),cldfrc(im,lm)
2529 integer i,L
2530 _RL getcon,cp,alhl,elocp,cpoel,t,p,s,qs,dqsdt,dq
2531 _RL k,krd,kmm,f
2532
2533 cp = getcon('CP')
2534 alhl = getcon('LATENT HEAT COND')
2535 cpoel = cp/alhl
2536 elocp = alhl/cp
2537
2538 do L=lm,2,-1
2539 do i=1,im
2540 dq = q(i,L)+cldwat(i,L)-q(i,L-1)-cldwat(i,L-1)
2541 if( dq.eq.0.0 ) dq = 1.0e-20
2542 k = 1.0 + cpoel*plke(i,L)*( th(i,L)-th(i,L-1) ) / dq
2543
2544 t = th(i,L)*plk(i,L)
2545 p = pl(i,L)
2546 call qsat ( t,p,qs,dqsdt,.true. )
2547
2548 krd = ( cpoel*t*(1+elocp*dqsdt) )/( 1 + 1.608*dqsdt*t )
2549
2550 kmm = ( 1+elocp*dqsdt )*( 1 + 0.392*cpoel*t )
2551 . / ( 2+(1+1.608*cpoel*t)*elocp*dqsdt )
2552
2553 s = ( (k-krd)/(kmm-krd) )
2554 f = 1.0 - min( 1.0, max(0.0,1.0-exp(-s)) )
2555
2556 cldfrc(i,L) = cldfrc(i,L)*f
2557 cldwat(i,L) = cldwat(i,L)*f
2558
2559 enddo
2560 enddo
2561
2562 return
2563 end
2564
2565 subroutine back2grd(gathered,indeces,scattered,irun)
2566 implicit none
2567 integer i,irun,indeces(irun)
2568 _RL gathered(irun),scattered(irun)
2569 _RL temp(irun)
2570 do i = 1,irun
2571 temp(indeces(i)) = gathered(i)
2572 enddo
2573 do i = 1,irun
2574 scattered(i) = temp(i)
2575 enddo
2576 return
2577 end

  ViewVC Help
Powered by ViewVC 1.1.22