UPP (develop)
Loading...
Searching...
No Matches
CALWXT_REVISED.f
1 SUBROUTINE calwxt_revised_post(T,Q,PMID,PINT,HTM,LMH,PREC,ZINT,IWX)
2!
3! FILE: CALWXT.f
4! WRITTEN: 11 NOVEMBER 1993, MICHAEL BALDWIN
5! REVISIONS:
6! 30 SEPT 1994-SETUP NEW DECISION TREE (M BALDWIN)
7! 12 JUNE 1998-CONVERSION TO 2-D (T BLACK)
8! 01-10-25 H CHUANG - MODIFIED TO PROCESS HYBRID MODEL OUTPUT
9! 02-01-15 MIKE BALDWIN - WRF VERSION
10! 05-07-07 BINBIN ZHOU - ADD PREC FOR RSM
11! 05-08-24 GEOFF MANIKIN - MODIFIED THE AREA REQUIREMENTS
12! TO MAKE AN ALTERNATE ALGORITHM
13! 19-10-30 Bo CUI - REMOVE "GOTO" STATEMENT
14! 21-10-31 JESSE MENG - 2D DECOMPOSITION
15!
16!
17! ROUTINE TO COMPUTE PRECIPITATION TYPE USING A DECISION TREE
18! APPROACH THAT USES VARIABLES SUCH AS INTEGRATED WET BULB TEMP
19! BELOW FREEZING AND LOWEST LAYER TEMPERATURE
20!
21! SEE BALDWIN AND CONTORNO PREPRINT FROM 13TH WEATHER ANALYSIS
22! AND FORECASTING CONFERENCE FOR MORE DETAILS
23! (OR BALDWIN ET AL, 10TH NWP CONFERENCE PREPRINT)
24!
25! SINCE THE ORIGINAL VERSION OF THE ALGORITHM HAS A HIGH BIAS
26! FOR FREEZING RAIN AND SLEET, THE GOAL IS TO BALANCE THAT BIAS
27! WITH A VERSION MORE LIKELY TO PREDICT SNOW
28!
29 use params_mod, only: h1m12, d00, d608, h1, rog
30 use ctlblk_mod, only: jsta, jend, modelname, pthresh, im, jsta_2l, jend_2u, lm,&
31 lp1, spval, ista, iend, ista_2l, iend_2u
32!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
33 implicit none
34!
35! LIST OF VARIABLES NEEDED
36! PARAMETERS:
37! D608,ROG,H1,D00
38!HC PARAMETER(D608=0.608,ROG=287.04/9.8,H1=1.0,D00=0.0)
39!
40! INPUT:
41! T,Q,PMID,HTM,LMH,PREC,ZINT
42 REAL,dimension(ista_2l:iend_2u,jsta_2l:jend_2u,LM), intent(in) :: T,Q,PMID,HTM
43 REAL,dimension(ista_2l:iend_2u,jsta_2l:jend_2u,LP1),intent(in) :: PINT,ZINT
44 REAL,dimension(ista_2l:iend_2u,jsta_2l:jend_2u), intent(in) :: LMH
45 REAL,dimension(ista_2l:iend_2u,jsta_2l:jend_2u), intent(in) :: PREC
46! OUTPUT:
47! IWX - INSTANTANEOUS WEATHER TYPE.
48! ACTS LIKE A 4 BIT BINARY
49! 1111 = RAIN/FREEZING RAIN/ICE PELLETS/SNOW
50! WHERE THE ONE'S DIGIT IS FOR SNOW
51! THE TWO'S DIGIT IS FOR ICE PELLETS
52! THE FOUR'S DIGIT IS FOR FREEZING RAIN
53! AND THE EIGHT'S DIGIT IS FOR RAIN
54 integer, DIMENSION(ista:iend,jsta:jend),intent(inout) :: IWX
55! INTERNAL:
56!
57 REAL, ALLOCATABLE :: TWET(:,:,:)
58 integer,DIMENSION(ista:iend,jsta:jend) :: KARR,LICEE
59 real, dimension(ista:iend,jsta:jend) :: TCOLD,TWARM
60!
61 integer I,J,L,LMHK,LICE,IFREL,IWRML,IFRZL
62 real PSFCK,TDCHK,A,TDKL,TDPRE,TLMHK,TWRMK,AREAS8,AREAP4,AREA1, &
63 surfw,surfc,dzkl,pintk1,pintk2,pm150,qkl,tkl,pkl,area0, &
64 areap0
65
66 logical :: jcontinue=.true.
67
68! SUBROUTINES CALLED:
69! WETBULB
70!
71!
72! INITIALIZE WEATHER TYPE ARRAY TO ZERO (IE, OFF).
73! WE DO THIS SINCE WE WANT IWX TO REPRESENT THE
74! INSTANTANEOUS WEATHER TYPE ON RETURN.
75!
76!
77! ALLOCATE LOCAL STORAGE
78!
79 ALLOCATE ( twet(ista_2l:iend_2u,jsta_2l:jend_2u,lm) )
80!
81!$omp parallel do
82 DO j=jsta,jend
83 DO i=ista,iend
84 iwx(i,j) = 0
85 ENDDO
86 ENDDO
87
88!
89!!$omp parallel do
90!!$omp& private(a,lmhk,pkl,psfck,qkl,tdchk,tdkl,tdpre,tkl)
91 DO 800 j=jsta,jend
92 DO 800 i=ista,iend
93 lmhk=nint(lmh(i,j))
94!
95! SKIP THIS POINT IF NO PRECIP THIS TIME STEP
96!
97 IF (prec(i,j)<=pthresh) cycle
98!
99! FIND COLDEST AND WARMEST TEMPS IN SATURATED LAYER BETWEEN
100! 70 MB ABOVE GROUND AND 500 MB
101! ALSO FIND HIGHEST SATURATED LAYER IN THAT RANGE
102!
103!meb
104 psfck=pint(i,j,lmhk+1)
105!meb
106 tdchk=2.0
107
108 jcontinue=.true.
109 do while (jcontinue)
110
111 760 tcold(i,j)=t(i,j,lmhk)
112 twarm(i,j)=t(i,j,lmhk)
113 licee(i,j)=lmhk
114!
115 DO 775 l=1,lmhk
116 qkl=q(i,j,l)
117 qkl=amax1(h1m12,qkl)
118 tkl=t(i,j,l)
119 pkl=pmid(i,j,l)
120!
121! SKIP PAST THIS IF THE LAYER IS NOT BETWEEN 70 MB ABOVE GROUND
122! AND 500 MB
123!
124 IF (pkl<50000.0.OR.pkl>psfck-7000.0) cycle
125 IF(qkl<spval)THEN
126 a=alog(qkl*pkl/(610.78*(0.378*qkl+0.622)))
127 tdkl=(237.3*a)/(17.269-a)+273.15
128 tdpre=tkl-tdkl
129 IF (tdpre<tdchk.AND.tkl<tcold(i,j)) tcold(i,j)=tkl
130 IF (tdpre<tdchk.AND.tkl>twarm(i,j)) twarm(i,j)=tkl
131 IF (tdpre<tdchk.AND.l<licee(i,j)) licee(i,j)=l
132 ENDIF
133 775 CONTINUE
134!
135! IF NO SAT LAYER AT DEW POINT DEP=TDCHK, INCREASE TDCHK
136! AND START AGAIN (BUT DON'T MAKE TDCHK > 6)
137!
138 IF (tcold(i,j)==t(i,j,lmhk).AND.tdchk<6.0) THEN
139 tdchk=tdchk+2.0
140 ELSE
141 jcontinue=.false.
142 ENDIF
143 enddo ! enddo jcontinue
144 800 CONTINUE
145!
146! LOWEST LAYER T
147!
148 DO 850 j=jsta,jend
149 DO 850 i=ista,iend
150 karr(i,j)=0
151 IF (prec(i,j)<=pthresh) cycle
152 lmhk=nint(lmh(i,j))
153 tlmhk=t(i,j,lmhk)
154!
155! DECISION TREE TIME
156!
157 IF (tcold(i,j)>269.15) THEN
158 IF (tlmhk<=273.15) THEN
159! TURN ON THE FLAG FOR
160! FREEZING RAIN = 4
161! IF ITS NOT ON ALREADY
162! IZR=MOD(IWX(I,J),8)/4
163! IF (IZR<1) IWX(I,J)=IWX(I,J)+4
164 iwx(i,j)=iwx(i,j)+4
165 cycle
166 ELSE
167! TURN ON THE FLAG FOR
168! RAIN = 8
169! IF ITS NOT ON ALREADY
170! IRAIN=IWX(I,J)/8
171! IF (IRAIN<1) IWX(I,J)=IWX(I,J)+8
172 iwx(i,j)=iwx(i,j)+8
173 cycle
174 ENDIF
175 ENDIF
176 karr(i,j)=1
177 850 CONTINUE
178!
179! COMPUTE WET BULB ONLY AT POINTS THAT NEED IT
180!
181 CALL wetbulb(t,q,pmid,htm,karr,twet)
182!
183!!$omp parallel do
184!!$omp& private(area1,areap4,areap0,areas8,dzkl,ifrzl,iwrml,lice,
185!!$omp& lmhk,pintk1,pintk2,pm150,psfck,surfc,surfw,
186!!$omp& tlmhk,twrmk)
187 DO 1900 j=jsta,jend
188 DO 1900 i=ista,iend
189 IF(karr(i,j)>0)THEN
190 lmhk=nint(lmh(i,j))
191 lice=licee(i,j)
192!meb
193 psfck=pint(i,j,lmhk+1)
194!meb
195 tlmhk=t(i,j,lmhk)
196 twrmk=twarm(i,j)
197!
198! TWET AREA VARIABLES
199! CALCULATE ONLY WHAT IS NEEDED
200! FROM GROUND TO 150 MB ABOVE SURFACE
201! FROM GROUND TO TCOLD LAYER
202! AND FROM GROUND TO 1ST LAYER WHERE WET BULB T < 0.0
203!
204! PINTK1 IS THE PRESSURE AT THE BOTTOM OF THE LAYER
205! PINTK2 IS THE PRESSURE AT THE TOP OF THE LAYER
206!
207! AREAP4 IS THE AREA OF TWET ABOVE -4 C BELOW HIGHEST SAT LYR
208! AREAP0 IS THE AREA OF TWET ABOVE 0 C BELOW HIGHEST SAT LYR
209!
210 areas8=d00
211 areap4=d00
212 areap0=d00
213 surfw =d00
214 surfc =d00
215
216!
217 DO 1945 l=lmhk,lice,-1
218 dzkl=zint(i,j,l)-zint(i,j,l+1)
219 area1=(twet(i,j,l)-269.15)*dzkl
220 area0=(twet(i,j,l)-273.15)*dzkl
221 IF (twet(i,j,l)>=269.15) areap4=areap4+area1
222 IF (twet(i,j,l)>=273.15) areap0=areap0+area0
223 1945 CONTINUE
224!
225! IF (AREAP4<3000.0) THEN
226! TURN ON THE FLAG FOR
227! SNOW = 1
228! IF ITS NOT ON ALREADY
229! ISNO=MOD(IWX(I,J),2)
230! IF (ISNO<1) IWX(I,J)=IWX(I,J)+1
231! IWX(I,J)=IWX(I,J)+1
232! GO TO 1900
233! ENDIF
234 IF (areap0<350.0) THEN
235! TURN ON THE FLAG FOR
236! SNOW = 1
237 iwx(i,j)=iwx(i,j)+1
238 cycle
239 ENDIF
240!
241! AREAS8 IS THE NET AREA OF TWET W.R.T. FREEZING IN LOWEST 150MB
242!
243 pintk1=psfck
244 pm150=psfck-15000.
245!
246 DO 1955 l=lmhk,1,-1
247 pintk2=pint(i,j,l)
248 IF(pintk1<pm150) THEN
249 pintk1=pintk2
250 ELSE
251 dzkl=zint(i,j,l)-zint(i,j,l+1)
252!
253! SUM PARTIAL LAYER IF IN 150 MB AGL LAYER
254!
255 IF(pintk2<pm150) &
256 dzkl=t(i,j,l)*(q(i,j,l)*d608+h1)*rog* &
257 alog(pintk1/pm150)
258 area1=(twet(i,j,l)-273.15)*dzkl
259 areas8=areas8+area1
260 pintk1=pintk2
261 ENDIF
262 1955 CONTINUE
263!
264! SURFW IS THE AREA OF TWET ABOVE FREEZING BETWEEN THE GROUND
265! AND THE FIRST LAYER ABOVE GROUND BELOW FREEZING
266! SURFC IS THE AREA OF TWET BELOW FREEZING BETWEEN THE GROUND
267! AND THE WARMEST SAT LAYER
268!
269 ifrzl=0
270 iwrml=0
271!
272 DO 2050 l=lmhk,1,-1
273 IF (ifrzl==0.AND.t(i,j,l)<273.15) ifrzl=1
274 IF (iwrml==0.AND.t(i,j,l)>=twrmk) iwrml=1
275!
276 IF (iwrml==0.OR.ifrzl==0) THEN
277 dzkl=zint(i,j,l)-zint(i,j,l+1)
278 area1=(twet(i,j,l)-273.15)*dzkl
279 IF(ifrzl==0.AND.twet(i,j,l)>=273.15)surfw=surfw+area1
280 IF(iwrml==0.AND.twet(i,j,l)<=273.15)surfc=surfc+area1
281 ENDIF
282 2050 CONTINUE
283 IF(surfc<-3000.0.OR. &
284 & (areas8<-3000.0.AND.surfw<50.0)) THEN
285! TURN ON THE FLAG FOR
286! ICE PELLETS = 2
287! IF ITS NOT ON ALREADY
288! IIP=MOD(IWX(I,J),4)/2
289! IF (IIP<1) IWX(I,J)=IWX(I,J)+2
290 iwx(i,j)=iwx(i,j)+2
291 cycle
292 ENDIF
293!
294 IF(tlmhk<273.15) THEN
295! TURN ON THE FLAG FOR
296! FREEZING RAIN = 4
297! IF ITS NOT ON ALREADY
298! IZR=MOD(IWX(K),8)/4
299! IF (IZR<1) IWX(K)=IWX(K)+4
300 iwx(i,j)=iwx(i,j)+4
301 ELSE
302! TURN ON THE FLAG FOR
303! RAIN = 8
304! IF ITS NOT ON ALREADY
305! IRAIN=IWX(K)/8
306! IF (IRAIN<1) IWX(K)=IWX(K)+8
307 iwx(i,j)=iwx(i,j)+8
308 ENDIF
309 ENDIF
310 1900 CONTINUE
311! print *, 'revised check ', IWX(500,800)
312
313!---------------------------------------------------------
314 DEALLOCATE (twet)
315
316
317 RETURN
318 END