9 public derive_fields, mixing_ratio
13 type :: precipitations_t
18 integer :: convection = 4
19 end type precipitations_t
20 type(precipitations_t
),
parameter :: precips = precipitations_t()
25 subroutine derive_fields(imp_physics,t, rh, pres, hgt, totalWater, totalCond,&
26 nz, topok, hprcp, hcprcp, cin, cape, &
27 ept, wbt, twp, pc, kx, lx, tott, prcptype)
43 integer,
intent(in) :: imp_physics
44 integer,
intent(in) :: nz, topok
45 real,
intent(in),
dimension(nz) :: t, rh, pres, hgt
46 real,
intent(in),
dimension(nz) :: totalwater,totalcond
47 real,
intent(in) :: hprcp, hcprcp, cin, cape
48 real,
intent(out) :: ept(nz), wbt(nz), twp(nz)
49 real,
intent(out) :: pc, kx, lx, tott
50 integer,
intent(out) :: prcptype
52 real,
allocatable :: td(:), tlc(:), wvm(:)
59 td(:) = getdewpointtemp(t(:), rh(:))
60 tlc(:) = get_tlcl(t(:), td(:))
61 wvm(:) = mixing_ratio(td(:), pres(:))
63 ept(:) = getthetae(pres(:), t(:), wvm(:), tlc(:))
64 wbt(:) = getthetaw(t(:), td(:), pres(:))
66 call calc_totalwaterpath(hgt, pres, t, totalwater, nz, twp)
69 pc = getprecipcond(totalcond, nz, topok)
72 call calc_indice(t, td, pres, wvm, nz, topok, kx, lx, tott)
74 prcptype=getpreciptype(imp_physics,pres,hgt,t,rh,wbt,nz,&
75 hprcp,hcprcp,pc,cin,cape,lx)
82 end subroutine derive_fields
87 elemental real function getdewpointtemp(t, rh)
89 real,
intent(in) :: t, rh
96 vapr = (max(1.0e-6,rh) / 100.0) * getvaporpres(t)
99 getdewpointtemp = 243.5 * rm / (17.67 - rm) + 273.15
102 end function getdewpointtemp
109 elemental real function getthetae(pres, t, wvm, tAtLCL)
111 real,
intent(in) :: pres, t, wvm, tatlcl
117 e = (2.0/7.0) * ( 1.0 - 0.28 * rmix * 0.001)
118 thtam = t * ((100000.0/pres)**e)
119 getthetae = thtam * exp( (3.376/tatlcl - 0.00254) * &
120 (rmix * (1. + (.81 * rmix * 0.001)) ))
123 end function getthetae
128 elemental real function getvaporpres(t)
130 real,
intent(in) :: t
136 getvaporpres = 611.2 * exp( 17.67*tc/(tc+243.5))
139 end function getvaporpres
144 elemental real function get_tlcl(t, td)
146 real,
intent(in) :: t, td
148 get_tlcl = 1. / (1./(td - 56.) + log(t/td)/800.0) + 56.0
151 end function get_tlcl
157 elemental real function mixing_ratio(td, pres)
159 real,
intent(in) :: td, pres
163 corr = 1.001 + ( (pres/100.0 - 100) /900.) * 0.0034
164 e = getvaporpres(td) * corr
165 mixing_ratio = 0.62197 * (e/(pres-e)) * 1000.0
168 end function mixing_ratio
174 elemental real function getthetaw(t, td, pres)
176 real,
intent(in) :: t, td, pres
178 real :: vl, cp, dt, top, bottom, twet, twetnew
191 twet = (top + bottom) / 2.0
192 rmixr = mixing_ratio(td, pres) / 1000.0
193 smixr = mixing_ratio(twet, pres) / 1000.0
195 twetnew = t - (vl/cp) * (smixr-rmixr)
196 if(twetnew < twet)
then
202 dt = abs(twet - twetnew)
211 end function getthetaw
215 real function getprecipcond(totalCond, nz, topoK)
217 real,
intent(in) :: totalcond(nz)
218 integer,
intent(in) :: nz, topok
224 do k = topok, topok-2, -1
225 getprecipcond = totalcond(k) + getprecipcond
229 end function getprecipcond
233 subroutine calc_indice(t, td, pres, wvm, nz, topoK, &
234 kindex, liftedindex, totaltotals)
236 integer,
intent(in) :: nz, topok
237 real,
intent(in) :: t(nz), td(nz), pres(nz), wvm(nz)
238 real,
intent(out) :: kindex, liftedindex, totaltotals
241 real :: t500hpa, t700hpa, t850hpa, dpt700hpa, dpt850hpa
243 real :: surfacetemp,surfacedewpttemp,surfacepressure,surfacemixr
244 real :: tempatlcl, theta, pressatlcl, thetaeatlcl, tempfromthetae, tem
258 if ((pres(k)- 50000.0 >= 0.) .and. (pres(k-1)- 50000.0 < 0.) )
then
259 if (abs(pres(k)- 50000.0) <= 0.1)
then
261 elseif (abs(pres(k-1)- 50000.0) <= 0.1)
then
264 t500hpa = t(k) - ((pres(k)-50000.0)/(pres(k)-pres(k-1))) * (t(k)-t(k-1))
269 if ((pres(k)- 70000.0 >= 0.) .and. (pres(k-1)- 70000.0 < 0.) )
then
270 if (abs(pres(k)- 70000.0) <= 0.1)
then
273 elseif (abs(pres(k-1)- 70000.0) <= 0.1)
then
277 tem = (pres(k)-70000.0)/(pres(k)-pres(k-1))
278 t700hpa = t(k) - tem * (t(k)-t(k-1))
279 dpt700hpa = td(k) - tem * (td(k)-td(k-1))
284 if ((pres(k)- 85000.0 >= 0.) .and. (pres(k-1)- 85000.0 < 0.) )
then
285 if (abs(pres(k)- 85000.0) <= 0.1)
then
288 elseif (abs(pres(k-1)- 85000.0) <= 0.1)
then
292 tem = (pres(k)-85000.0)/(pres(k)-pres(k-1))
293 t850hpa = t(k) - tem * (t(k)-t(k-1))
294 dpt850hpa = td(k) - tem * (td(k)-td(k-1))
301 kindex = (t850hpa-t500hpa) + (dpt850hpa-273.15) - (t700hpa-dpt700hpa)
304 totaltotals = t850hpa + dpt850hpa - t500hpa * 2.0
308 surfacetemp = t(topok)
309 surfacedewpttemp = td(topok)
310 surfacepressure = pres(topok)
311 surfacemixr = wvm(topok)
313 tempatlcl = get_tlcl(surfacetemp, surfacedewpttemp)
314 theta = surfacetemp * (100000.0/surfacepressure)**0.286
315 pressatlcl = 100000.0 * (tempatlcl / theta)**3.4965
316 thetaeatlcl = getthetae(pressatlcl,tempatlcl,surfacemixr,tempatlcl)
318 tempfromthetae = gettemperature(thetaeatlcl, 50000.0)
320 liftedindex = t500hpa - tempfromthetae
323 end subroutine calc_indice
326 real function gettemperature(thetaE, pres)
328 real,
intent(in) :: thetae, pres
330 real :: guess, w1, w2, thetae1, thetae2, cor
333 guess = (thetae - 0.5 * (max(thetae - 270.0, 0.0))**1.05) * &
337 w1 = calcrsubs(pres, guess)
338 w2 = calcrsubs(pres, guess + 1.0)
339 thetae1 = getthetae(pres, guess, w1, guess)
340 thetae2 = getthetae(pres, guess + 1.0, w2, guess + 1.0)
341 cor = (thetae - thetae1)/(thetae2 - thetae1)
344 if (abs(cor) < 0.01)
then
349 gettemperature = guess
352 end function gettemperature
355 elemental real function calcrsubs(pres, temp)
357 real,
intent(in) :: pres, temp
361 esubs = calcesubs(temp)
362 calcrsubs = (0.622*esubs)/(pres - esubs)
365 end function calcrsubs
368 elemental real function calcesubs(temp)
370 real,
intent(in) :: temp
373 real,
parameter :: coeffs(9) = &
374 (/ 610.5851, 44.40316, 1.430341, 0.02641412, 2.995057e-4,&
375 2.031998e-6, 6.936113e-9, 2.564861e-12, -3.704404e-14 /)
378 x= max(-80.0, temp - 273.15)
380 calcesubs = coeffs(1) + x*(coeffs(2) + x*(coeffs(3) + x*(coeffs(4) + &
381 x*(coeffs(5) + x*(coeffs(6) + x*(coeffs(7) + &
382 x*(coeffs(8) + x*coeffs(9))))))))
384 end function calcesubs
387 subroutine calc_totalwaterpath(hgt, pres, t, totalWater, nz, twp)
389 integer,
intent(in) :: nz
390 real,
intent(in) :: hgt(nz), pres(nz), t(nz), totalwater(nz)
391 real,
intent(out):: twp(nz)
393 real :: condensateintegrated
394 integer :: topidx, baseidx
404 lp200:
do m = k, 2, -1
405 if (totalwater(m) > 0.001)
then
413 if (totalwater(m) > 0.001)
then
421 if(topidx == -1 .or. baseidx == -1) cycle
422 condensateintegrated = 0.0
423 lp400:
do m = topidx, baseidx
424 if (t(m) <= 273.15)
then
425 condensateintegrated = condensateintegrated + &
426 ( ((totalwater(m)*pres(m)) / (287.05 * t(m)) ) * &
427 (hgt(m-1) - hgt(m) + 0.001))
431 twp(k) = condensateintegrated
436 end subroutine calc_totalwaterpath
439 integer function getpreciptype(imp_physics,pres, hgt, t, rh, wbt, nz, &
440 hprcp, hcprcp, pc, cin, cape, lx)
442 integer,
intent(in) :: imp_physics
443 integer,
intent(in) :: nz
444 real,
intent(in) :: pres(nz), hgt(nz), t(nz), rh(nz), wbt(nz)
445 real,
intent(in) :: hprcp,hcprcp, pc, cin, cape, lx
447 integer,
allocatable :: wxtype(:)
449 integer :: idxmelt, idxrefz, iceidx
450 real :: coldtemp, tcoldarea, wetbuldarea
452 real :: precip, precip_standard
461 if(imp_physics == 98 .or. imp_physics == 99)
then
463 precip_standard =0.045
464 else if(imp_physics == 11 .or. imp_physics == 8)
then
466 precip_standard = 0.01
470 if_conv:
if ( hcprcp >= 1.0 .and. cape >= 1000.0 .and. cin > -100. &
472 getpreciptype = precips% CONVECTION
477 if( wbt(k) > 273.15)
then
481 if (t(m) < 273.15)
then
492 lp200:
do k = nz, 1, -1
496 if_pres_wbt:
if ( pres(k) >= 15000. .or. &
497 (wbt(k) >= 200. .and. wbt(k) <= 1000.))
then
505 if( (pres(m)>=30000. .or. (hgt(m)>hgt(nz)+700.))&
506 .and. (rh(m) > 90.) .and. (t(m) < coldtemp))
then
513 if_iceidx:
if (iceidx /= k)
then
519 if (wbt(m) >= 273.15)
then
520 tcoldarea = tcoldarea + (wbt(m)-273.15) * &
529 if ( (hgt(m) <= hgt(nz)+1500) .and. &
531 wetbuldarea = wetbuldarea + (wbt(m)-273.15) * &
537 if (coldtemp > 265.15)
then
538 wxtype(k) = precips% OTHER
539 else if (tcoldarea < 350.)
then
540 if(t(k) <= 273.15)
then
541 wxtype(k) = precips% SNOW
543 wxtype(k) = precips% RAIN
545 else if (wetbuldarea <= -250.)
then
546 wxtype(k) = precips% OTHER
547 else if(t(k) <= 273.15)
then
548 wxtype(k) = precips% OTHER
550 wxtype(k) = precips% RAIN
554 wxtype(k) = precips% NONE
560 getpreciptype = precips% NONE
563 if(precip >= precip_standard)
then
565 if (wxtype(k) > getpreciptype)
then
566 getpreciptype = wxtype(k)
575 end function getpreciptype
577 end module derivedfields
584 use derivedfields
, only : mixing_ratio
589 public calc_cloudlayers
592 integer,
parameter :: maxlayers = 30
599 real,
allocatable :: layerq(:)
601 integer :: topidx(maxlayers)
602 integer :: baseidx(maxlayers)
603 real :: ctt(maxlayers)
609 subroutine calc_cloudlayers(rh,t,pres,ept,vv, nz, topoK, xlat, xlon,&
612 integer,
intent(in) :: nz, topok
613 real,
intent(in) :: rh(nz),t(nz),pres(nz), ept(nz), vv(nz)
614 real,
intent(in) :: xlat,xlon
615 integer,
intent(out) :: region
616 type(clouds_t
),
intent(inout) :: clouds
619 integer :: in_cld,cur_base, num_lyr
621 integer :: k, kk, n, kstart
624 if (abs(xlat)<23.5)
then
627 elseif ( (abs(xlat)>=23.5).and.(abs(xlat)<66))
then
643 if((rh(1) >= t_rh) .and. (rh(2) >= t_rh))
then
644 num_lyr = num_lyr + 1
647 clouds%topIdx(num_lyr) = 1
648 clouds%ctt(num_lyr) = t(1)
652 if ((rh(kk-1)>=t_rh).and.(rh(kk)<t_rh))
then
653 if ((rh(kk+1)<t_rh).and.(in_cld==1))
then
654 clouds%baseIdx(num_lyr) = kk
661 kstart = cur_base + 1
663 do k = kstart, topok-1
664 if (cur_base<=k)
then
665 if ((rh(k-1)<t_rh).and.(in_cld==0))
then
666 if ((rh(k)>=t_rh).and.(rh(k+1)>=t_rh))
then
667 num_lyr = num_lyr + 1
670 clouds%topIdx(num_lyr) = k
671 clouds%ctt(num_lyr) = t(k)
674 do kk=clouds%topIdx(num_lyr),topok-1
675 if ((rh(kk-1)>=t_rh).and.(rh(kk)<t_rh))
then
676 if ((rh(kk+1)<t_rh).and.(in_cld==1))
then
677 clouds%baseIdx(num_lyr) = kk
691 clouds%baseIdx(num_lyr) = topok
694 clouds%nLayers = num_lyr
696 clouds%wmnIdx = getwarmnoseidx(t, nz)
701 clouds%avv = getaveragevertvel(t, vv, nz, &
702 clouds%topIdx(num_lyr), clouds%baseIdx(num_lyr))
708 call calc_layerq(t, rh, pres, ept, nz, clouds)
711 end subroutine calc_cloudlayers
726 integer function getwarmnoseidx(t, nz)
728 integer,
intent(in) :: nz
729 real,
intent(in) :: t(nz)
731 logical :: abovefreezing
736 abovefreezing = .false.
744 if(t(k) <= 273.15)
then
745 if (abovefreezing)
then
746 getwarmnoseidx = tmpindex
751 abovefreezing = .false.
753 if(.not. abovefreezing) tmpindex = k
754 abovefreezing = .true.
760 end function getwarmnoseidx
767 real function getaveragevertvel(t,vv,nz, topIdx_lowest,baseIdx_lowest)
769 integer,
intent(in) :: nz
770 real,
intent(in) :: t(nz), vv(nz)
771 integer,
intent(in) :: topidx_lowest,baseidx_lowest
773 integer :: numvertvel, k, start_base
781 if(baseidx_lowest == nz)
then
784 start_base = baseidx_lowest
787 do k = start_base, topidx_lowest, -1
788 if ( t(k) <= 260.15 .and. t(k) >= 257.15)
then
789 sumvertvel = sumvertvel + vv(k)
790 numvertvel = numvertvel + 1
791 else if (t(k) < 257.15)
then
792 sumvertvel = sumvertvel + vv(k) + vv(k+1)
798 if (numvertvel == 0)
then
799 getaveragevertvel = 0.0
801 getaveragevertvel = sumvertvel / numvertvel
805 end function getaveragevertvel
808 subroutine calc_layerq(t, rh, pres, ept, nz, clouds)
810 integer,
intent(in) :: nz
811 real,
intent(in) :: t(nz), rh(nz), pres(nz), ept(nz)
812 type(clouds_t
),
intent(inout) :: clouds
815 real :: mean_rh, te_map, rh_map, adj_q
816 real :: totalq, testthetae, q_top, q_base, q_layer, delta_te, sum_rh
817 integer :: base_k, test_k, num_layers
819 integer :: k, m, n, kk
821 clouds%layerQ(:) = 0.0
824 lp_nlayers:
do n = 1, clouds%nLayers
825 lp_k_inlayer:
do k = clouds%topIdx(n), clouds%baseIdx(n)-1
831 base_k = clouds%baseIdx(n)
833 if((testthetae-ept(m)>4.) .and. (clouds%baseIdx(n)<=m))
then
840 lp_kk:
do kk = k, base_k-1
841 q_top = mixing_ratio(t(kk), pres(kk))
842 q_base = mixing_ratio(t(kk+1), pres(kk+1))
843 q_layer = q_base - q_top
846 q_layer = q_layer * (pres(kk)/(287.05 * t(kk)))
848 if (q_layer < 0.0) q_layer = 0.0
850 delta_te = testthetae - ept(kk+1)
857 num_layers = num_layers + 1
858 sum_rh = sum_rh + rh(m)
861 if (num_layers == 0)
then
864 mean_rh = sum_rh / num_layers
867 te_map = dq_delta_te_map(delta_te)
868 rh_map = dq_rh_map(mean_rh)
869 adj_q = q_layer * te_map * rh_map
871 totalq = totalq + adj_q
874 clouds%layerQ(k) = totalq
879 end subroutine calc_layerq
883 real function dq_rh_map(rh)
885 real,
intent(in) :: rh
889 else if( rh >= 100.)
then
892 dq_rh_map = (rh-70.)/30.
896 end function dq_rh_map
900 real function dq_delta_te_map(delTE)
902 real,
intent(in) :: delte
904 if (delte <= 0.0)
then
905 dq_delta_te_map = 1.0
906 else if (delte <= 4.0)
then
907 dq_delta_te_map = (4. - delte) / 4.0
909 dq_delta_te_map = 0.0
913 end function dq_delta_te_map
915 end module cloudlayers
921 module icingpotential
923 use ctlblk_mod
, only: me
925 use derivedfields
, only : precips
926 use cloudlayers
, only : clouds_t
936 subroutine icing_pot(hgt, rh, t, liqCond, vv, nz, clouds, ice_pot)
938 integer,
intent(in) :: nz
939 real,
intent(in) :: hgt(nz), rh(nz), t(nz), liqcond(nz), vv(nz)
940 type(clouds_t
),
intent(in) :: clouds
941 real,
intent(out) :: ice_pot(nz)
943 real :: ctt, slw, slw_fac, vv_fac
947 real :: mapt,maprh,mapctt,mapvv,mapslw
952 lp_n:
do n = 1, clouds%nLayers
957 lp_k:
do k = clouds%topIdx(n), clouds%baseIdx(n)
960 if ( ctt>=259.15 .and. ctt<=273.15 )
then
966 ice_pot(k) = tmap(t(k)) * rh_map(rh(k)) * ctt_map(ctt)
969 if (vv_map(vv(k))>=0.0)
then
970 vv_fac = (1.0-ice_pot(k)) * (0.25*vv_map(vv(k)))
972 vv_fac = ice_pot(k) * (0.25*vv_map(vv(k)))
976 slw_fac = (1.0 - ice_pot(k))*(0.4*slw_map(slw))
979 if (ice_pot(k)>0.001)
then
980 ice_pot(k) = ice_pot(k) + vv_fac + slw_fac
984 if (ice_pot(k)>1.0)
then
1005 end subroutine icing_pot
1008 real function tmap(temp)
1010 real,
intent(in) :: temp
1012 if((temp>=248.15).and.(temp<=263.15))
then
1013 tmap=((temp-248.15)/15.0)**1.5
1014 elseif((temp>263.15).and.(temp<268.15))
then
1016 elseif((temp>268.15).and.(temp<273.15))
then
1017 tmap=((273.15 - temp)/5.0)**0.75
1027 real function rh_map(rh)
1029 real,
intent(in) :: rh
1033 elseif ( (rh<=95.0).and.(rh>=50.0) )
then
1034 rh_map=((20./9.) * ((rh/100.0)-0.5))**3.0
1044 real function ctt_map(ctt)
1046 real,
intent(in) :: ctt
1048 if((ctt>=261.0).and.(ctt<280.))
then
1050 elseif((ctt>223.0).and.(ctt<261.0 ))
then
1051 ctt_map = 0.2 + 0.8*(((ctt - 223.0)/38.0)**2.0)
1052 elseif(ctt<223.0)
then
1059 end function ctt_map
1064 real function vv_map(vv)
1066 real,
intent(in) :: vv
1070 elseif (vv<-0.5)
then
1073 vv_map = -1.0 * (vv/0.5)
1082 real function slw_map(slw)
1084 real,
intent(in) :: slw
1088 elseif (slw<=0.001)
then
1095 end function slw_map
1097 end module icingpotential
1109 integer :: no_precipitaion = 0
1110 integer :: precipitaion_below_warmnose = 1
1111 integer :: precipitaion_above_warmnose = 2
1112 integer :: all_snow = 3
1113 integer :: cold_rain = 4
1114 integer :: warm_precipitaion = 5
1115 integer :: freezing_precipitaion = 6
1116 integer :: convection = 7
1117 end type scenarios_t
1118 type(scenarios_t
),
parameter :: scenarios = scenarios_t()
1126 real function twp_map(v, scenario)
1128 real,
intent(in) :: v
1129 integer,
intent(in) :: scenario
1131 select case (scenario)
1132 case(scenarios%PRECIPITAION_ABOVE_WARMNOSE, scenarios%ALL_SNOW, &
1133 scenarios%COLD_RAIN, scenarios%FREEZING_PRECIPITAION)
1137 else if( v < 1000.)
then
1138 twp_map = v / 1000.0
1142 case( scenarios%WARM_PRECIPITAION)
1146 else if( v < 500.)
then
1153 end function twp_map
1156 real function t_map(v, scenario)
1158 real,
intent(in) :: v
1159 integer,
intent(in) :: scenario
1161 select case (scenario)
1162 case(scenarios%PRECIPITAION_BELOW_WARMNOSE)
1164 if(v <= 269.15)
then
1166 elseif( v <= 273.15)
then
1167 t_map = (273.15 - v ) / 4.0
1175 if(v <= 248.15)
then
1177 elseif( v <= 248.65)
then
1178 t_map = (49.162215 - 0.1961*v) * 2.
1179 elseif( v <= 249.15)
then
1180 t_map = (20.617195 - 0.0813*v) * 2.
1181 elseif(v <= 251.15)
then
1182 t_map = (52.02265 - 0.203*v) / 2.0
1183 elseif(v <= 253.15)
then
1184 t_map = (36.14997 - 0.1398*v) / 2.0
1185 elseif(v <= 255.15)
then
1186 t_map = (29.51744 - 0.1136*v) / 2.0
1187 elseif(v <= 257.15)
then
1188 t_map = (25.613645-0.0983*v) / 2.0
1189 elseif(v <= 259.15)
then
1190 t_map = (22.91357 - 0.0878*v) / 2.0
1191 elseif( v <= 261.15)
then
1192 t_map = (20.918115 - 0.0801*v) / 2.
1201 real function prcpcondensate_map(v, scenario)
1203 real,
intent(in) :: v
1204 integer,
intent(in) :: scenario
1205 prcpcondensate_map = 0.0
1206 select case (scenario)
1207 case(scenarios%NO_PRECIPITAION)
1209 case(scenarios%PRECIPITAION_BELOW_WARMNOSE, scenarios%COLD_RAIN, &
1210 scenarios%PRECIPITAION_ABOVE_WARMNOSE)
1213 prcpcondensate_map = 0.
1214 elseif(v <= 0.2)
then
1215 prcpcondensate_map = (v - 0.05) / 0.15
1217 prcpcondensate_map = 1.
1219 case(scenarios%ALL_SNOW)
1222 prcpcondensate_map = 0.
1223 elseif(v <= 0.25)
then
1224 prcpcondensate_map = (v - 0.05) / 0.2
1226 prcpcondensate_map = 1.0
1228 case(scenarios%WARM_PRECIPITAION, scenarios%FREEZING_PRECIPITAION)
1231 prcpcondensate_map = 0.
1232 elseif(v <= 0.15)
then
1233 prcpcondensate_map = (.5*v - 0.025) / 0.1
1235 prcpcondensate_map = 0.5
1239 end function prcpcondensate_map
1242 real function deltaz_map(v, scenario)
1244 real,
intent(in) :: v
1245 integer,
intent(in) :: scenario
1247 select case (scenario)
1248 case (scenarios%NO_PRECIPITAION, scenarios%WARM_PRECIPITAION)
1252 else if(v <= 1828.8)
then
1253 deltaz_map = v /1828.8
1257 case(scenarios%PRECIPITAION_BELOW_WARMNOSE)
1264 case(scenarios%ALL_SNOW)
1268 else if(v <= 2743.2)
then
1269 deltaz_map = (v - 914.4) / 1828.8
1273 case(scenarios%COLD_RAIN, scenarios%FREEZING_PRECIPITAION)
1277 else if(v <= 4267.2)
then
1278 deltaz_map = (v - 1524.) / 2743.2
1284 end function deltaz_map
1292 real function ctt_map(v)
1294 real,
intent(in) :: v
1295 if(v <= 223.15)
then
1297 elseif(v <= 233.15)
then
1298 ctt_map = (20.36251 - 0.0554*v) / 10.
1299 elseif(v <= 243.15)
then
1300 ctt_map = (46.19553 - 0.1662*v) / 10.
1301 elseif(v <= 253.15)
then
1302 ctt_map = (73.13655 - 0.277*v) / 10.
1303 elseif(v <= 261.15)
then
1304 ctt_map = (78.71061 - 0.3014*v) / 8.
1308 end function ctt_map
1311 real function vv_map(v)
1313 real,
intent(in) :: v
1316 else if(v <= 0.)
then
1326 real function cldtopdist_map(v)
1328 real,
intent(in) :: v
1329 if( v <= 609.6)
then
1330 cldtopdist_map = 1.0
1331 elseif( v <= 3048.)
then
1332 cldtopdist_map = (3048. - v) / 2438.4
1336 end function cldtopdist_map
1341 real function cldbasedist_map(v)
1343 real,
intent(in) :: v
1344 if( v <= 304.8)
then
1345 cldbasedist_map = 1.0
1346 elseif( v <= 1524.)
then
1347 cldbasedist_map = (1524. - v) / 1219.2
1349 cldbasedist_map = 0.
1351 end function cldbasedist_map
1354 real function deltaq_map(v)
1356 real,
intent(in) :: v
1359 elseif( v <= 1.0)
then
1364 end function deltaq_map
1366 real function moisture_map_cond(rh, liqCond, iceCond, pres, t)
1368 real,
intent(in) :: rh, liqcond, icecond, pres, t
1369 real :: liq_m3, ice_m3, rhfactor, liqfactor,icefactor
1371 liq_m3 = (liqcond * pres) / (287.05 * t)
1372 ice_m3 = (icecond * pres) / (287.05 * t)
1373 rhfactor = 0.5 * rh_map(rh)
1374 liqfactor = 0.6* condensate_map(liq_m3)
1375 icefactor = 0.4 * condensate_map(ice_m3)
1376 moisture_map_cond = rhfactor + (1.0 - rhfactor) * (liqfactor + icefactor)
1378 end function moisture_map_cond
1381 real function moisture_map_cwat(rh, cwat, pres, t)
1383 real,
intent(in) :: rh, cwat, pres, t
1384 real :: cwat_m3, rhfactor, cwatfactor
1386 cwat_m3 = (cwat * pres) / (287.05 * t)
1387 rhfactor = 0.5 * rh_map(rh)
1388 cwatfactor = condensate_map(cwat_m3)
1389 moisture_map_cwat = rhfactor + (1.0 - rhfactor) * cwatfactor
1391 end function moisture_map_cwat
1395 real function rh_map(v)
1397 real,
intent(in) :: v
1400 else if(v <= 100)
then
1401 rh_map = (v - 70.) / 30.
1409 real function condensate_map(v)
1411 real,
intent(in) :: v
1414 elseif( v <= 0.2)
then
1415 condensate_map = (v - 0.004) / 0.196
1419 end function condensate_map
1428 real function convect_t_map(v)
1430 real,
intent(in) :: v
1431 if(v <= 243.15)
then
1433 else if(v <= 265.15)
then
1434 convect_t_map = (v -243.15)/22.
1435 else if(v <= 269.15)
then
1437 else if(v <= 270.15)
then
1438 convect_t_map = -0.13*v + 35.9895
1439 else if(v <= 271.15)
then
1440 convect_t_map = -0.16*v + 44.094
1441 else if(v <= 272.15)
then
1442 convect_t_map = -0.21*v + 57.6515
1443 else if(v <= 273.15)
then
1444 convect_t_map = -0.50*v + 136.575
1448 end function convect_t_map
1452 real function convect_qpf_map(v)
1454 real,
intent(in) :: v
1457 else if(v <= 3.0)
then
1458 convect_qpf_map = (v - 1.) / 2.
1460 convect_qpf_map = 1.
1462 end function convect_qpf_map
1465 real function convect_cape_map(v)
1467 real,
intent(in) :: v
1469 if (v <= 1000.0)
then
1470 convect_cape_map = 0.0
1471 else if(v <= 2500.)
then
1472 convect_cape_map = (v - 1000.0) / 1400.
1474 convect_cape_map = 1.0
1476 end function convect_cape_map
1480 real function convect_liftedidx_map(v)
1482 real,
intent(in) :: v
1484 convect_liftedidx_map = 1.0
1485 else if(v <= 0.)
then
1486 convect_liftedidx_map = -0.1 * v
1488 convect_liftedidx_map = 0.
1490 end function convect_liftedidx_map
1494 real function convect_kidx_map(v)
1496 real,
intent(in) :: v
1498 convect_kidx_map = 0
1499 else if(v <= 40.0)
then
1500 convect_kidx_map = (v - 20.) / 20.
1502 convect_kidx_map = 1.
1504 end function convect_kidx_map
1507 real function convect_totals_map(v)
1509 real,
intent(in) :: v
1511 convect_totals_map = 0
1512 else if( v <= 55.0)
then
1513 convect_totals_map = (v - 20) / 35.
1515 convect_totals_map = 1.0
1517 end function convect_totals_map
1519 end module severitymaps
1524 module icingseverity
1525 use derivedfields
, only : precips
1526 use cloudlayers
, only : clouds_t
1534 real,
parameter :: cold_prcp_t_thld = 261.15
1538 real :: ctt, layerq, avv
1539 integer :: nlayers, topidx, baseidx, wmnidx
1541 real :: deltaz, cloudtopdist
1543 logical :: lowestcloud
1551 subroutine icing_sev(imp_physics,hgt, rh, t, pres, vv, liqCond, iceCond, twp, &
1552 ice_pot, nz, hcprcp, cape, lx, kx, tott, pc, prcptype, clouds, &
1555 integer,
intent(in) :: imp_physics
1556 integer,
intent(in) :: nz
1557 real,
intent(in) :: hgt(nz), rh(nz), t(nz), pres(nz), vv(nz)
1558 real,
intent(in) :: liqcond(nz), icecond(nz), twp(nz), ice_pot(nz)
1559 real,
intent(in) :: hcprcp, cape, lx, kx, tott, pc
1560 integer,
intent(in) :: prcptype
1561 type(clouds_t
),
intent(in) :: clouds
1562 real,
intent(out) :: iseverity(nz)
1571 lowestcloud = .false.
1574 lp_n:
do n = 1, clouds%nLayers
1577 wmnidx = clouds%wmnIdx
1579 if (n == clouds%nLayers) lowestcloud = .true.
1582 topidx = clouds%topIdx(n)
1583 baseidx = clouds%baseIdx(n)
1585 lp_k:
do k = topidx, baseidx
1586 if (ice_pot(k) < 0.001) cycle
1591 layerq = clouds%layerQ(k)
1593 deltaz = hgt(k) - hgt(baseidx)
1594 cloudtopdist = hgt(topidx) - hgt(k)
1597 if(imp_physics == 98 .or. imp_physics == 99)
then
1598 moistint = moisture_map_cwat(rh(k), liqcond(k)+icecond(k), pres(k), t(k))
1600 moistint = moisture_map_cond(rh(k), liqcond(k),icecond(k), pres(k), t(k))
1603 if(isconvection(prcptype))
then
1604 call convectionscenario &
1605 (t(k), hcprcp, cape, lx, kx, tott, severity)
1606 elseif(isclassicprcpblwwmn(prcptype, k))
then
1607 call classicprcpblwwmnscenario &
1608 (t(k), avv, pc, ice_pot(k), severity)
1609 elseif(isclassicprcpabvwmn(prcptype, k))
then
1610 call classicprcpabvwmnscenario &
1611 (twp(k), vv(k), t(k), pc, ice_pot(k),severity)
1612 elseif(isnoprecip(prcptype))
then
1613 call noprcpscenario &
1614 (vv(k), t(k), pc, ice_pot(k), severity)
1615 elseif(issnow(prcptype))
then
1617 (twp(k), vv(k), t(k), pc, ice_pot(k), severity)
1618 elseif(iscoldrain(prcptype))
then
1619 call coldrainscenario &
1620 (twp(k), vv(k), t(k), pc, ice_pot(k), severity)
1621 elseif(iswarmprecip(prcptype))
then
1622 call warmprecipscenario &
1623 (twp(k), vv(k), t(k), pc, ice_pot(k), severity)
1624 elseif(isfreezingprecip(prcptype))
then
1625 call freezingprecipscenario &
1626 (twp(k), vv(k), t(k), pc, ice_pot(k), severity)
1662 iseverity(k)=min(1., severity)
1663 iseverity(k)=max(0., severity)
1669 end subroutine icing_sev
1672 logical function isconvection(prcpType)
1674 integer,
intent(in) :: prcptype
1676 isconvection = prcptype == precips%CONVECTION
1679 end function isconvection
1682 subroutine convectionscenario(t, hcprcp, cape, lx, kx, tott, severity)
1684 real,
intent(in) :: t, hcprcp, cape, lx, kx, tott
1685 real,
intent(out) :: severity
1689 real :: tint, prcpint, capeint, lxint, kxint,tottint
1691 real :: weights(6) = (/ 4.0, 3.0, 5.0, 5.0, 2.0, 2.0 /)
1693 sumweight = sum(weights)
1695 tint = convect_t_map(t) ** 0.5
1697 prcpint = convect_qpf_map(hcprcp)
1698 capeint = convect_cape_map(cape)
1699 lxint = convect_liftedidx_map(lx)
1700 kxint = convect_kidx_map(kx)
1701 tottint = convect_totals_map(tott)
1703 scenario = scenarios%CONVECTION
1705 severity = (weights(1) * tint + weights(2) * prcpint + &
1706 weights(3) * capeint + weights(4) * lxint + &
1707 weights(5) * kxint + weights(6) * tottint) / &
1710 end subroutine convectionscenario
1714 logical function isclassicprcpblwwmn(prcpType, k)
1716 integer,
intent(in) :: prcptype, k
1719 if ((prcptype /= precips%NONE .and. prcptype /= precips%SNOW) .and.&
1720 (wmnidx > 0 .and. k <= wmnidx .and. topidx > wmnidx) .and. &
1721 (ctt < cold_prcp_t_thld))
then
1722 isclassicprcpblwwmn = .true.
1724 isclassicprcpblwwmn = .false.
1728 end function isclassicprcpblwwmn
1731 subroutine classicprcpblwwmnscenario(t, vv, pc, ice_pot, severity)
1733 real,
intent(in) :: t, pc, vv, ice_pot
1734 real,
intent(out) :: severity
1738 real :: tint, pcint, dzint
1739 real :: cttint, vvint
1743 real :: weights(7) = (/ 3.0, 4.0, 3.0, 3.0, 3.5, 2.5, 2.0 /)
1745 sumweight = sum(weights)
1747 scenario = scenarios%PRECIPITAION_BELOW_WARMNOSE
1749 tint = t_map(t, scenario)
1750 pcint = prcpcondensate_map(pc, scenario)
1751 dzint = deltaz_map(deltaz, scenario)
1752 cttint = ctt_map(ctt)
1755 severity = (weights(1) * tint + weights(2) * pcint + &
1756 weights(3) * dzint + weights(4) * cttint + &
1757 weights(5) * moistint + weights(6) * vvint + &
1758 weights(7) * ice_pot) / sumweight
1761 end subroutine classicprcpblwwmnscenario
1766 logical function isclassicprcpabvwmn(prcpType, k)
1768 integer,
intent(in) :: prcptype, k
1771 if((prcptype /= precips%NONE .and. prcptype /= precips%SNOW) .and.&
1772 (wmnidx > 0 .and. k > wmnidx .and. topidx > wmnidx) .and. &
1773 (ctt < cold_prcp_t_thld))
then
1774 isclassicprcpabvwmn = .true.
1776 isclassicprcpabvwmn = .false.
1780 end function isclassicprcpabvwmn
1784 subroutine classicprcpabvwmnscenario(twp, vv, t, pc, ice_pot, severity)
1786 real,
intent(in) :: twp, vv, t, pc, ice_pot
1787 real,
intent(out) :: severity
1790 real :: twpint, vvint
1791 real :: tempadj, cttadj, pcadj
1795 real :: weights(4) = (/ 3.0, 3.5, 4.0, 2.0 /)
1797 sumweight = sum(weights)
1799 scenario = scenarios%PRECIPITAION_ABOVE_WARMNOSE
1801 twpint = twp_map(twp, scenario)
1808 call cal_dampingfactors(scenario, t, pc, tempadj, cttadj, pcadj)
1810 severity = (weights(1) * twpint + weights(2) * vvint + &
1811 weights(3) * moistint + weights(4) * ice_pot) / &
1812 (sumweight + 6.5*tempadj + 3.0*cttadj + 3.0*pcadj)
1815 end subroutine classicprcpabvwmnscenario
1820 logical function isnoprecip(prcpType)
1822 integer,
intent(in) :: prcptype
1824 isnoprecip = prcptype == precips%NONE
1827 end function isnoprecip
1832 subroutine noprcpscenario(vv, t, pc, ice_pot, severity)
1834 real,
intent(in) :: vv, t, pc, ice_pot
1835 real,
intent(out) :: severity
1839 real :: dqint, dzint, vvint
1840 real :: tempadj, cttadj, pcadj
1844 real :: weights(5) = (/ 3.0, 3.5, 4.0, 4.0, 3.0 /)
1846 sumweight = sum(weights)
1848 scenario = scenarios%NO_PRECIPITAION
1850 dqint = deltaq_map(layerq)
1851 dzint = deltaz_map(deltaz, scenario)
1858 call cal_dampingfactors(scenario, t, pc, tempadj, cttadj, pcadj)
1860 severity = (weights(1) * dqint + weights(2) * dzint + &
1861 weights(3) * vvint + weights(4) * moistint + &
1862 weights(5) * ice_pot) / &
1863 (sumweight + 9.0*tempadj + 4.5*cttadj)
1866 end subroutine noprcpscenario
1871 logical function issnow(prcpType)
1873 integer,
intent(in) :: prcptype
1875 issnow = prcptype == precips%SNOW
1882 subroutine snowscenario(twp, vv, t, pc, ice_pot, severity)
1884 real,
intent(in) :: twp, vv, t, pc, ice_pot
1885 real,
intent(out) :: severity
1889 real :: twpint, dzint, vvint
1890 real :: tempadj, cttadj, pcadj
1894 real :: weights(5) = (/ 3.0, 3.5, 3.5, 4.0, 3.0 /)
1896 sumweight = sum(weights)
1898 scenario = scenarios%ALL_SNOW
1900 twpint = twp_map(twp, scenario)
1901 dzint = deltaz_map(deltaz, scenario)
1908 call cal_dampingfactors(scenario, t, pc, tempadj, cttadj, pcadj)
1910 severity = (weights(1) * twpint + weights(2) * dzint + &
1911 weights(3) * vvint + weights(4) * moistint + &
1912 weights(5) * ice_pot) / &
1913 (sumweight + 9.0*tempadj + 4.0*cttadj + 4.0*pcadj)
1916 end subroutine snowscenario
1920 logical function iscoldrain(prcpType)
1922 integer,
intent(in) :: prcptype
1925 iscoldrain = prcptype == precips%RAIN .and. ctt < cold_prcp_t_thld
1928 end function iscoldrain
1932 subroutine coldrainscenario(twp, vv, t, pc, ice_pot, severity)
1934 real,
intent(in) :: twp, vv, t, pc, ice_pot
1935 real,
intent(out) :: severity
1939 real :: twpint, dzint, vvint
1940 real :: tempadj, cttadj, pcadj
1944 real :: weights(5) = (/ 3.0, 3.5, 3.5, 4.0, 2.0 /)
1946 sumweight = sum(weights)
1948 scenario = scenarios%COLD_RAIN
1950 twpint = twp_map(twp, scenario)
1951 dzint = deltaz_map(deltaz, scenario)
1958 call cal_dampingfactors(scenario, t, pc, tempadj, cttadj, pcadj)
1960 severity = (weights(1) * twpint + weights(2) * dzint + &
1961 weights(3) * vvint + weights(4) * moistint + &
1962 weights(5) * ice_pot) / &
1963 (sumweight + 8.0*tempadj + 4.0*cttadj + 4.0*pcadj)
1966 end subroutine coldrainscenario
1974 logical function iswarmprecip(prcpType)
1976 integer,
intent(in) :: prcptype
1979 if((prcptype == precips%RAIN .or. prcptype == precips%OTHER) .and. &
1980 (ctt >= cold_prcp_t_thld))
then
1981 iswarmprecip = .true.
1982 elseif((prcptype /= precips%NONE .and. prcptype /= precips%SNOW).and.&
1983 (wmnidx > 0 .and. topidx <= wmnidx) .and. &
1984 (ctt >= cold_prcp_t_thld))
then
1985 iswarmprecip = .true.
1987 iswarmprecip = .false.
1991 end function iswarmprecip
1995 subroutine warmprecipscenario(twp, vv, t, pc, ice_pot, severity)
1997 real,
intent(in) :: twp, vv, t, pc, ice_pot
1998 real,
intent(out) :: severity
2002 real :: twpint, dzint, vvint, dqint
2003 real :: tempadj, cttadj, pcadj
2007 real :: weights(6) = (/ 3.0, 3.0, 3.0, 3.5, 4.0, 2.0 /)
2009 sumweight = sum(weights)
2011 scenario = scenarios%WARM_PRECIPITAION
2013 twpint = twp_map(twp, scenario)
2014 dzint = deltaz_map(deltaz, scenario)
2016 dqint = deltaq_map(layerq)
2022 call cal_dampingfactors(scenario, t, pc, tempadj, cttadj, pcadj)
2024 severity = (weights(1) * twpint + weights(2) * dzint + &
2025 weights(3) * vvint + weights(4) * moistint + &
2026 weights(5) * ice_pot + weights(6) * dqint) / &
2027 (sumweight + 9.0*tempadj + 4.5*cttadj + 4.5*pcadj)
2030 end subroutine warmprecipscenario
2034 logical function isfreezingprecip(prcpType)
2036 integer,
intent(in) :: prcptype
2039 if (prcptype == precips%OTHER .and. ctt < cold_prcp_t_thld)
then
2040 isfreezingprecip = .true.
2042 isfreezingprecip = .false.
2046 end function isfreezingprecip
2051 subroutine freezingprecipscenario(twp, vv, t, pc, ice_pot, severity)
2053 real,
intent(in) :: twp, vv, t, pc, ice_pot
2054 real,
intent(out) :: severity
2058 real :: twpint, dzint, vvint
2059 real :: tempadj, cttadj, pcadj
2063 real :: weights(5) = (/ 3.0, 3.5, 3.5, 4.0, 2.0 /)
2065 sumweight = sum(weights)
2067 scenario = scenarios%FREEZING_PRECIPITAION
2069 twpint = twp_map(twp, scenario)
2070 dzint = deltaz_map(deltaz, scenario)
2077 call cal_dampingfactors(scenario, t, pc, tempadj, cttadj, pcadj)
2079 severity = (weights(1) * twpint + weights(2) * dzint + &
2080 weights(3) * vvint + weights(4) * moistint + &
2081 weights(5) * ice_pot) / &
2082 (sumweight + 8.0*tempadj + 4.0*cttadj + 4.0*pcadj)
2085 end subroutine freezingprecipscenario
2103 subroutine cal_dampingfactors(scenario, t, pc, tempAdj,cttAdj,condAdj)
2105 integer,
intent(in) :: scenario
2106 real,
intent(in) :: t, pc
2107 real,
intent(out) :: tempadj, cttadj, condadj
2110 tempadj = t_map(t, scenario)
2112 cttadj = ctt_map(ctt) * cldtopdist_map(cloudtopdist)
2114 if (lowestcloud)
then
2115 condadj = cldbasedist_map(deltaz)
2117 condadj = condadj * prcpcondensate_map(pc, scenario)
2121 end subroutine cal_dampingfactors
2124 end module icingseverity
2129 subroutine icing_algo(i,j,pres,temp,rh,hgt,omega,wh,&
2130 q,cwat,qqw,qqi,qqr,qqs,qqg,&
2131 nz,xlat,xlon,xalt,prate,cprate,&
2132 cape,cin, ice_pot, ice_sev)
2133 use ctlblk_mod
, only: imp_physics, spval, dtq2,me
2134 use physcons_post,only: g => con_g, fv => con_fvirt, rd => con_rd
2136 use derivedfields
, only : derive_fields
2137 use cloudlayers
, only : calc_cloudlayers, clouds_t
2138 use icingpotential
, only : icing_pot
2139 use icingseverity
, only : icing_sev
2143 integer,
external :: gettopok
2178 integer,
intent(in) :: i,j, nz
2179 real,
intent(in),
dimension(nz) :: pres,temp,rh,hgt,omega,wh
2181 real,
intent(in),
dimension(nz) :: q,cwat,qqw,qqi,qqr,qqs,qqg
2182 real,
intent(in) :: xlat, xlon, xalt
2183 real,
intent(in) :: prate, cprate
2184 real,
intent(in) :: cape, cin
2185 real,
intent(out) :: ice_pot(nz), ice_sev(nz)
2187 real :: hprcp, hcprcp
2188 integer :: topok, region, prcptype
2189 real,
allocatable :: vv(:), ept(:), wbt(:), twp(:)
2190 real,
allocatable :: icecond(:),liqcond(:), totwater(:),totcond(:)
2191 real :: pc, kx, lx, tott
2192 type(clouds_t
) :: clouds
2202 allocate(icecond(nz),liqcond(nz))
2203 allocate(totwater(nz),totcond(nz))
2205 allocate(clouds%layerQ(nz))
2219 topok = gettopok(hgt,xalt,nz)
2222 hprcp = prate * 1000. / dtq2 * 3600.
2223 hcprcp = cprate * 1000. / dtq2 * 3600.
2225 if(imp_physics == 98 .or. imp_physics == 99)
then
2231 if(cwat(k) == spval)
then
2238 if(temp(k) >=259.15)
then
2239 liqcond(k) = cwat(k) * 1000.
2243 icecond(k) = cwat(k) * 1000.
2246 totwater(k) = cwat(k) * 1000.
2248 totcond(k) = cwat(k) * 1000.
2250 else if(imp_physics == 11 .or. imp_physics == 8)
then
2255 if(vv(k) /= spval)
then
2256 tv = temp(k)*(1.+q(k)*fv)
2257 vv(k)=-vv(k)*g*(pres(k)/(rd*tv))
2260 if(qqw(k) /= spval .and. qqr(k) /= spval)
then
2261 liqcond(k) = (qqw(k) + qqr(k)) * 1000.
2266 if(qqi(k) /= spval .and. qqs(k) /= spval .and. qqg(k) /= spval)
then
2267 icecond(k) = (qqi(k) + qqs(k) + qqg(k)) * 1000.
2272 if(liqcond(k) /= spval .and. icecond(k) /= spval)
then
2273 totwater(k) = liqcond(k) + icecond(k)
2278 if(qqr(k) /= spval .and. qqs(k) /= spval .and. qqg(k) /= spval)
then
2279 totcond(k) = (qqr(k) + qqs(k) + qqg(k)) * 1000.
2285 print *,
"This schema is not supported. imp_physics=", imp_physics
2289 call derive_fields(imp_physics,temp, rh, pres, hgt, totwater,totcond, &
2290 nz, topok, hprcp, hcprcp, cin, cape,&
2291 ept, wbt, twp, pc, kx, lx, tott, prcptype)
2293 call calc_cloudlayers(rh, temp, pres, ept, vv, nz,topok, xlat,xlon,&
2301 call icing_pot(hgt, rh, temp, liqcond, vv, nz, clouds, ice_pot)
2304 call icing_sev(imp_physics, hgt, rh, temp, pres, vv, liqcond, icecond, twp, &
2305 ice_pot, nz, hcprcp, cape, lx, kx, tott, pc, prcptype, clouds, &
2325 deallocate(icecond,liqcond)
2326 deallocate(totwater,totcond)
2327 deallocate(clouds%layerQ)
2330 end subroutine icing_algo
2333 integer function gettopok(hgt, alt, nz)
2335 real,
intent(in) :: hgt(nz)
2336 real,
intent(in) :: alt
2337 integer,
intent(in) :: nz
2341 if(hgt(nz) >= alt)
then
2345 if ((hgt(k-1) > alt) .and. (hgt(k) <= alt))
then
2353 end function gettopok