/[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.30 - (show annotations) (download)
Wed Mar 9 23:25:18 2005 UTC (19 years, 4 months ago) by molod
Branch: MAIN
Changes since 1.29: +1 -29 lines
Get rid of UDIAG fills for moist u and v tendency

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

  ViewVC Help
Powered by ViewVC 1.1.22