/[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.33 - (show annotations) (download)
Mon May 23 02:09:58 2005 UTC (19 years, 1 month ago) by jmc
Branch: MAIN
Changes since 1.32: +3 -3 lines
fix 2 diagnostics_fill calls.

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

  ViewVC Help
Powered by ViewVC 1.1.22