65 SUBROUTINE calllws(U,V,H,LLWS)
68 USE vrbls2d, only: fis, u10, v10
70 use ctlblk_mod
, only: jsta, jend, im, jm, lsm, spval
76 REAL,
DIMENSION(IM,JM,LSM),
INTENT(IN) :: u,v,h
77 REAL,
DIMENSION(IM,JM),
INTENT(INOUT) :: llws
78 REAL :: z1,z2,hz1,dh,u2,v2,w2,rt
89 z1 = 10.0 + fis(i,j)*gi
91 IF(z1<h(i,j,lsm))
THEN
95 IF(z1>=h(i,j,lp).AND.z1<h(i,j,lp-1))
THEN
101 hz1 = h(i,j,k1-1) - z1
105 IF((hz1+10)>609.6)
THEN
106 u2= u10(i,j) + (u(i,j,k1-1)-u10(i,j))*599.6/hz1
107 v2= v10(i,j) + (v(i,j,k1-1)-v10(i,j))*599.6/hz1
108 z2= fis(i,j)*gi + 609.6
111 dh=dh+(h(i,j,lp-1) - h(i,j,lp))
112 IF((dh+hz1+10)>609.6)
THEN
114 rt=(z2-h(i,j,lp))/(h(i,j,lp-1)-h(i,j,lp))
115 u2=u(i,j,lp)+rt*(u(i,j,lp-1)-u(i,j,lp))
116 v2=v(i,j,lp)+rt*(v(i,j,lp-1)-v(i,j,lp))
125 if(u10(i,j)<spval.and.v10(i,j)<spval) &
126 llws(i,j)=sqrt((u2-u10(i,j))**2+(v2-v10(i,j))**2)/ &
160 SUBROUTINE calicing (T1,RH,OMGA, ICING)
161 use ctlblk_mod
, only: jsta, jend, im, spval
167 REAL,
DIMENSION(IM,jsta:jend),
INTENT(IN) :: t1,rh,omga
168 REAL,
DIMENSION(IM,jsta:jend),
INTENT(INOUT) :: icing
175 IF(omga(i,j)<spval.AND.t1(i,j)<spval.AND.rh(i,j)<spval)
THEN
176 IF(omga(i,j) < 0.0 .AND. &
177 (t1(i,j) <= 273.0 .AND. t1(i,j) >= 251.0) &
178 .AND. rh(i,j) >= 70.0)
THEN
219 SUBROUTINE calcat(U,V,H,U_OLD,V_OLD,H_OLD,CAT)
220 use masks, only: dx, dy
221 use ctlblk_mod
, only: spval, jsta_2l, jend_2u, jsta_m, jend_m, &
223 use gridspec_mod
, only: gridtype
231 REAL,
DIMENSION(IM,jsta_2l:jend_2u),
INTENT(IN) :: u,v,h, &
234 REAL,
DIMENSION(IM,jsta_2l:jend_2u),
INTENT(INOUT) :: cat
236 REAL dsh, dst, def, cvg, vws, trbindx
237 INTEGER ihe(jm),ihw(jm)
239 integer istart,istop,jstart,jstop
240 real vws1,vws2,vws3,vws4
247 IF(gridtype ==
'A')
THEN
254 ELSE IF(gridtype==
'E')
THEN
261 ELSE IF(gridtype==
'B')
THEN
269 print*,
'no gridtype specified, exit calcat comp'
281 DO 100 j=jstart,jstop
284 IF(gridtype==
'B')
THEN
285 IF(u(i,j)<spval.and.u(i,j-1)<spval.and.u(i-1,j)<spval.and.u(i-1,j-1)<spval.and.&
286 v(i,j)<spval.and.v(i,j-1)<spval.and.v(i-1,j)<spval.and.v(i-1,j-1)<spval)
THEN
288 dsh=(0.5*(v(i,j)+v(i,j-1))-0.5*(v(i-1,j)+v(i-1,j-1)))*10000./dx(i,j) &
289 +(0.5*(u(i,j)+u(i-1,j))-0.5*(u(i,j-1)+u(i-1,j-1)))*10000./dy(i,j)
291 dst =(0.5*(u(i,j)+u(i,j-1))-0.5*(u(i-1,j)+u(i-1,j-1)))*10000./dx(i,j) &
292 -(0.5*(v(i,j)+v(i-1,j))-0.5*(v(i,j-1)+v(i-1,j-1)))*10000./dy(i,j)
293 def = sqrt(dsh*dsh + dst*dst)
296 cvg = -((0.5*(u(i,j)+u(i,j-1))-0.5*(u(i-1,j)+u(i-1,j-1)))*10000./dx(i,j) &
297 +(0.5*(v(i,j)+v(i-1,j))-0.5*(v(i,j-1)+v(i-1,j-1)))*10000./dy(i,j))
303 IF(u(i,j+1)<spval.and.u(i,j-1)<spval.and.u(i+ihe(j),j)<spval.and.u(i+ihw(j),j)<spval.and.&
304 v(i,j+1)<spval.and.v(i,j-1)<spval.and.v(i+ihe(j),j)<spval.and.v(i+ihw(j),j)<spval)
THEN
306 dsh = (v(i+ihe(j),j) - v(i+ihw(j),j))*10000./(2*dx(i,j)) &
307 + (u(i,j+1) - u(i,j-1))*10000./(2*dy(i,j))
310 dst = (u(i+ihe(j),j) - u(i+ihw(j),j))*10000./(2*dx(i,j)) &
311 - (v(i,j+1) - v(i,j-1))*10000./(2*dy(i,j))
313 def = sqrt(dsh*dsh + dst*dst)
316 cvg = -( (u(i+ihe(j),j) - u(i+ihw(j),j))*10000./(2*dx(i,j)) &
317 +(v(i,j+1) - v(i,j-1))*10000./(2*dy(i,j)) )
324 IF(gridtype ==
'A')
THEN
326 IF(u_old(i,j)<spval.and.u(i,j)<spval.and.&
327 v_old(i,j)<spval.and.v(i,j)<spval.and.&
328 h_old(i,j)<spval.and.h(i,j)<spval)
THEN
329 vws = ( sqrt(u_old(i,j)**2+v_old(i,j)**2 ) - &
330 sqrt(u(i,j)**2+v(i,j)**2 ) ) * &
331 1000.0/(h_old(i,j) - h(i,j))
335 else IF(gridtype ==
'E')
THEN
337 IF(u_old(i+ihe(j),j)<spval.and.u(i+ihe(j),j)<spval.and.&
338 v_old(i+ihe(j),j)<spval.and.v(i+ihe(j),j)<spval)
THEN
340 vws1 = ( sqrt(u_old(i+ihe(j),j)**2+v_old(i+ihe(j),j)**2 ) -&
341 sqrt(u(i+ihe(j),j)**2+v(i+ihe(j),j)**2 ) )
346 IF(u_old(i+ihw(j),j)<spval.and.u(i+ihw(j),j)<spval.and.&
347 v_old(i+ihw(j),j)<spval.and.v(i+ihw(j),j)<spval)
THEN
348 vws2 = ( sqrt(u_old(i+ihw(j),j)**2+v_old(i+ihw(j),j)**2 ) -&
349 sqrt(u(i+ihw(j),j)**2+v(i+ihw(j),j)**2 ) )
354 IF(u_old(i,j-1)<spval.and.u(i,j-1)<spval.and.&
355 v_old(i,j-1)<spval.and.v(i,j-1)<spval)
THEN
356 vws3 = ( sqrt(u_old(i,j-1)**2+v_old(i,j-1)**2 ) - &
357 sqrt(u(i,j-1)**2+v(i,j-1)**2 ) )
362 IF(u_old(i,j+1)<spval.and.u(i,j+1)<spval.and.&
363 v_old(i,j+1)<spval.and.v(i,j+1)<spval)
THEN
364 vws4 = ( sqrt(u_old(i,j+1)**2+v_old(i,j+1)**2 ) - &
365 sqrt(u(i,j+1)**2+v(i,j+1)**2 ) )
370 IF(vws1<spval.and.vws2<spval.and.vws3<spval.and.vws4<spval.and.&
371 h_old(i,j)<spval.and.h(i,j)<spval)
THEN
372 vws=1000.0*(vws1+vws2+vws3+vws4)/4.0/(h_old(i,j) - h(i,j))
376 ELSE IF(gridtype ==
'B')
THEN
377 IF(u_old(i+ihe(j),j)<spval.and.u(i+ihe(j),j)<spval.and.&
378 v_old(i+ihe(j),j)<spval.and.v(i+ihe(j),j)<spval)
THEN
379 vws1 = ( sqrt(u_old(i+ihe(j),j)**2+v_old(i+ihe(j),j)**2 ) -&
380 sqrt(u(i+ihe(j),j)**2+v(i+ihe(j),j)**2 ) )
385 IF(u_old(i+ihw(j),j)<spval.and.u(i+ihw(j),j)<spval.and.&
386 v_old(i+ihw(j),j)<spval.and.v(i+ihw(j),j)<spval)
THEN
387 vws2 = ( sqrt(u_old(i+ihw(j),j)**2+v_old(i+ihw(j),j)**2 ) -&
388 sqrt(u(i+ihw(j),j)**2+v(i+ihw(j),j)**2 ) )
393 IF(u_old(i,j-1)<spval.and.u(i,j-1)<spval.and.&
394 v_old(i,j-1)<spval.and.v(i,j-1)<spval)
THEN
395 vws3 = ( sqrt(u_old(i,j-1)**2+v_old(i,j-1)**2 ) - &
396 sqrt(u(i,j-1)**2+v(i,j-1)**2 ) )
401 IF(u_old(i-1,j-1)<spval.and.u(i-1,j-1)<spval.and.&
402 v_old(i-1,j-1)<spval.and.v(i-1,j-1)<spval)
THEN
403 vws4 = ( sqrt(u_old(i-1,j-1)**2+v_old(i-1,j-1)**2 ) - &
404 sqrt(u(i-1,j-1)**2+v(i-1,j-1)**2 ) )
409 IF(vws1<spval.and.vws2<spval.and.vws3<spval.and.vws4<spval.and.&
410 h_old(i,j)<spval.and.h(i,j)<spval)
THEN
411 vws=1000.0*(vws1+vws2+vws3+vws4)/4.0/(h_old(i,j) - h(i,j))
417 IF(vws<spval.and.def<spval.and.cvg<spval)
THEN
418 trbindx = abs(vws)*(def + abs(cvg))
422 ELSE IF(trbindx<=8.)
THEN
424 ELSE IF(trbindx<=12.)
THEN
451 SUBROUTINE calceiling (CLDZ,TCLD,CEILING)
454 use ctlblk_mod
, only: jsta, jend, spval, im, modelname
460 REAL,
DIMENSION(IM,jsta:jend),
INTENT(IN) :: cldz, tcld
461 REAL,
DIMENSION(IM,jsta:jend),
INTENT(INOUT) :: ceiling
468 IF(abs(tcld(i,j)-spval) <= small)
THEN
470 ELSE IF(tcld(i,j) >= 50.)
THEN
471 if(modelname ==
'RAPR')
then
472 ceiling(i,j) = cldz(i,j) - fis(i,j)*gi
474 ceiling(i,j) = cldz(i,j)
477 ceiling(i,j) = 20000.0
480 IF(ceiling(i,j) < 0.0) ceiling(i,j)=20000.0
505 SUBROUTINE calfltcnd (CEILING,FLTCND)
507 use ctlblk_mod
, only: jsta, jend, im, spval
513 REAL,
DIMENSION(IM,jsta:jend),
INTENT(IN) :: ceiling
514 REAL,
DIMENSION(IM,jsta:jend),
INTENT(INOUT) :: fltcnd
524 IF(ceiling(i,j)<spval.and.vis(i,j)<spval)
THEN
525 ceil = ceiling(i,j) * 3.2808
526 visi = vis(i,j) / 1609.0
528 IF(ceil<500.0 .OR. visi<1.0 )
THEN
531 ELSE IF( (ceil>=500.AND.ceil<1000.0) .OR. &
532 (visi>=1.0.AND.visi<3.0) )
THEN
535 ELSE IF( (ceil>=1000.AND.ceil<=3000.0) .OR. &
536 (visi>=3.0.AND.visi<=5.0) )
THEN
539 ELSE IF( ceil>3000.0 .OR. visi>5.0)
THEN