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