UPP  001
All Data Structures Files Functions Pages
AVIATION.f
Go to the documentation of this file.
1 
3 !
65  SUBROUTINE calllws(U,V,H,LLWS)
66 
67 !
68  USE vrbls2d, only: fis, u10, v10
69  use params_mod, only: gi
70  use ctlblk_mod, only: jsta, jend, im, jm, lsm, spval
71 !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72  implicit none
73 !
74 ! DECLARE VARIABLES.
75 !
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
79  INTEGER :: k1,k2
80  integer i,j,lp
81 
82 !***************************************************************
83 !
84 !
85 
86  DO 100 j=jsta,jend
87  DO i=1,im
88 
89  z1 = 10.0 + fis(i,j)*gi !Height of 10m levels geographic height (from sea level)
90 
91  IF(z1<h(i,j,lsm)) THEN !First search location of 10m wind level
92  k1 = lsm + 1 !to see it is in which pressure levels
93  ELSE
94  DO lp = lsm,2,-1 !If not found, keep searching upward
95  IF(z1>=h(i,j,lp).AND.z1<h(i,j,lp-1)) THEN
96  k1 = lp
97  END IF
98  END DO
99  END IF
100 
101  hz1 = h(i,j,k1-1) - z1 !Distance between K1-1 and 10m level
102 
103  dh = 0.0
104 
105  IF((hz1+10)>609.6) THEN !Then, search 2000ft(609.6m) location
106  u2= u10(i,j) + (u(i,j,k1-1)-u10(i,j))*599.6/hz1 !found it between K1-1 and K1, then linear
107  v2= v10(i,j) + (v(i,j,k1-1)-v10(i,j))*599.6/hz1 !interpolate to get wind at 2000ft U2,V2
108  z2= fis(i,j)*gi + 609.6
109  ELSE !otherwise, keep on search upward
110  DO lp = k1-1,2,-1
111  dh=dh+(h(i,j,lp-1) - h(i,j,lp))
112  IF((dh+hz1+10)>609.6) THEN !found the 2000ft level
113  z2=fis(i,j)*gi+609.6
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))
117  k2=lp
118  exit
119  END IF
120  END DO
121  END IF
122 
123 !computer vector difference
124  llws(i,j) = spval
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)/ &
127  609.6 * 1.943*609.6 !unit: knot/2000ft
128  ENDDO
129 
130 100 CONTINUE
131 
132  RETURN
133  END
134 
160  SUBROUTINE calicing (T1,RH,OMGA, ICING)
161  use ctlblk_mod, only: jsta, jend, im, spval
162 !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163  implicit none
164 !
165 ! DECLARE VARIABLES.
166 !
167  REAL, DIMENSION(IM,jsta:jend), INTENT(IN) :: t1,rh,omga
168  REAL, DIMENSION(IM,jsta:jend), INTENT(INOUT) :: icing
169  integer i,j
170 !***************************************************************
171 !
172 !
173  DO j=jsta,jend
174  DO i=1,im
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
179 
180  icing(i,j) = 1.0
181  ELSE
182  icing(i,j) = 0.0
183  END IF
184  ELSE
185  icing(i,j) = spval
186  ENDIF
187  ENDDO
188  ENDDO
189 
190  RETURN
191  END
192 
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, &
222  im, jm
223  use gridspec_mod, only: gridtype
224 !
225 !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226  implicit none
227 
228 !
229 ! DECLARE VARIABLES.
230 !
231  REAL,DIMENSION(IM,jsta_2l:jend_2u),INTENT(IN) :: u,v,h, &
232  u_old,v_old,h_old
233 ! INTEGER,INTENT(IN) :: L
234  REAL,DIMENSION(IM,jsta_2l:jend_2u),INTENT(INOUT) :: cat
235 
236  REAL dsh, dst, def, cvg, vws, trbindx
237  INTEGER ihe(jm),ihw(jm)
238  integer i,j
239  integer istart,istop,jstart,jstop
240  real vws1,vws2,vws3,vws4
241 
242 !***************************************************************
243 !
244 !
245  cat=spval
246  DO j=jsta_2l,jend_2u
247  IF(gridtype == 'A')THEN
248  ihw(j)=-1
249  ihe(j)=1
250  istart=2
251  istop=im-1
252  jstart=jsta_m
253  jstop=jend_m
254  ELSE IF(gridtype=='E')THEN
255  ihw(j)=-mod(j,2)
256  ihe(j)=ihw(j)+1
257  istart=2
258  istop=im-1
259  jstart=jsta_m
260  jstop=jend_m
261  ELSE IF(gridtype=='B')THEN
262  ihw(j)=-1
263  ihe(j)=0
264  istart=2
265  istop=im-1
266  jstart=jsta_m
267  jstop=jend_m
268  ELSE
269  print*,'no gridtype specified, exit calcat comp'
270  return
271  END IF
272  ENDDO
273 
274  call exch_f(u)
275  call exch_f(v)
276  call exch_f(u_old)
277  call exch_f(v_old)
278  call exch_f(h)
279  call exch_f(h_old)
280 
281  DO 100 j=jstart,jstop
282  DO i=istart,istop
283 !
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
287 !dsh=dv/dx+du/dy
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)
290 !dst=du/dx-dv/dy
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)
294 
295 !cvg=-(du/dx+dv/dy)
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))
298  ELSE
299  def = spval
300  cvg = spval
301  ENDIF
302  ELSE
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
305 !dsh=dv/dx+du/dy
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))
308 
309 !dst=du/dx-dv/dy
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))
312 
313  def = sqrt(dsh*dsh + dst*dst)
314 
315 !cvg=-(du/dx+dv/dy)
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)) )
318  ELSE
319  def = spval
320  cvg = spval
321  ENDIF
322  END IF
323 
324  IF(gridtype == 'A')THEN
325 !vws=d|U|/dz
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))
332  ELSE
333  vws = spval
334  ENDIF
335  else IF(gridtype == 'E')THEN
336 !vws=d|U|/dz
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
339 
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 ) )
342  ELSE
343  vws1 = spval
344  ENDIF
345 !vws=d|U|/dz
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 ) )
350  ELSE
351  vws2 = spval
352  ENDIF
353 !vws=d|U|/dz
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 ) )
358  ELSE
359  vws3 = spval
360  ENDIF
361 !vws=d|U|/dz
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 ) )
366  ELSE
367  vws4 = spval
368  ENDIF
369 
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))
373  ELSE
374  vws = spval
375  ENDIF
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 ) )
381  ELSE
382  vws1 = spval
383  ENDIF
384 !vws=d|U|/dz
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 ) )
389  ELSE
390  vws2 = spval
391  ENDIF
392 !vws=d|U|/dz
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 ) )
397  ELSE
398  vws3 = spval
399  ENDIF
400 !vws=d|U|/dz
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 ) )
405  ELSE
406  vws4 = spval
407  ENDIF
408 
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))
412  ELSE
413  vws=spval
414  ENDIF
415  END IF
416 
417  IF(vws<spval.and.def<spval.and.cvg<spval)THEN
418  trbindx = abs(vws)*(def + abs(cvg))
419 
420  IF(trbindx<=4.) THEN
421  cat(i,j) = 0.0
422  ELSE IF(trbindx<=8.) THEN
423  cat(i,j)=1.0
424  ELSE IF(trbindx<=12.) THEN
425  cat(i,j)=2.0
426  ELSE
427  cat(i,j)=3.0
428  END IF
429  ELSE
430  cat(i,j)=spval
431  ENDIF
432  ENDDO
433 
434 100 CONTINUE
435 
436  RETURN
437  END
438 
451  SUBROUTINE calceiling (CLDZ,TCLD,CEILING)
452  USE vrbls2d, only: fis
453  use params_mod, only: small, gi
454  use ctlblk_mod, only: jsta, jend, spval, im, modelname
455 !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
456  implicit none
457 !
458 ! DECLARE VARIABLES.
459 !
460  REAL, DIMENSION(IM,jsta:jend), INTENT(IN) :: cldz, tcld
461  REAL, DIMENSION(IM,jsta:jend), INTENT(INOUT) :: ceiling
462  integer i,j
463 !***************************************************************
464 !
465 !
466  DO j=jsta,jend
467  DO i=1,im
468  IF(abs(tcld(i,j)-spval) <= small) THEN
469  ceiling(i,j)=spval
470  ELSE IF(tcld(i,j) >= 50.) THEN
471  if(modelname == 'RAPR')then
472  ceiling(i,j) = cldz(i,j) - fis(i,j)*gi
473  else
474  ceiling(i,j) = cldz(i,j) ! for RAP/HRRR - FIS(I,J)*GI
475  endif
476  ELSE
477  ceiling(i,j) = 20000.0
478  END IF
479 
480  IF(ceiling(i,j) < 0.0) ceiling(i,j)=20000.0
481 
482  ENDDO
483  ENDDO
484 
485  RETURN
486  END
487 
505  SUBROUTINE calfltcnd (CEILING,FLTCND)
506  use vrbls2d, only: vis
507  use ctlblk_mod, only: jsta, jend, im, spval
508 !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
509  implicit none
510 !
511 ! DECLARE VARIABLES.
512 !
513  REAL, DIMENSION(IM,jsta:jend), INTENT(IN) :: ceiling
514  REAL, DIMENSION(IM,jsta:jend), INTENT(INOUT) :: fltcnd
515  REAL ceil,visi
516  integer i,j
517 !
518 !***************************************************************
519 !
520 !
521  DO j=jsta,jend
522  DO i=1,im
523 
524  IF(ceiling(i,j)<spval.and.vis(i,j)<spval)THEN
525  ceil = ceiling(i,j) * 3.2808 !from m -> feet
526  visi = vis(i,j) / 1609.0 !from m -> miles
527 
528  IF(ceil<500.0 .OR. visi<1.0 ) THEN
529  fltcnd(i,j) = 1.0
530 
531  ELSE IF( (ceil>=500.AND.ceil<1000.0) .OR. &
532  (visi>=1.0.AND.visi<3.0) ) THEN
533  fltcnd(i,j) = 2.0
534 
535  ELSE IF( (ceil>=1000.AND.ceil<=3000.0) .OR. &
536  (visi>=3.0.AND.visi<=5.0) ) THEN
537  fltcnd(i,j) = 3.0
538 
539  ELSE IF( ceil>3000.0 .OR. visi>5.0) THEN
540  fltcnd(i,j) = 4.0
541 
542  END IF
543  ELSE
544  fltcnd(i,j) = spval
545  ENDIF
546  ENDDO
547  ENDDO
548 
549  RETURN
550  END
Definition: MASKS_mod.f:1
subroutine exch_f(a)
Definition: EXCH.f:65