/[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.38 - (show annotations) (download)
Thu Apr 2 21:31:05 2009 UTC (15 years, 3 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f, checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63g, checkpoint63a, checkpoint63b, checkpoint63c, checkpoint64, checkpoint65, checkpoint62, checkpoint63, checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65o, checkpoint65l, checkpoint65m, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65g, checkpoint65d, checkpoint65e, checkpoint62c, checkpoint62b, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62w, checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62z, checkpoint62y, checkpoint62x, checkpoint61n, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61q, checkpoint61z, checkpoint61x, checkpoint61y, HEAD
Changes since 1.37: +11 -5 lines
forgot this one from previous check-in.

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

  ViewVC Help
Powered by ViewVC 1.1.22