CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
calpreciptype.f90
1
3
7contains
8
10 subroutine calpreciptype(kdt,nrcm,im,ix,lm,lp1,randomno, &
11 xlat,xlon, &
12 gt0,gq0,prsl,prsi,prec, & !input
13 phii,tskin, & !input
14 domr,domzr,domip,doms) !output
15
16!$$$ subprogram documentation block
17! . . .
18! subprogram: calpreciptype compute dominant precip type
19! prgrmmr: chuang org: w/np2 date: 2008-05-28
20!
21!
22! abstract:
23! this routine computes precipitation type.
24! . it is adopted from post but was made into a column to used by gfs model
25!
26! --------------------------------------------------------------------
27 use funcphys, only : fpvs,ftdp,fpkap,ftlcl,stma,fthe
28 use physcons
29 use machine , only : kind_phys
30!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
31 implicit none
32!
33 real(kind=kind_phys), parameter :: pthresh = 0.0, oneog = 1.0/con_g
34 integer,parameter :: nalg = 5
35!
36! declare variables.
37!
38 integer,intent(in) :: kdt,nrcm,im,ix,lm,lp1
39 real(kind=kind_phys),intent(in) :: xlat(im),xlon(im)
40 real(kind=kind_phys),intent(in) :: randomno(ix,nrcm)
41 real(kind=kind_phys),dimension(im), intent(in) :: prec,tskin
42 real(kind=kind_phys),dimension(ix,lm), intent(in) :: gt0,gq0,prsl
43 real(kind=kind_phys),dimension(ix,lp1),intent(in) :: prsi,phii
44 real(kind=kind_phys),dimension(im), intent(out) :: domr,domzr,domip,doms
45
46 integer, dimension(nalg) :: sleet,rain,freezr,snow
47 real(kind=kind_phys),dimension(lm) :: t,q,pmid
48 real(kind=kind_phys),dimension(lp1) :: pint,zint
49 real(kind=kind_phys), allocatable :: twet(:),rh(:),td(:)
50!
51 integer i,iwx,isno,iip,izr,irain,k,k1
52 real(kind=kind_phys) es,qc,pv,tdpd,pr,tr,pk,tlcl,thelcl,qwet, &
53 time_vert,time_ncep,time_ramer,time_bourg,time_revised,&
54 time_dominant,btim,timef,ranl(2)
55
56!
57! computes wet bulb here since two algorithms use it
58! lp1=lm+1
59! convert geopotential to height
60! do l=1,lp1
61! zint(l)=zint(l)/con_g
62! end do
63! don't forget to flip 3d arrays around because gfs counts from bottom up
64
65 allocate ( twet(lm),rh(lm),td(lm) )
66
67! print*,'debug calpreciptype: ', im,lm,lp1,nrcm
68
69! time_vert = 0.
70! time_ncep = 0.
71! time_ramer = 0.
72! time_bourg = 0.
73! time_revised = 0.
74
75 do i=1,im
76 if (prec(i) > pthresh) then
77 do k=1,lm
78 k1 = lm-k+1
79 t(k1) = gt0(i,k)
80 q(k1) = gq0(i,k)
81 pmid(k1) = prsl(i,k) ! pressure in pascals
82!
84!
85 pv = pmid(k1)*q(k1)/(con_eps-con_epsm1*q(k1))
86 td(k1) = ftdp(pv)
87 tdpd = t(k1)-td(k1)
88! if (pmid(k1) >= 50000.) then ! only compute twet below 500mb to save time
89 if (tdpd > 0.) then
90 pr = pmid(k1)
91 tr = t(k1)
92 pk = fpkap(pr)
93 tlcl = ftlcl(tr,tdpd)
94 thelcl = fthe(tlcl,pk*tlcl/tr)
95 call stma(thelcl,pk,twet(k1),qwet)
96 else
97 twet(k1) = t(k1)
98 endif
99! endif
100 es = min(fpvs(t(k1)), pmid(k1))
101 qc = con_eps*es / (pmid(k1)+con_epsm1*es)
102 rh(k1) = max(con_epsq,q(k1)) / qc
103
104 k1 = lp1-k+1
105 pint(k1) = prsi(i,k)
106 zint(k1) = phii(i,k) * oneog
107
108 enddo
109 pint(1) = prsi(i,lp1)
110 zint(1) = phii(i,lp1) * oneog
111
112!-------------------------------------------------------------------------------
113! if(kdt>15.and.kdt<20) time_vert = time_vert + (timef() - btim)
114! debug print statement
115! if (abs(xlon(i)*57.29578-114.0) .lt. 0.2 .and. &
116! abs(xlat(i)*57.29578-40.0) .lt. 0.2)then
117! print*,'debug in calpreciptype: i,im,lm,lp1,xlon,xlat,prec,tskin,nrcm,randomno', &
118! i,im,lm,lp1,xlon(i)*57.29578,xlat(i)*57.29578,prec(i),tskin(i),, &
119! nrcm,randomno(i,1:nrcm)
120! do l=1,lm
121! print*,'debug in calpreciptype: l,t,q,p,pint,z,twet', &
122! l,t(l),q(l), &
123! pmid(l),pint(l),zint(l),twet(l)
124! end do
125! print*,'debug in calpreciptype: lp1,pint,z ', lp1,pint(lp1),zint(lp1)
126! end if
127! end debug print statement
128! call wetbulb(lm,con_rocp,con_epsq,t,q,pmid,twet)
129! if(kdt>10.and.kdt<20)btim = timef()
130!-------------------------------------------------------------------------------
131!
132! instantaneous precipitation type.
133
135 call calwxt(lm,lp1,t,q,pmid,pint,con_fvirt,con_rog,con_epsq,zint,iwx,twet)
136 snow(1) = mod(iwx,2)
137 sleet(1) = mod(iwx,4)/2
138 freezr(1) = mod(iwx,8)/4
139 rain(1) = iwx/8
140
141! dominant precipitation type
142
143!gsm if dominant precip type is requested, 4 more algorithms
144!gsm will be called. the tallies are then summed in calwxt_dominant
145
146! ramer algorithm
147! allocate ( rh(lm),td(lm) )
148! do l=1,lm
149!hc: use rh and td consistent with gfs ice physics
150! es=fpvs(t(l))
151! es=min(es,pmid(l))
152! qc=con_eps*es/(pmid(l)+con_epsm1*es)
153! rh(l)=max(con_epsq,q(l))/qc
154! pv = pmid(l)*q(l)/(con_eps-con_epsm1*q(l))
155! td(l)=ftdp(pv)
156! end do
157! if(kdt>10.and.kdt<20)btim = timef()
158
159! write(0,*)' i=',i,' lm=',lm,' lp1=',lp1,' t=',t(1),q(1),pmid(1) &
160! &,' pint=',pint(1),' prec=',prec(i),' pthresh=',pthresh
161
163 call calwxt_ramer(lm,lp1,t,q,pmid,rh,td,pint,iwx)
164
165!
166 snow(2) = mod(iwx,2)
167 sleet(2) = mod(iwx,4)/2
168 freezr(2) = mod(iwx,8)/4
169 rain(2) = iwx/8
170
171! bourgouin algorithm
172! iseed=44641*(int(sdat(1)-1)*24*31+int(sdat(2))*24+ihrst)+ &
173! & mod(ifhr*60+ifmin,44641)+4357
174
176 ranl = randomno(i,1:2)
177 call calwxt_bourg(lm,lp1,ranl,con_g,t,q,pmid,pint,zint(1),iwx)
178
179!
180 snow(3) = mod(iwx,2)
181 sleet(3) = mod(iwx,4)/2
182 freezr(3) = mod(iwx,8)/4
183 rain(3) = iwx/8
184!
185! revised ncep algorithm
186!
188 call calwxt_revised(lm,lp1,t,q,pmid,pint, &
189 con_fvirt,con_rog,con_epsq,zint,twet,iwx)
190!
191 snow(4) = mod(iwx,2)
192 sleet(4) = mod(iwx,4)/2
193 freezr(4) = mod(iwx,8)/4
194 rain(4) = iwx/8
195!
196! explicit algorithm (under 18 not admitted without parent or guardian)
197
198 snow(5) = 0
199 sleet(5) = 0
200 freezr(5) = 0
201 rain(5) = 0
202!
203 call calwxt_dominant(nalg,rain(1),freezr(1),sleet(1), &
204 snow(1),domr(i),domzr(i),domip(i),doms(i))
205
206 else ! prec < pthresh
207 domr(i) = 0.
208 domzr(i) = 0.
209 domip(i) = 0.
210 doms(i) = 0.
211 end if
212 enddo ! end loop for i
213
214 deallocate (twet,rh,td)
215 return
216 end
217
221 subroutine calwxt(lm,lp1,t,q,pmid,pint, &
222 d608,rog,epsq,zint,iwx,twet)
223 use machine , only : kind_phys
224!
225! file: calwxt.f
226! written: 11 november 1993, michael baldwin
227! revisions:
228! 30 sept 1994-setup new decision tree (m baldwin)
229! 12 june 1998-conversion to 2-d (t black)
230! 01-10-25 h chuang - modified to process hybrid model output
231! 02-01-15 mike baldwin - wrf version
232!
233!
234! routine to compute precipitation type using a decision tree
235! approach that uses variables such as integrated wet bulb temp
236! below freezing and lowest layer temperature
237!
238! see baldwin and contorno preprint from 13th weather analysis
239! and forecasting conference for more details
240! (or baldwin et al, 10th nwp conference preprint)
241!
242!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
243 implicit none
244!
245! input:
246! t,q,pmid,htm,lmh,zint
247!
248 integer,intent(in) :: lm,lp1
249 real(kind=kind_phys),dimension(lm),intent(in) :: t,q,pmid,twet
250 real(kind=kind_phys),dimension(lp1),intent(in) :: zint,pint
251 integer,intent(out) :: iwx
252 real(kind=kind_phys),intent(in) :: d608,rog,epsq
253
254
255! output:
256! iwx - instantaneous weather type.
257! acts like a 4 bit binary
258! 1111 = rain/freezing rain/ice pellets/snow
259! where the one's digit is for snow
260! the two's digit is for ice pellets
261! the four's digit is for freezing rain
262! and the eight's digit is for rain
263!
264! internal:
265!
266! real(kind=kind_phys), allocatable :: twet(:)
267 real(kind=kind_phys), parameter :: d00=0.0
268 integer karr,licee
269 real(kind=kind_phys) tcold,twarm
270
271! subroutines called:
272! wetbulb
273!
274!
275! initialize weather type array to zero (ie, off).
276! we do this since we want iwx to represent the
277! instantaneous weather type on return.
278!
279!
280! allocate local storage
281!
282
283 integer l,lice,iwrml,ifrzl
284 real(kind=kind_phys) psfck,tdchk,a,tdkl,tdpre,tlmhk,twrmk,areas8,areap4, &
285 surfw,surfc,dzkl,area1,pintk1,pintk2,pm150,pkl,tkl,qkl
286
287! allocate ( twet(lm) )
288!
289 iwx = 0
290!
291! find coldest and warmest temps in saturated layer between
292! 70 mb above ground and 500 mb
293! also find highest saturated layer in that range
294!
295!meb
296 psfck = pint(lm+1)
297!meb
298 tdchk = 2.0
299 760 tcold = t(lm)
300 twarm = t(lm)
301 licee = lm
302!
303 do l=1,lm
304 qkl = q(l)
305 qkl = max(epsq,qkl)
306 tkl = t(l)
307 pkl = pmid(l)
308!
309! skip past this if the layer is not between 70 mb above ground and 500 mb
310!
311 if (pkl < 50000.0 .or. pkl > psfck-7000.0) cycle
312 a = log(qkl*pkl/(6.1078*(0.378*qkl+0.622)))
313 tdkl = (237.3*a) / (17.269-a) + 273.15
314 tdpre = tkl - tdkl
315 if (tdpre < tdchk .and. tkl < tcold) tcold = tkl
316 if (tdpre < tdchk .and. tkl > twarm) twarm = tkl
317 if (tdpre < tdchk .and. l < licee) licee = l
318 enddo
319!
320! if no sat layer at dew point dep=tdchk, increase tdchk
321! and start again (but don't make tdchk > 6)
322!
323 if (tcold == t(lm) .and. tdchk < 6.0) then
324 tdchk = tdchk + 2.0
325 goto 760
326 endif
327!
328! lowest layer t
329!
330 karr = 0
331 tlmhk = t(lm)
332!
333! decision tree time
334!
335 if (tcold > 269.15) then
336 if (tlmhk <= 273.15) then
337
338! turn on the flag for freezing rain = 4 if its not on already
339! izr=mod(iwx(i,j),8)/4
340! if (izr.lt.1) iwx(i,j)=iwx(i,j)+4
341
342 iwx = iwx + 4
343 goto 850
344 else
345! turn on the flag for rain = 8
346! if its not on already
347! irain=iwx(i,j)/8
348! if (irain.lt.1) iwx(i,j)=iwx(i,j)+8
349
350 iwx = iwx + 8
351 goto 850
352 endif
353 endif
354 karr = 1
355 850 continue
356!
357! compute wet bulb only at points that need it
358!
359! call wetbulb(lm,t,q,pmid,karr,twet)
360! call wetfrzlvl(twet,zwet)
361!
362 if (karr > 0) then
363 lice=licee
364!meb
365 psfck = pint(lm+1)
366!meb
367 tlmhk = t(lm)
368 twrmk = twarm
369!
370! twet area variables calculate only what is needed
371! from ground to 150 mb above surface from ground to tcold layer
372! and from ground to 1st layer where wet bulb t < 0.0
373!
374! pintk1 is the pressure at the bottom of the layer
375! pintk2 is the pressure at the top of the layer
376!
377! areap4 is the area of twet above -4 c below highest sat lyr
378!
379 areas8 = d00
380 areap4 = d00
381 surfw = d00
382 surfc = d00
383!
384 do l=lm,lice,-1
385 area1 = (twet(l)-269.15) * (zint(l)-zint(l+1))
386 if (twet(l) >= 269.15) areap4 = areap4 + area1
387 enddo
388!
389 if (areap4 < 3000.0) then
390! turn on the flag for snow = 1
391! if its not on already
392! isno=mod(iwx(i,j),2)
393! if (isno.lt.1) iwx(i,j)=iwx(i,j)+1
394
395 iwx = iwx + 1
396 return
397 endif
398!
399! areas8 is the net area of twet w.r.t. freezing in lowest 150mb
400!
401 pintk1 = psfck
402 pm150 = psfck - 15000.
403!
404 do l=lm,1,-1
405 pintk2 = pint(l)
406 if (pintk1 >= pm150) then
407 dzkl = zint(l)-zint(l+1)
408! sum partial layer if in 150 mb agl layer
409 if (pintk2 < pm150) &
410 dzkl = t(l)*(q(l)*d608+1.0)*rog*log(pintk1/pm150)
411 area1 = (twet(l)-273.15)*dzkl
412 areas8 = areas8 + area1
413 endif
414 pintk1 = pintk2
415 enddo
416!
417! surfw is the area of twet above freezing between the ground
418! and the first layer above ground below freezing
419! surfc is the area of twet below freezing between the ground
420! and the warmest sat layer
421!
422 ifrzl = 0
423 iwrml = 0
424!
425 do l=lm,1,-1
426 if (ifrzl == 0 .and. t(l) < 273.15) ifrzl = 1
427 if (iwrml == 0 .and. t(l) >= twrmk) iwrml = 1
428!
429 if (iwrml == 0 .or. ifrzl == 0) then
430! if(pmid(l) < 50000.)print*,'need twet above 500mb'
431 dzkl = zint(l)-zint(l+1)
432 area1 = (twet(l)-273.15)*dzkl
433 if(ifrzl == 0 .and. twet(l) >= 273.15) surfw = surfw + area1
434 if(iwrml == 0 .and. twet(l) <= 273.15) surfc = surfc + area1
435 endif
436 enddo
437 if(surfc < -3000.0 .or. (areas8 < -3000.0 .and. surfw < 50.0)) then
438! turn on the flag for ice pellets = 2 if its not on already
439! iip=mod(iwx(i,j),4)/2
440! if (iip.lt.1) iwx(i,j)=iwx(i,j)+2
441 iwx = iwx + 2
442!
443 elseif(tlmhk < 273.15) then
444! turn on the flag for freezing rain = 4 if its not on already
445! izr=mod(iwx(k),8)/4
446! if (izr.lt.1) iwx(k)=iwx(k)+4
447 iwx = iwx + 4
448 else
449! turn on the flag for rain = 8 if its not on already
450! irain=iwx(k)/8
451! if (irain.lt.1) iwx(k)=iwx(k)+8
452 iwx = iwx + 8
453 endif
454 endif
455!---------------------------------------------------------
456! deallocate (twet)
457
458 return
459 end
460!
461!
462!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
463!
464! dophase is a subroutine written and provided by jim ramer at noaa/fsl
465!
466! ramer, j, 1993: an empirical technique for diagnosing precipitation
467! type from model output. preprints, 5th conf. on aviation
468! weather systems, vienna, va, amer. meteor. soc., 227-230.
469!
470! code adapted for wrf post 24 august 2005 g manikin
471!ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
474 subroutine calwxt_ramer(lm,lp1,t,q,pmid,rh,td,pint,ptyp)
475
476! subroutine dophase(pq, ! input pressure sounding mb
477! + t, ! input temperature sounding k
478! + pmid, ! input pressure
479! + pint, ! input interface pressure
480! + q, ! input spec humidityfraction
481! + lmh, ! input number of levels in sounding
482! + ptyp) ! output(2) phase 2=rain, 3=frzg, 4=solid,
483! 6=ip jc 9/16/99
484! use params_mod
485! use ctlblk_mod
486!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
487 use machine , only : kind_phys
488 implicit none
489!
490 real(kind=kind_phys),parameter :: twice=266.55,rhprcp=0.80,deltag=1.02, &
491 & emelt=0.045,rlim=0.04,slim=0.85
492 real(kind=kind_phys),parameter :: twmelt=273.15,tz=273.15,efac=1.0 ! specify in params now
493!
494 integer*4 i, k1, lll, k2, toodry
495!
496 real(kind=kind_phys) xxx ,mye, icefrac
497 integer, intent(in) :: lm,lp1
498 real(kind=kind_phys),dimension(lm), intent(in) :: t,q,pmid,rh,td
499 real(kind=kind_phys),dimension(lp1),intent(in) :: pint
500 integer, intent(out) :: ptyp
501!
502 real(kind=kind_phys),dimension(lm) :: tq,pq,rhq,twq
503!
504 integer j,l,lev,ii
505 real(kind=kind_phys) rhmax,twmax,ptop,dpdrh,twtop,rhtop,wgt1,wgt2, &
506 rhavg,dtavg,dpk,ptw,pbot
507! real(kind=kind_phys) b,qtmp,rate,qc
508!
509!
510! initialize.
511 icefrac = -9999.
512!
513
514 ptyp = 0
515 do l = 1,lm
516 lev = lp1 - l
517! p(l)=pmid(l)
518! qc=pq0/p(l) * exp(a2*(t(l)-a3)/(t(l)-a4))
519!gsm forcing q (qtmp) to be positive to deal with negative q values
520! causing problems later in this subroutine
521! qtmp=max(h1m12,q(l))
522! rhqtmp(lev)=qtmp/qc
523 rhq(lev) = rh(l)
524 pq(lev) = pmid(l) * 0.01
525 tq(lev) = t(l)
526 enddo
527
528
529!
530!cc rate restriction removed by john cortinas 3/16/99
531!
532! construct wet-bulb sounding, locate generating level.
533 twmax = -999.0
534 rhmax = 0.0
535 k1 = 0 ! top of precip generating layer
536 k2 = 0 ! layer of maximum rh
537!
538 if (rhq(1) < rhprcp) then
539 toodry = 1
540 else
541 toodry = 0
542 end if
543!
544 pbot = pq(1)
545! nq=lm
546 do l = 1, lm
547! xxx = tdofesat(esat(tq(l))*rhq(l))
548! xxx = td(l) !hc: use td consistent with gfs ice physics
549 xxx = td(lp1-l) !hc: use td consistent with gfs ice physics
550 if (xxx < -500.) return
551 twq(l) = xmytw(tq(l),xxx,pq(l))
552 twmax = max(twq(l),twmax)
553 if (pq(l) >= 400.0) then
554 if (rhq(l) > rhmax) then
555 rhmax = rhq(l)
556 k2 = l
557 end if
558!
559 if (l /= 1) then
560 if (rhq(l) >= rhprcp .or. toodry == 0) then
561 if (toodry /= 0) then
562 dpdrh = log(pq(l)/pq(l-1)) / (rhq(l)-rhq(l-1))
563 pbot = exp(log(pq(l))+(rhprcp-rhq(l))*dpdrh)
564!
565 ptw = pq(l)
566 toodry = 0
567 else if (rhq(l)>= rhprcp) then
568 ptw = pq(l)
569 else
570 toodry = 1
571 dpdrh = log(pq(l)/pq(l-1)) / (rhq(l)-rhq(l-1))
572 ptw = exp(log(pq(l))+(rhprcp-rhq(l))*dpdrh)
573
574!lin dpdrh = (pq(i)-pq(i-1))/(rhq(i)-rhq(i-1))
575!lin ptw = pq(i)+(rhprcp-rhq(i))*dpdrh
576!
577 end if
578!
579 if (pbot/ptw >= deltag) then
580!lin if (pbot-ptw.lt.deltag) goto 2003
581 k1 = l
582 ptop = ptw
583 end if
584 end if
585 end if
586 end if
587 enddo
588!
589! gross checks for liquid and solid precip which dont require generating level.
590!
591 if (twq(1) >= 273.15+2.0) then
592 ptyp = 8 ! liquid
593 icefrac = 0.0
594 return
595 end if
596!
597 if (twmax <= twice) then
598 icefrac = 1.0
599 ptyp = 1 ! solid
600 return
601 end if
602!
603! check to see if we had no success with locating a generating level.
604!
605 if (k1 == 0) return
606!
607 if (ptop == pq(k1)) then
608 twtop = twq(k1)
609 rhtop = rhq(k1)
610 k2 = k1
611 k1 = k1 - 1
612 else
613 k2 = k1
614 k1 = k1 - 1
615 wgt1 = log(ptop/pq(k2)) / log(pq(k1)/pq(k2))
616 wgt2 = 1.0 - wgt1
617 twtop = twq(k1) * wgt1 + twq(k2) * wgt2
618 rhtop = rhq(k1) * wgt1 + rhq(k2) * wgt2
619 end if
620!
621! calculate temp and wet-bulb ranges below precip generating level.
622 do l = 1, k1
623 twmax = max(twq(l),twmax)
624 enddo
625!
626! gross check for solid precip, initialize ice fraction.
627! if (i.eq.1.and.j.eq.1) write (*,*) 'twmax=',twmax,twice,'twtop=',twtop
628
629 if (twtop <= twice) then
630 icefrac = 1.0
631 if (twmax <= twmelt) then ! gross check for solid precip.
632 ptyp = 1 ! solid precip
633 return
634 end if
635 lll = 0
636 else
637 icefrac = 0.0
638 lll = 1
639 end if
640!
641! loop downward through sounding from highest precip generating level.
642 30 continue
643!
644 if (icefrac >= 1.0) then ! starting as all ice
645 if (twq(k1) < twmelt) go to 40 ! cannot commence melting
646 if (twq(k1) == twtop) go to 40 ! both equal twmelt, nothing h
647 wgt1 = (twmelt-twq(k1)) / (twtop-twq(k1))
648 rhavg = rhq(k1) + wgt1 * (rhtop-rhq(k1)) * 0.5
649 dtavg = (twmelt-twq(k1)) * 0.5
650 dpk = wgt1 * log(pq(k1)/ptop) !lin dpk=wgt1*(pq(k1)-ptop)
651! mye=emelt*(1.0-(1.0-rhavg)*efac)
652 mye = emelt * rhavg ** efac
653 icefrac = icefrac + dpk * dtavg / mye
654 else if (icefrac <= 0.0) then ! starting as all liquid
655 lll = 1
656! goto 1020
657 if (twq(k1) > twice) go to 40 ! cannot commence freezing
658 if (twq(k1) == twtop) then
659 wgt1 = 0.5
660 else
661 wgt1 = (twice-twq(k1)) / (twtop-twq(k1))
662 end if
663 rhavg = rhq(k1) + wgt1 * (rhtop-rhq(k1)) * 0.5
664 dtavg = twmelt - (twq(k1)+twice) * 0.5
665 dpk = wgt1 * log(pq(k1)/ptop) !lin dpk=wgt1*(pq(k1)-ptop)
666! mye = emelt*(1.0-(1.0-rhavg)*efac)
667 mye = emelt * rhavg ** efac
668 icefrac = icefrac + dpk * dtavg / mye
669 else if ((twq(k1) <= twmelt).and.(twq(k1) < twmelt)) then ! mix
670 rhavg = (rhq(k1)+rhtop) * 0.5
671 dtavg = twmelt - (twq(k1)+twtop) * 0.5
672 dpk = log(pq(k1)/ptop) !lin dpk=pq(k1)-ptop
673! mye = emelt*(1.0-(1.0-rhavg)*efac)
674 mye = emelt * rhavg ** efac
675 icefrac = icefrac + dpk * dtavg / mye
676 else ! mix where tw curve crosses twmelt in layer
677 if (twq(k1) == twtop) go to 40 ! both equal twmelt, nothing h
678 wgt1 = (twmelt-twq(k1)) / (twtop-twq(k1))
679 wgt2 = 1.0 - wgt1
680 rhavg = rhtop + wgt2 * (rhq(k1)-rhtop) * 0.5
681 dtavg = (twmelt-twtop) * 0.5
682 dpk = wgt2 * log(pq(k1)/ptop) !lin dpk=wgt2*(pq(k1)-ptop)
683! mye = emelt*(1.0-(1.0-rhavg)*efac)
684 mye = emelt * rhavg ** efac
685 icefrac = icefrac + dpk * dtavg / mye
686 icefrac = min(1.0,max(icefrac,0.0))
687 if (icefrac <= 0.0) then
688! goto 1020
689 if (twq(k1) > twice) go to 40 ! cannot commence freezin
690 wgt1 = (twice-twq(k1)) / (twtop-twq(k1))
691 dtavg = twmelt - (twq(k1)+twice) * 0.5
692 else
693 dtavg = (twmelt-twq(k1)) * 0.5
694 end if
695 rhavg = rhq(k1) + wgt1 * (rhtop-rhq(k1)) * 0.5
696 dpk = wgt1 * log(pq(k1)/ptop) !lin dpk=wgt1*(pq(k1)-ptop)
697! mye = emelt*(1.0-(1.0-rhavg)*efac)
698 mye = emelt * rhavg ** efac
699 icefrac = icefrac + dpk * dtavg / mye
700 end if
701!
702 icefrac = min(1.0,max(icefrac,0.0))
703
704! if (i.eq.1.and.j.eq.1) write (*,*) 'new icefrac:', icefrac, icefrac
705!
706! get next level down if there is one, loop back.
707 40 continue
708 if (k1 > 1) then
709 twtop = twq(k1)
710 ptop = pq(k1)
711 rhtop = rhq(k1)
712 k1 = k1 - 1
713 go to 30
714 end if
715!
716! determine precip type based on snow fraction and surface wet-bulb.
717!
718 if (icefrac >= slim) then
719 if (lll /= 0) then
720 ptyp = 2 ! ice pellets jc 9/16/99
721 else
722 ptyp = 1 ! snow
723 end if
724 else if (icefrac <= rlim) then
725 if (twq(1).lt.tz) then
726 ptyp = 4 ! freezing precip
727 else
728 ptyp = 8 ! rain
729 end if
730 else
731 if (twq(1) < tz) then
732!gsm not sure what to do when 'mix' is predicted; in previous
733!gsm versions of this code for which i had to have an answer,
734!gsm i chose sleet. here, though, since we have 4 other
735!gsm algorithms to provide an answer, i will not declare a
736!gsm type from the ramer in this situation and allow the
737!gsm other algorithms to make the call.
738
739 ptyp = 0 ! don't know
740! ptyp = 5 ! mix
741 else
742! ptyp = 5 ! mix
743 ptyp = 0 ! don't know
744 end if
745 end if
746
747 return
748!
749 end
750!
751!
752!--------------------------------------------------------------------------
753 function xmytw(t,td,p)
754!
755 use machine , only : kind_phys
756 implicit none
757!
758 integer*4 cflag, l
759 real(kind=kind_phys) f, c0, c1, c2, k, kd, kw, ew, t, td, p, ed, fp, s, &
760 & de, xmytw
761 data f, c0, c1, c2 /0.0006355, 26.66082, 0.0091379024, 6106.3960/
762!
763 xmytw = (t+td) / 2
764 if (td >= t) return
765!
766 if (t < 100.0) then
767 k = t + 273.15
768 kd = td + 273.15
769 if (kd >= k) return
770 cflag = 1
771 else
772 k = t
773 kd = td
774 cflag = 0
775 end if
776!
777 ed = c0 - c1 * kd - c2 / kd
778 if (ed < -14.0 .or. ed > 7.0) return
779 ed = exp(ed)
780 ew = c0 - c1 * k - c2 / k
781 if (ew < -14.0 .or. ew > 7.0) return
782 ew = exp(ew)
783 fp = p * f
784 s = (ew-ed) / (k-kd)
785 kw = (k*fp+kd*s) / (fp+s)
786!
787 do l = 1, 5
788 ew = c0 - c1 * kw - c2 / kw
789 if (ew < -14.0 .or. ew > 7.0) return
790 ew = exp(ew)
791 de = fp * (k-kw) + ed - ew
792 if (abs(de/ew) < 1e-5) exit
793 s = ew * (c1-c2/(kw*kw)) - fp
794 kw = kw - de / s
795 enddo
796!
797! print *, 'kw ', kw
798 if (cflag /= 0) then
799 xmytw = kw - 273.15
800 else
801 xmytw = kw
802 end if
803!
804 return
805 end
806!
807!
808!$$$ subprogram documentation block
809!
810! subprogram: calwxt_bourg calculate precipitation type (bourgouin)
811! prgmmr: baldwin org: np22 date: 1999-07-06
812!
813! abstract: this routine computes precipitation type
814! using a decision tree approach that uses the so-called
815! "energy method" of bourgouin of aes (canada) 1992
816!
817! program history log:
818! 1999-07-06 m baldwin
819! 1999-09-20 m baldwin make more consistent with bourgouin (1992)
820! 2005-08-24 g manikin added to wrf post
821! 2007-06-19 m iredell mersenne twister, best practices
822! 2008-03-03 g manikin added checks to prevent stratospheric warming
823! episodes from being seen as "warm" layers
824! impacting precip type
825!
826! usage: call calwxt_bourg(im,jm,jsta_2l,jend_2u,jsta,jend,lm,lp1, &
827! & iseed,g, &
828! & t,q,pmid,pint,lmh,zint,ptype)
829! input argument list:
830! im integer i dimension
831! jm integer j dimension
832! jsta_2l integer j dimension start point (including haloes)
833! jend_2u integer j dimension end point (including haloes)
834! jsta integer j dimension start point (excluding haloes)
835! jend integer j dimension end point (excluding haloes)
836! lm integer k dimension
837! lp1 integer k dimension plus 1
838! iseed integer random number seed
839! g real gravity (m/s**2)
840! pthresh real precipitation threshold (m)
841! t real(im,jsta_2l:jend_2u,lm) mid layer temp (k)
842! q real(im,jsta_2l:jend_2u,lm) specific humidity (kg/kg)
843! pmid real(im,jsta_2l:jend_2u,lm) mid layer pressure (pa)
844! pint real(im,jsta_2l:jend_2u,lp1) interface pressure (pa)
845! lmh real(im,jsta_2l:jend_2u) max number of layers
846! zint real(im,jsta_2l:jend_2u,lp1) interface height (m)
847! output argument list:
848! ptype real(im,jm) instantaneous weather type ()
849! acts like a 4 bit binary
850! 1111 = rain/freezing rain/ice pellets/snow
851! where the one's digit is for snow
852! the two's digit is for ice pellets
853! the four's digit is for freezing rain
854! and the eight's digit is for rain
855! in other words...
856! ptype=1 snow
857! ptype=2 ice pellets/mix with ice pellets
858! ptype=4 freezing rain/mix with freezing rain
859! ptype=8 rain
860!
861! modules used:
862! mersenne_twister pseudo-random number generator
863!
864! subprograms called:
865! random_number pseudo-random number generator
866!
867! attributes:
868! language: fortran 90
869!
870! remarks: vertical order of arrays must be layer 1 = top
871! and layer lmh = bottom
872!
873!$$$
877!of aes (canada) 1992
878 subroutine calwxt_bourg(lm,lp1,rn,g,t,q,pmid,pint,zint,ptype)
879 use machine , only : kind_phys
880 implicit none
881!
882! input:
883 integer,intent(in) :: lm,lp1
884 real(kind=kind_phys),intent(in) :: g,rn(2)
885 real(kind=kind_phys),intent(in), dimension(lm) :: t, q, pmid
886 real(kind=kind_phys),intent(in), dimension(lp1) :: pint, zint
887!
888! output:
889 integer, intent(out) :: ptype
890!
891 integer ifrzl,iwrml,l,lhiwrm
892 real(kind=kind_phys) pintk1,areane,tlmhk,areape,pintk2,surfw,area1,dzkl,psfck,r1,r2
893!
894! initialize weather type array to zero (ie, off).
895! we do this since we want ptype to represent the
896! instantaneous weather type on return.
897!
898 ptype = 0
899 psfck = pint(lm+1)
900
901! find the depth of the warm layer based at the surface
902! this will be the cut off point between computing
903! the surface based warm air and the warm air aloft
904!
905! lowest layer t
906!
907 tlmhk = t(lm)
908 iwrml = lm + 1
909 if (tlmhk >= 273.15) then
910 do l = lm, 2, -1
911 if (t(l) >= 273.15 .and. t(l-1) < 273.15 .and. &
912 & iwrml == lm+1) iwrml = l
913 end do
914 end if
915!
916! now find the highest above freezing level
917!
918 lhiwrm = lm + 1
919 do l = lm, 1, -1
920! gsm added 250 mb check to prevent stratospheric warming situations
921! from counting as warm layers aloft
922 if (t(l) >= 273.15 .and. pmid(l) > 25000.) lhiwrm = l
923 end do
924
925! energy variables
926! surfw is the positive energy between the ground
927! and the first sub-freezing layer above ground
928! areane is the negative energy between the ground
929! and the highest layer above ground
930! that is above freezing
931! areape is the positive energy "aloft"
932! which is the warm energy not based at the ground
933! (the total warm energy = surfw + areape)
934!
935! pintk1 is the pressure at the bottom of the layer
936! pintk2 is the pressure at the top of the layer
937! dzkl is the thickness of the layer
938! ifrzl is a flag that tells us if we have hit
939! a below freezing layer
940!
941 pintk1 = psfck
942 ifrzl = 0
943 areane = 0.0
944 areape = 0.0
945 surfw = 0.0
946
947 do l = lm, 1, -1
948 if (ifrzl == 0 .and. t(l) <= 273.15) ifrzl = 1
949 pintk2 = pint(l)
950 dzkl = zint(l)-zint(l+1)
951 if (t(l) >= 273.15 .and. pmid(l) > 25000.) then
952 area1 = log(t(l)/273.15) * g * dzkl
953 if (l < iwrml) then
954 areape = areape + area1
955 else
956 surfw = surfw + area1
957 endif
958 elseif (l > lhiwrm) then
959 area1 = log(t(l)/273.15) * g * dzkl
960 areane = areane + abs(area1)
961 endif
962 pintk1 = pintk2
963 enddo
964
965!
966! decision tree time
967!
968 if (areape < 2.0) then ! very little or no positive energy aloft, check for
969 ! positive energy just above the surface to determine rain vs. snow
970 if (surfw < 5.6) then ! not enough positive energy just above the surface snow = 1
971 ptype = 1
972 else if (surfw > 13.2) then ! enough positive energy just above the surface rain = 8
973 ptype = 8
974 else ! transition zone, assume equally likely rain/snow
975 ! picking a random number, if <=0.5 snow
976 r1 = rn(1)
977 if (r1 <= 0.5) then ! snow = 1
978 ptype = 1
979 else ! rain = 8
980 ptype = 8
981 end if
982 end if
983!
984 else ! some positive energy aloft, check for enough negative energy
985 ! to freeze and make ice pellets to determine ip vs. zr
986
987 if (areane > 66.0+0.66*areape) then
988! enough negative area to make ip,
989! now need to check if there is enough positive energy
990! just above the surface to melt ip to make rain
991 if (surfw < 5.6) then ! not enough energy at the surface to melt ip ice pellets = 2
992 ptype = 2
993 elseif (surfw > 13.2) then ! enough energy at the surface to melt ip rain = 8
994 ptype = 8
995 else ! transition zone, assume equally likely ip/rain picking a random number, if <=0.5 ip
996 r1 = rn(1)
997 if (r1 <= 0.5) then ! ice pellets = 2
998 ptype = 2
999 else ! rain = 8
1000 ptype = 8
1001 end if
1002 end if
1003 elseif (areane < 46.0+0.66*areape) then
1004! not enough negative energy to refreeze, check surface temp to determine rain vs. zr
1005 if (tlmhk < 273.15) then ! freezing rain = 4
1006 ptype = 4
1007 else ! rain = 8
1008 ptype = 8
1009 end if
1010 else
1011! transition zone, assume equally likely ip/zr picking a random number, if <=0.5 ip
1012 r1 = rn(1)
1013 if (r1 <= 0.5) then
1014! still need to check positive energy just above the surface to melt ip vs. rain
1015 if (surfw < 5.6) then ! ice pellets = 2
1016 ptype = 2
1017 else if (surfw > 13.2) then ! rain = 8
1018 ptype = 8
1019 else
1020! transition zone, assume equally likely ip/rain picking a random number, if <=0.5 ip
1021 r2 = rn(2)
1022 if (r2 <= 0.5) then ! ice pellets = 2
1023 ptype = 2
1024 else ! rain = 8
1025 ptype = 8
1026 end if
1027 end if
1028 else
1029! not enough negative energy to refreeze, check surface temp to determine rain vs. zr
1030 if (tlmhk < 273.15) then ! freezing rain = 4
1031 ptype = 4
1032 else ! rain = 8
1033 ptype = 8
1034 end if
1035 end if
1036 end if
1037 end if
1038!
1039 return
1040 end
1041!
1048 subroutine calwxt_revised(lm,lp1,t,q,pmid,pint, &
1049 d608,rog,epsq,zint,twet,iwx)
1050!
1051! file: calwxt.f
1052! written: 11 november 1993, michael baldwin
1053! revisions:
1054! 30 sept 1994-setup new decision tree (m baldwin)
1055! 12 june 1998-conversion to 2-d (t black)
1056! 01-10-25 h chuang - modified to process hybrid model output
1057! 02-01-15 mike baldwin - wrf version
1058! 05-07-07 binbin zhou - add prec for rsm
1059! 05-08-24 geoff manikin - modified the area requirements
1060! to make an alternate algorithm
1061!
1062!
1063! routine to compute precipitation type using a decision tree
1064! approach that uses variables such as integrated wet bulb temp
1065! below freezing and lowest layer temperature
1066!
1067! see baldwin and contorno preprint from 13th weather analysis
1068! and forecasting conference for more details
1069! (or baldwin et al, 10th nwp conference preprint)
1070!
1071! since the original version of the algorithm has a high bias
1072! for freezing rain and sleet, the goal is to balance that bias
1073! with a version more likely to predict snow
1074!
1075! use params_mod
1076! use ctlblk_mod
1077!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1078 use machine , only : kind_phys
1079 implicit none
1080!
1081! list of variables needed
1082! parameters:
1083! d608,rog,h1,d00
1084!hc parameter(d608=0.608,rog=287.04/9.8,h1=1.0,d00=0.0)
1085!
1086! input:
1087! t,q,pmid,htm,lmh,zint
1088
1089 integer,intent(in) :: lm,lp1
1090 real(kind=kind_phys),dimension(lm),intent(in) :: t,q,pmid,twet
1091 real(kind=kind_phys),dimension(lp1),intent(in) :: pint,zint
1092 real(kind=kind_phys),intent(in) :: d608,rog,epsq
1093! output:
1094! iwx - instantaneous weather type.
1095! acts like a 4 bit binary
1096! 1111 = rain/freezing rain/ice pellets/snow
1097! where the one's digit is for snow
1098! the two's digit is for ice pellets
1099! the four's digit is for freezing rain
1100! and the eight's digit is for rain
1101 integer, intent(out) :: iwx
1102! internal:
1103!
1104 real(kind=kind_phys), parameter :: d00=0.0
1105 integer karr,licee
1106 real(kind=kind_phys) tcold,twarm
1107!
1108 integer l,lmhk,lice,iwrml,ifrzl
1109 real(kind=kind_phys) psfck,tdchk,a,tdkl,tdpre,tlmhk,twrmk,areas8,areap4,area1, &
1110 surfw,surfc,dzkl,pintk1,pintk2,pm150,qkl,tkl,pkl,area0,areap0
1111
1112! subroutines called:
1113! wetbulb
1114!
1115!
1116! initialize weather type array to zero (ie, off).
1117! we do this since we want iwx to represent the
1118! instantaneous weather type on return.
1119!
1120!
1121! allocate local storage
1122!
1123!
1124 iwx = 0
1125 lmhk=lm
1126!
1127! find coldest and warmest temps in saturated layer between
1128! 70 mb above ground and 500 mb
1129! also find highest saturated layer in that range
1130!
1131!meb
1132 psfck = pint(lp1)
1133!meb
1134 tdchk = 2.0
1135 760 tcold = t(lmhk)
1136 twarm = t(lmhk)
1137 licee = lmhk
1138!
1139 do l=1,lmhk
1140 qkl = q(l)
1141 qkl = max(epsq,qkl)
1142 tkl = t(l)
1143 pkl = pmid(l)
1144!
1145! skip past this if the layer is not between 70 mb above ground
1146! and 500 mb
1147!
1148 if (pkl < 50000.0 .or. pkl > psfck-7000.0) cycle
1149 a = log(qkl*pkl/(6.1078*(0.378*qkl+0.622)))
1150 tdkl = (237.3*a)/(17.269-a)+273.15
1151 tdpre = tkl-tdkl
1152 if (tdpre < tdchk .and. tkl < tcold) tcold = tkl
1153 if (tdpre < tdchk .and. tkl > twarm) twarm = tkl
1154 if (tdpre < tdchk .and. l < licee) licee = l
1155 enddo
1156!
1157! if no sat layer at dew point dep=tdchk, increase tdchk
1158! and start again (but don't make tdchk > 6)
1159!
1160 if (tcold == t(lmhk) .and. tdchk < 6.0) then
1161 tdchk = tdchk + 2.0
1162 goto 760
1163 endif
1164!
1165! lowest layer t
1166!
1167 karr = 0
1168 lmhk = lm
1169 tlmhk = t(lmhk)
1170!
1171! decision tree time
1172!
1173 if (tcold > 269.15) then
1174 if (tlmhk <= 273.15) then
1175! turn on the flag for freezing rain = 4 if its not on already
1176! izr=mod(iwx,8)/4
1177! if (izr.lt.1) iwx=iwx+4
1178 iwx = iwx + 4
1179 goto 850
1180 else
1181! turn on the flag for rain = 8 if its not on already
1182! irain=iwx/8
1183! if (irain.lt.1) iwx=iwx+8
1184 iwx = iwx + 8
1185 goto 850
1186 endif
1187 endif
1188 karr = 1
1189 850 continue
1190!
1191 if (karr > 0)then
1192 lmhk = lm
1193 lice = licee
1194!meb
1195 psfck = pint(lp1)
1196!meb
1197 tlmhk = t(lmhk)
1198 twrmk = twarm
1199!
1200! twet area variables
1201! calculate only what is needed
1202! from ground to 150 mb above surface
1203! from ground to tcold layer
1204! and from ground to 1st layer where wet bulb t < 0.0
1205!
1206! pintk1 is the pressure at the bottom of the layer
1207! pintk2 is the pressure at the top of the layer
1208!
1209! areap4 is the area of twet above -4 c below highest sat lyr
1210! areap0 is the area of twet above 0 c below highest sat lyr
1211!
1212 areas8 = d00
1213 areap4 = d00
1214 areap0 = d00
1215 surfw = d00
1216 surfc = d00
1217
1218!
1219 do l=lmhk,lice,-1
1220 dzkl = zint(l)-zint(l+1)
1221 area1 = (twet(l)-269.15)*dzkl
1222 area0 = (twet(l)-273.15)*dzkl
1223 if (twet(l) >= 269.15) areap4 = areap4 + area1
1224 if (twet(l) >= 273.15) areap0 = areap0 + area0
1225 enddo
1226!
1227! if (areap4.lt.3000.0) then turn on the flag for snow = 1 if its not on already
1228! isno=mod(iwx,2)
1229! if (isno.lt.1) iwx=iwx+1
1230! iwx=iwx+1
1231! go to 1900
1232! endif
1233 if (areap0 < 350.0) then ! turn on the flag for snow = 1
1234 iwx = iwx + 1
1235 return
1236 endif
1237!
1238! areas8 is the net area of twet w.r.t. freezing in lowest 150mb
1239!
1240 pintk1 = psfck
1241 pm150 = psfck - 15000.
1242!
1243 do l=lmhk,1,-1
1244 pintk2 = pint(l)
1245 if(pintk1 >= pm150) then
1246 dzkl = zint(l)-zint(l+1)
1247!
1248! sum partial layer if in 150 mb agl layer
1249!
1250 if(pintk2 < pm150) dzkl = t(l)*(q(l)*d608+1.0)*rog* &
1251 log(pintk1/pm150)
1252 area1 = (twet(l)-273.15)*dzkl
1253 areas8 = areas8 + area1
1254 endif
1255 pintk1=pintk2
1256 enddo
1257!
1258! surfw is the area of twet above freezing between the ground
1259! and the first layer above ground below freezing
1260! surfc is the area of twet below freezing between the ground
1261! and the warmest sat layer
1262!
1263 ifrzl = 0
1264 iwrml = 0
1265!
1266 do l=lmhk,1,-1
1267 if (ifrzl == 0 .and. t(l) < 273.15) ifrzl = 1
1268 if (iwrml == 0 .and. t(l) >= twrmk) iwrml = 1
1269!
1270 if (iwrml == 0 .or. ifrzl == 0) then
1271! if(pmid(l) .lt. 50000.)print*,'twet needed above 500mb'
1272 dzkl = zint(l)-zint(l+1)
1273 area1 = (twet(l)-273.15)*dzkl
1274 if(ifrzl == 0 .and. twet(l) >= 273.15) surfw = surfw + area1
1275 if(iwrml == 0 .and. twet(l) <= 273.15) surfc = surfc + area1
1276 endif
1277 enddo
1278 if (surfc < -3000.0 .or. &
1279 & (areas8 < -3000.0 .and. surfw < 50.0)) then
1280! turn on the flag for ice pellets = 2 if its not on already
1281! iip=mod(iwx,4)/2
1282! if (iip.lt.1) iwx=iwx+2
1283 iwx = iwx + 2
1284 return
1285 endif
1286!
1287 if (tlmhk < 273.15) then
1288! turn on the flag for freezing rain = 4 if its not on already
1289! izr=mod(iwx(k),8)/4
1290! if (izr.lt.1) iwx(k)=iwx(k)+4
1291 iwx = iwx + 4
1292 else
1293! turn on the flag for rain = 8 if its not on already
1294! irain=iwx(k)/8
1295! if (irain.lt.1) iwx(k)=iwx(k)+8
1296 iwx = iwx + 8
1297 endif
1298! print *, 'revised check ', iwx(500,800)
1299 endif
1300
1301 return
1302 end
1303!
1308 subroutine calwxt_dominant(nalg,rain,freezr,sleet,snow, &
1309 & domr,domzr,domip,doms)
1310!
1311! written: 24 august 2005, g manikin
1312!
1313! this routine takes the precip type solutions from different
1314! algorithms and sums them up to give a dominant type
1315!
1316!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1317 use machine , only : kind_phys
1318 implicit none
1319!
1320! input:
1321 integer,intent(in) :: nalg
1322 real(kind=kind_phys),intent(out) :: doms,domr,domzr,domip
1323 integer,dimension(nalg),intent(in) :: rain,snow,sleet,freezr
1324 integer l
1325 real(kind=kind_phys) totsn,totip,totr,totzr
1326!--------------------------------------------------------------------------
1327! print* , 'into dominant'
1328 domr = 0.
1329 doms = 0.
1330 domzr = 0.
1331 domip = 0.
1332!
1333 totsn = 0
1334 totip = 0
1335 totr = 0
1336 totzr = 0
1338 do l = 1, nalg
1339 if (rain(l) > 0) then
1340 totr = totr + 1
1341 elseif (snow(l) > 0) then
1342 totsn = totsn + 1
1343 elseif (sleet(l) > 0) then
1344 totip = totip + 1
1345 elseif (freezr(l) > 0) then
1346 totzr = totzr + 1
1347 endif
1348 enddo
1349
1352 if (totsn > totip) then
1353 if (totsn > totzr) then
1354 if (totsn >= totr) then
1355 doms = 1
1356 else
1357 domr = 1
1358 endif
1359 elseif (totzr >= totr) then
1360 domzr = 1
1361 else
1362 domr = 1
1363 endif
1364 else if (totip > totzr) then
1365 if (totip >= totr) then
1366 domip = 1
1367 else
1368 domr = 1
1369 endif
1370 else if (totzr >= totr) then
1371 domzr = 1
1372 else
1373 domr = 1
1374 endif
1375!
1376 return
1377 end
1378end module calpreciptype_mod
This module defines four algorithms that are called to calculate dominant precipitation type,...