CCPP SciDoc v7.0.x  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
sfc_diff.f
1
5
7 module sfc_diff
8
9 use machine , only : kind_phys
10
11 implicit none
12
13 public :: sfc_diff_run
14 public :: stability
15
16 private
17
18 real (kind=kind_phys), parameter :: ca=0.4_kind_phys ! ca - von karman constant
19
20 contains
21
57 subroutine sfc_diff_run (im,rvrdm1,eps,epsm1,grav, & !intent(in)
58 & ps,t1,q1,z1,garea,wind, & !intent(in)
59 & prsl1,prslki,prsik1,prslk1, & !intent(in)
60 & sigmaf,vegtype,shdmax,ivegsrc, & !intent(in)
61 & z0pert,ztpert, & ! mg, sfc-perts !intent(in)
62 & flag_iter,redrag, & !intent(in)
63 & flag_lakefreeze,lakefrac,fice, & !intent(in)
64 & u10m,v10m,sfc_z0_type, & !hafs,z0 type !intent(in)
65 & u1,v1,usfco,vsfco,use_oceanuv, &
66 & wet,dry,icy, & !intent(in)
67 & thsfc_loc, & !intent(in)
68 & tskin_wat, tskin_lnd, tskin_ice, & !intent(in)
69 & tsurf_wat, tsurf_lnd, tsurf_ice, & !intent(in)
70 & z0rl_wat, z0rl_lnd, z0rl_ice, & !intent(inout)
71 & z0rl_wav, & !intent(inout)
72 & ustar_wat, ustar_lnd, ustar_ice, & !intent(inout)
73 & cm_wat, cm_lnd, cm_ice, & !intent(inout)
74 & ch_wat, ch_lnd, ch_ice, & !intent(inout)
75 & rb_wat, rb_lnd, rb_ice, & !intent(inout)
76 & stress_wat,stress_lnd,stress_ice, & !intent(inout)
77 & fm_wat, fm_lnd, fm_ice, & !intent(inout)
78 & fh_wat, fh_lnd, fh_ice, & !intent(inout)
79 & fm10_wat, fm10_lnd, fm10_ice, & !intent(inout)
80 & fh2_wat, fh2_lnd, fh2_ice, & !intent(inout)
81 & ztmax_wat, ztmax_lnd, ztmax_ice, & !intent(inout)
82 & zvfun, & !intent(out)
83 & errmsg, errflg) !intent(out)
84!
85 implicit none
86!
87 integer, parameter :: kp = kind_phys
88 integer, intent(in) :: im, ivegsrc
89 integer, intent(in) :: sfc_z0_type ! option for calculating surface roughness length over ocean
90 logical, intent(in) :: use_oceanuv ! option for including ocean current in the computation of flux
91
92 integer, dimension(:), intent(in) :: vegtype
93
94 logical, intent(in) :: redrag ! reduced drag coeff. flag for high wind over sea (j.han)
95 logical, dimension(:), intent(in) :: flag_iter, dry, icy
96 logical, dimension(:), intent(in) :: flag_lakefreeze
97 logical, dimension(:), intent(inout) :: wet
98
99 logical, intent(in) :: thsfc_loc ! Flag for reference pressure in theta calculation
100
101 real(kind=kind_phys), dimension(:), intent(in) :: u10m,v10m
102 real(kind=kind_phys), dimension(:), intent(in) :: u1,v1
103 real(kind=kind_phys), dimension(:), intent(in) :: usfco,vsfco
104 real(kind=kind_phys), intent(in) :: rvrdm1, eps, epsm1, grav
105 real(kind=kind_phys), dimension(:), intent(in) :: &
106 & ps,t1,q1,z1,garea,prsl1,prslki,prsik1,prslk1, &
107 & wind,sigmaf,shdmax, &
108 & z0pert,ztpert ! mg, sfc-perts
109 real(kind=kind_phys), dimension(:), intent(in) :: lakefrac
110 real(kind=kind_phys), dimension(:), intent(in) :: fice
111 real(kind=kind_phys), dimension(:), intent(in) :: &
112 & tskin_wat, tskin_lnd, tskin_ice, &
113 & tsurf_wat, tsurf_lnd, tsurf_ice
114
115 real(kind=kind_phys), dimension(:), intent(in) :: z0rl_wav
116 real(kind=kind_phys), dimension(:), intent(inout) :: &
117 & z0rl_wat, z0rl_lnd, z0rl_ice, &
118 & ustar_wat, ustar_lnd, ustar_ice, &
119 & cm_wat, cm_lnd, cm_ice, &
120 & ch_wat, ch_lnd, ch_ice, &
121 & rb_wat, rb_lnd, rb_ice, &
122 & stress_wat,stress_lnd,stress_ice, &
123 & fm_wat, fm_lnd, fm_ice, &
124 & fh_wat, fh_lnd, fh_ice, &
125 & fm10_wat, fm10_lnd, fm10_ice, &
126 & fh2_wat, fh2_lnd, fh2_ice, &
127 & ztmax_wat, ztmax_lnd, ztmax_ice
128 real(kind=kind_phys), dimension(:), intent(out) :: zvfun
129!
130 character(len=*), intent(out) :: errmsg
131 integer, intent(out) :: errflg
132!
133! locals
134!
135 integer i
136 real(kind=kind_phys) :: windrel
137!
138 real(kind=kind_phys) :: rat, tv1, thv1, restar, wind10m,
139 & czilc, tem1, tem2, virtfac
140!
141
142 real(kind=kind_phys) :: tvs, z0, z0max, ztmax, gdx
143!
144 real(kind=kind_phys), parameter :: z0lo=0.1, z0up=1.0
145!
146 real(kind=kind_phys), parameter ::
147 & one=1.0_kp, zero=0.0_kp, half=0.5_kp, qmin=1.0e-8_kp
148 &, charnock=.018_kp, z0s_max=.317e-2_kp &! a limiting value at high winds over sea
149 &, zmin=1.0e-6_kp &
150 &, vis=1.4e-5_kp, rnu=1.51e-5_kp, visi=one/vis &
151 &, log01=log(0.01_kp), log05=log(0.05_kp), log07=log(0.07_kp)
152
153! parameter (charnock=.014,ca=.4)!c ca is the von karman constant
154! parameter (alpha=5.,a0=-3.975,a1=12.32,b1=-7.755,b2=6.041)
155! parameter (a0p=-7.941,a1p=24.75,b1p=-8.705,b2p=7.899,vis=1.4e-5)
156
157! real(kind=kind_phys) aa1,bb1,bb2,cc,cc1,cc2,arnu
158! parameter (aa1=-1.076,bb1=.7045,cc1=-.05808)
159! parameter (bb2=-.1954,cc2=.009999)
160! parameter (arnu=.135*rnu)
161!
162! z0s_max=.196e-2 for u10_crit=25 m/s
163! z0s_max=.317e-2 for u10_crit=30 m/s
164! z0s_max=.479e-2 for u10_crit=35 m/s
165!
166! mbek -- toga-coare flux algorithm
167! parameter (rnu=1.51e-5,arnu=0.11*rnu)
168
169! Initialize CCPP error handling variables
170 errmsg = ''
171 errflg = 0
172
173! initialize variables. all units are supposedly m.k.s. unless specified
174! ps is in pascals, wind is wind speed,
175! surface roughness length is converted to m from cm
176!
177! write(0,*)'in sfc_diff, sfc_z0_type=',sfc_z0_type
178
179 do i=1,im
180 if(flag_iter(i) .or. flag_lakefreeze(i)) then
181
182 ! Need to initialize ztmax arrays
183 ztmax_lnd(i) = 1. ! log(1) = 0
184 ztmax_ice(i) = 1. ! log(1) = 0
185 ztmax_wat(i) = 1. ! log(1) = 0
186
187 virtfac = one + rvrdm1 * max(q1(i),qmin)
188
189 tv1 = t1(i) * virtfac ! Virtual temperature in middle of lowest layer
190 if(thsfc_loc) then ! Use local potential temperature
191 thv1 = t1(i) * prslki(i) * virtfac
192 else ! Use potential temperature reference to 1000 hPa
193 thv1 = t1(i) / prslk1(i) * virtfac
194 endif
195
196 zvfun(i) = zero
197 gdx = sqrt(garea(i))
198
199! compute stability dependent exchange coefficients
200! this portion of the code is presently suppressed
201!
202 if (dry(i)) then ! Some land
203
204 if(thsfc_loc) then ! Use local potential temperature
205 tvs = half * (tsurf_lnd(i)+tskin_lnd(i)) * virtfac
206 else ! Use potential temperature referenced to 1000 hPa
207 tvs = half * (tsurf_lnd(i)+tskin_lnd(i))/prsik1(i)
208 & * virtfac
209 endif
210
211 z0max = max(zmin, min(0.01_kp * z0rl_lnd(i), z1(i)))
212!** xubin's new z0 over land
213 tem1 = one - shdmax(i)
214 tem2 = tem1 * tem1
215 tem1 = one - tem2
216
217 if( ivegsrc == 1 ) then
218
219 if (vegtype(i) == 10) then
220 z0max = exp( tem2*log01 + tem1*log07 )
221 elseif (vegtype(i) == 6) then
222 z0max = exp( tem2*log01 + tem1*log05 )
223 elseif (vegtype(i) == 7) then
224! z0max = exp( tem2*log01 + tem1*log01 )
225 z0max = 0.01_kp
226 elseif (vegtype(i) == 16) then
227! z0max = exp( tem2*log01 + tem1*log01 )
228 z0max = 0.01_kp
229 else
230 z0max = exp( tem2*log01 + tem1*log(z0max) )
231 endif
232
233 elseif (ivegsrc == 2 ) then
234
235 if (vegtype(i) == 7) then
236 z0max = exp( tem2*log01 + tem1*log07 )
237 elseif (vegtype(i) == 8) then
238 z0max = exp( tem2*log01 + tem1*log05 )
239 elseif (vegtype(i) == 9) then
240! z0max = exp( tem2*log01 + tem1*log01 )
241 z0max = 0.01_kp
242 elseif (vegtype(i) == 11) then
243! z0max = exp( tem2*log01 + tem1*log01 )
244 z0max = 0.01_kp
245 else
246 z0max = exp( tem2*log01 + tem1*log(z0max) )
247 endif
248
249 endif
250! mg, sfc-perts: add surface perturbations to z0max over land
251 if (z0pert(i) /= zero ) then
252 z0max = z0max * (10.0_kp**z0pert(i))
253 endif
254
255 z0max = max(z0max, zmin)
256
257!! czilc = 10.0 ** (- (0.40/0.07) * z0) ! fei's canopy height dependance of czil
258! czilc = 0.8_kp
259!
260! tem1 = 1.0_kp - sigmaf(i)
261! ztmax_lnd(i) = z0max*exp( - tem1*tem1
262! & * czilc*ca*sqrt(ustar_lnd(i)*(0.01/1.5e-05)))
263!
264 czilc = 10.0_kp ** (- 4.0_kp * z0max) ! Trier et al. (2011,WAF)
265 czilc = max(min(czilc, 0.8_kp), 0.08_kp)
266 tem1 = 1.0_kp - sigmaf(i)
267 czilc = czilc * tem1 * tem1
268 ztmax_lnd(i) = z0max * exp( - czilc * ca
269 & * 258.2_kp * sqrt(ustar_lnd(i)*z0max) )
270!
271! mg, sfc-perts: add surface perturbations to ztmax/z0max ratio over land
272 if (ztpert(i) /= zero) then
273 ztmax_lnd(i) = ztmax_lnd(i) * (10.0_kp**ztpert(i))
274 endif
275 ztmax_lnd(i) = max(ztmax_lnd(i), zmin)
276!
277! compute a function of surface roughness & green vegetation fraction (zvfun)
278!
279 tem1 = (z0max - z0lo) / (z0up - z0lo)
280 tem1 = min(max(tem1, zero), 1.0_kp)
281 tem2 = max(sigmaf(i), 0.1_kp)
282 zvfun(i) = sqrt(tem1 * tem2)
283!
284 call stability
285! --- inputs:
286 & (z1(i), zvfun(i), gdx, tv1, thv1, wind(i),
287 & z0max, ztmax_lnd(i), tvs, grav, thsfc_loc,
288! --- outputs:
289 & rb_lnd(i), fm_lnd(i), fh_lnd(i), fm10_lnd(i), fh2_lnd(i),
290 & cm_lnd(i), ch_lnd(i), stress_lnd(i), ustar_lnd(i))
291 endif ! Dry points
292
293 if (icy(i)) then ! Some ice
294
295 zvfun(i) = zero
296
297 if(thsfc_loc) then ! Use local potential temperature
298 tvs = half * (tsurf_ice(i)+tskin_ice(i)) * virtfac
299 else ! Use potential temperature referenced to 1000 hPa
300 tvs = half * (tsurf_ice(i)+tskin_ice(i))/prsik1(i)
301 & * virtfac
302 endif
303
304 z0max = max(zmin, min(0.01_kp * z0rl_ice(i), z1(i)))
305!** xubin's new z0 over land and sea ice
306 tem1 = one - shdmax(i)
307 tem2 = tem1 * tem1
308 tem1 = one - tem2
309
310! Removed the following lines by W. Zheng, for effective z0m (z0max) is applied only
311! for land.
312!wz if( ivegsrc == 1 ) then
313!wz
314!wz z0max = exp( tem2*log01 + tem1*log(z0max) )
315!wz elseif (ivegsrc == 2 ) then
316!wz z0max = exp( tem2*log01 + tem1*log(z0max) )
317!wz endif
318
319 z0max = max(z0max, zmin)
320
321!! czilc = 10.0 ** (- (0.40/0.07) * z0) ! fei's canopy height
322! dependance of czil
323! czilc = 0.8_kp
324!
325! tem1 = 1.0_kp - sigmaf(i)
326! ztmax_ice(i) = z0max*exp( - tem1*tem1
327! & * czilc*ca*sqrt(ustar_ice(i)*(0.01/1.5e-05)))
328!
329 czilc = 10.0_kp ** (- 4.0_kp * z0max)
330 czilc = max(min(czilc, 0.8_kp), 0.08_kp)
331 tem1 = 1.0_kp - sigmaf(i)
332 czilc = czilc * tem1 * tem1
333 ztmax_ice(i) = z0max * exp( - czilc * ca
334 & * 258.2_kp * sqrt(ustar_ice(i)*z0max) )
335!
336 ztmax_ice(i) = max(ztmax_ice(i), 1.0e-6)
337!
338 call stability
339! --- inputs:
340 & (z1(i), zvfun(i), gdx, tv1, thv1, wind(i),
341 & z0max, ztmax_ice(i), tvs, grav, thsfc_loc,
342! --- outputs:
343 & rb_ice(i), fm_ice(i), fh_ice(i), fm10_ice(i), fh2_ice(i),
344 & cm_ice(i), ch_ice(i), stress_ice(i), ustar_ice(i))
345 endif ! Icy points
346
347! BWG: Everything from here to end of subroutine was after
348! the stuff now put into "stability"
349
350 if (wet(i)) then ! Some open ocean
351
352 zvfun(i) = zero
353
354 if(thsfc_loc) then ! Use local potential temperature
355 tvs = half * (tsurf_wat(i)+tskin_wat(i)) * virtfac
356 else
357 tvs = half * (tsurf_wat(i)+tskin_wat(i))/prsik1(i)
358 & * virtfac
359 endif
360
361 if (use_oceanuv) then
362 wind10m=sqrt((u10m(i)-usfco(i))**2+(v10m(i)-vsfco(i))**2)
363 windrel=sqrt((u1(i)-usfco(i))**2+(v1(i)-vsfco(i))**2)
364 else
365 wind10m=sqrt(u10m(i)*u10m(i)+v10m(i)*v10m(i))
366 windrel=wind(i)
367 endif
368
369 z0 = 0.01_kp * z0rl_wat(i)
370 z0max = max(zmin, min(z0,z1(i)))
371!
372!** test xubin's new z0
373
374! ztmax = z0max
375
376 restar = max(ustar_wat(i)*z0max*visi, 0.000001_kp)
377
378! restar = log(restar)
379! restar = min(restar,5.)
380! restar = max(restar,-5.)
381! rat = aa1 + (bb1 + cc1*restar) * restar
382! rat = rat / (1. + (bb2 + cc2*restar) * restar))
383! rat taken from zeng, zhao and dickinson 1997
384
385 rat = min(7.0_kp, 2.67_kp * sqrt(sqrt(restar)) - 2.57_kp)
386 ztmax_wat(i) = max(z0max * exp(-rat), zmin)
387!
388 if (sfc_z0_type == 6) then
389 call znot_t_v6(wind10m, ztmax_wat(i)) ! 10-m wind,m/s, ztmax(m)
390 else if (sfc_z0_type == 7) then
391 call znot_t_v7(wind10m, ztmax_wat(i)) ! 10-m wind,m/s, ztmax(m)
392 else if (sfc_z0_type > 0) then
393 write(0,*)'no option for sfc_z0_type=',sfc_z0_type
394 errflg = 1
395 errmsg = 'ERROR(sfc_diff_run): no option for sfc_z0_type'
396 return
397 endif
398!
399 call stability
400! --- inputs:
401 & (z1(i), zvfun(i), gdx, tv1, thv1, windrel,
402 & z0max, ztmax_wat(i), tvs, grav, thsfc_loc,
403! --- outputs:
404 & rb_wat(i), fm_wat(i), fh_wat(i), fm10_wat(i), fh2_wat(i),
405 & cm_wat(i), ch_wat(i), stress_wat(i), ustar_wat(i))
406!
407! update z0 over ocean
408!
409 if ((sfc_z0_type == -1) .and.
410 & (lakefrac(i) == 0.0 .and. fice(i) == 0.0) .and.
411 & (z0rl_wav(i)>1.0e-7_kp .and. z0rl_wav(i)<0.1_kp)) then
412 ! using wave model derived momentum roughness
413 tem1 = 0.11 * vis / ustar_wat(i)
414 z0 = tem1 + 0.01_kp * z0rl_wav(i)
415
416 if (redrag) then
417 z0rl_wat(i) = 100.0_kp * max(min(z0,z0s_max),1.0e-7_kp)
418 else
419 z0rl_wat(i) = 100.0_kp * max(min(z0,0.1_kp), 1.e-7_kp)
420 endif
421
422 elseif ((sfc_z0_type == 0) .or.
423 & ((sfc_z0_type == -1) .and.
424 & (z0rl_wav(i)<=1.0e-7_kp .or. z0rl_wav(i)>=0.1_kp))) then
425! z0 = (charnock / grav) * ustar_wat(i) * ustar_wat(i)
426 tem1 = 0.11 * vis / ustar_wat(i)
427 z0 = tem1 + (charnock/grav)*ustar_wat(i)*ustar_wat(i)
428
429
430! mbek -- toga-coare flux algorithm
431! z0 = (charnock / grav) * ustar(i)*ustar(i) + arnu/ustar(i)
432! new implementation of z0
433! cc = ustar(i) * z0 / rnu
434! pp = cc / (1. + cc)
435! ff = grav * arnu / (charnock * ustar(i) ** 3)
436! z0 = arnu / (ustar(i) * ff ** pp)
437
438 if (redrag) then
439 z0rl_wat(i) = 100.0_kp * max(min(z0,z0s_max),1.0e-7_kp)
440 else
441 z0rl_wat(i) = 100.0_kp * max(min(z0,0.1_kp), 1.e-7_kp)
442 endif
443
444 elseif (sfc_z0_type == 6) then ! wang
445 call znot_m_v6(wind10m, z0) ! wind, m/s, z0, m
446 z0rl_wat(i) = 100.0_kp * z0 ! cm
447 elseif (sfc_z0_type == 7) then ! wang
448 call znot_m_v7(wind10m, z0) ! wind, m/s, z0, m
449 z0rl_wat(i) = 100.0_kp * z0 ! cm
450 else
451 z0rl_wat(i) = 1.0e-4_kp
452 endif
453!
454 endif ! end of if(open ocean)
455!
456 endif ! end of if(flagiter) loop
457 enddo
458
459 end subroutine sfc_diff_run
460
461!----------------------------------------
463 subroutine stability &
464! --- inputs:
465 & ( z1, zvfun, gdx, tv1, thv1, wind, z0max, ztmax, tvs, grav, &
466 & thsfc_loc, &
467! --- outputs:
468 & rb, fm, fh, fm10, fh2, cm, ch, stress, ustar)
469!-----
470
471 integer, parameter :: kp = kind_phys
472! --- inputs:
473 real(kind=kind_phys), intent(in) :: &
474 & z1, zvfun, gdx, tv1, thv1, wind, z0max, ztmax, tvs, grav
475 logical, intent(in) :: thsfc_loc
476
477! --- outputs:
478 real(kind=kind_phys), intent(out) :: &
479 & rb, fm, fh, fm10, fh2, cm, ch, stress, ustar
480
481! --- locals:
482 real(kind=kind_phys), parameter :: alpha=5.0_kp, a0=-3.975_kp &
483 &, a1=12.32_kp, alpha4=4.0_kp*alpha &
484 &, b1=-7.755_kp, b2=6.041_kp &
485 &, xkrefsqr=0.3_kp, xkmin=0.05_kp &
486 &, xkgdx=3000.0_kp &
487 &, a0p=-7.941_kp, a1p=24.75_kp, b1p=-8.705_kp, b2p=7.899_kp&
488 &, zolmin=-10.0_kp, zero=0.0_kp, one=1.0_kp
489
490 real(kind=kind_phys) aa, aa0, bb, bb0, dtv, adtv,
491 & hl1, hl12, pm, ph, pm10, ph2,
492 & z1i,
493 & fms, fhs, hl0, hl0inf, hlinf,
494 & hl110, hlt, hltinf, olinf,
495 & tem1, tem2, zolmax
496
497 real(kind=kind_phys) xkzo
498
499 z1i = one / z1
500
501!
502! set background diffusivities with one for gdx >= xkgdx and
503! as a function of horizontal grid size for gdx < xkgdx
504! (i.e., gdx/xkgdx for gdx < xkgdx)
505!
506 if(gdx >= xkgdx) then
507 xkzo = one
508 else
509 xkzo = gdx / xkgdx
510 endif
511
512 tem1 = tv1 - tvs
513 if(tem1 > zero) then
514 tem2 = xkzo * zvfun
515 xkzo = min(max(tem2, xkmin), xkzo)
516 endif
517
518 zolmax = xkrefsqr / sqrt(xkzo)
519
520! compute stability indices (rb and hlinf)
521
522 dtv = thv1 - tvs
523 adtv = max(abs(dtv),0.001_kp)
524 dtv = sign(1.0_kp,dtv) * adtv
525
526 if(thsfc_loc) then ! Use local potential temperature
527 rb = max(-5000.0_kp, (grav+grav) * dtv * z1
528 & / ((thv1 + tvs) * wind * wind))
529 else ! Use potential temperature referenced to 1000 hPa
530 rb = max(-5000.0_kp, grav * dtv * z1
531 & / (tv1 * wind * wind))
532 endif
533
534 tem1 = one / z0max
535 tem2 = one / ztmax
536 fm = log((z0max+z1) * tem1)
537 fh = log((ztmax+z1) * tem2)
538 fm10 = log((z0max+10.0_kp) * tem1)
539 fh2 = log((ztmax+2.0_kp) * tem2)
540 hlinf = rb * fm * fm / fh
541 hlinf = min(max(hlinf,zolmin),zolmax)
542!
543! stable case
544!
545 if (dtv >= zero) then
546 hl1 = hlinf
547 if(hlinf > 0.25_kp) then
548 tem1 = hlinf * z1i
549 hl0inf = z0max * tem1
550 hltinf = ztmax * tem1
551 aa = sqrt(one + alpha4 * hlinf)
552 aa0 = sqrt(one + alpha4 * hl0inf)
553 bb = aa
554 bb0 = sqrt(one + alpha4 * hltinf)
555 pm = aa0 - aa + log( (aa + one)/(aa0 + one) )
556 ph = bb0 - bb + log( (bb + one)/(bb0 + one) )
557 fms = fm - pm
558 fhs = fh - ph
559 hl1 = fms * fms * rb / fhs
560 hl1 = min(hl1, zolmax)
561 endif
562!
563! second iteration
564!
565 tem1 = hl1 * z1i
566 hl0 = z0max * tem1
567 hlt = ztmax * tem1
568 aa = sqrt(one + alpha4 * hl1)
569 aa0 = sqrt(one + alpha4 * hl0)
570 bb = aa
571 bb0 = sqrt(one + alpha4 * hlt)
572 pm = aa0 - aa + log( (one+aa)/(one+aa0) )
573 ph = bb0 - bb + log( (one+bb)/(one+bb0) )
574 hl110 = hl1 * 10.0_kp * z1i
575 aa = sqrt(one + alpha4 * hl110)
576 pm10 = aa0 - aa + log( (one+aa)/(one+aa0) )
577 hl12 = (hl1+hl1) * z1i
578! aa = sqrt(one + alpha4 * hl12)
579 bb = sqrt(one + alpha4 * hl12)
580 ph2 = bb0 - bb + log( (one+bb)/(one+bb0) )
581!
582! unstable case - check for unphysical obukhov length
583!
584 else ! dtv < 0 case
585 olinf = z1 / hlinf
586 tem1 = 50.0_kp * z0max
587 if(abs(olinf) <= tem1) then
588 hlinf = -z1 / tem1
589 hlinf = max(hlinf, zolmin)
590 endif
591!
592! get pm and ph
593!
594 if (hlinf >= -0.5_kp) then
595 hl1 = hlinf
596 pm = (a0 + a1*hl1) * hl1 / (one+ (b1+b2*hl1) *hl1)
597 ph = (a0p + a1p*hl1) * hl1 / (one+ (b1p+b2p*hl1)*hl1)
598 hl110 = hl1 * 10.0_kp * z1i
599 pm10 = (a0 + a1*hl110) * hl110/(one+(b1+b2*hl110)*hl110)
600 hl12 = (hl1+hl1) * z1i
601 ph2 = (a0p + a1p*hl12) * hl12/(one+(b1p+b2p*hl12)*hl12)
602 else ! hlinf < 0.05
603 hl1 = -hlinf
604 tem1 = one / sqrt(hl1)
605 pm = log(hl1) + 2.0_kp * sqrt(tem1) - .8776_kp
606 ph = log(hl1) + 0.5_kp * tem1 + 1.386_kp
607! pm = log(hl1) + 2.0 * hl1 ** (-.25) - .8776
608! ph = log(hl1) + 0.5 * hl1 ** (-.5) + 1.386
609 hl110 = hl1 * 10.0_kp * z1i
610 pm10 = log(hl110) + 2.0_kp/sqrt(sqrt(hl110)) - 0.8776_kp
611! pm10 = log(hl110) + 2. * hl110 ** (-.25) - .8776
612 hl12 = (hl1+hl1) * z1i
613 ph2 = log(hl12) + 0.5_kp / sqrt(hl12) + 1.386_kp
614! ph2 = log(hl12) + .5 * hl12 ** (-.5) + 1.386
615 endif
616
617 endif ! end of if (dtv >= 0 ) then loop
618!
619! finish the exchange coefficient computation to provide fm and fh
620!
621 fm = fm - pm
622 fh = fh - ph
623 fm10 = fm10 - pm10
624 fh2 = fh2 - ph2
625 cm = ca * ca / (fm * fm)
626 ch = ca * ca / (fm * fh)
627 tem1 = 0.00001_kp/z1
628 cm = max(cm, tem1)
629 ch = max(ch, tem1)
630 stress = cm * wind * wind
631 ustar = sqrt(stress)
632
633!.................................
634 end subroutine stability
635!---------------------------------
636
637
640 SUBROUTINE znot_m_v6(uref, znotm)
641 use machine , only : kind_phys
642 IMPLICIT NONE
643! Calculate areodynamical roughness over water with input 10-m wind
644! For low-to-moderate winds, try to match the Cd-U10 relationship from COARE V3.5 (Edson et al. 2013)
645! For high winds, try to fit available observational data
646!
647! Bin Liu, NOAA/NCEP/EMC 2017
648!
649! uref(m/s) : wind speed at 10-m height
650! znotm(meter): areodynamical roughness scale over water
651!
652
653 REAL(kind=kind_phys), INTENT(IN) :: uref
654 REAL(kind=kind_phys), INTENT(OUT):: znotm
655 real(kind=kind_phys), parameter :: p13 = -1.296521881682694e-02,
656 & p12 = 2.855780863283819e-01, p11 = -1.597898515251717e+00,
657 & p10 = -8.396975715683501e+00,
658
659 & p25 = 3.790846746036765e-10, p24 = 3.281964357650687e-09,
660 & p23 = 1.962282433562894e-07, p22 = -1.240239171056262e-06,
661 & p21 = 1.739759082358234e-07, p20 = 2.147264020369413e-05,
662
663 & p35 = 1.840430200185075e-07, p34 = -2.793849676757154e-05,
664 & p33 = 1.735308193700643e-03, p32 = -6.139315534216305e-02,
665 & p31 = 1.255457892775006e+00, p30 = -1.663993561652530e+01,
666
667 & p40 = 4.579369142033410e-04
668
669
670 if (uref >= 0.0 .and. uref <= 6.5 ) then
671 znotm = exp(p10 + uref * (p11 + uref * (p12 + uref*p13)))
672 elseif (uref > 6.5 .and. uref <= 15.7) then
673 znotm = p20 + uref * (p21 + uref * (p22 + uref * (p23
674 & + uref * (p24 + uref * p25))))
675 elseif (uref > 15.7 .and. uref <= 53.0) then
676 znotm = exp( p30 + uref * (p31 + uref * (p32 + uref * (p33
677 & + uref * (p34 + uref * p35)))))
678 elseif ( uref > 53.0) then
679 znotm = p40
680 else
681 print*, 'Wrong input uref value:',uref
682 endif
683
684 END SUBROUTINE znot_m_v6
685
691!
692!! uref(m/s) : wind speed at 10-m height
693!! znott(meter): scalar roughness scale over water
694 SUBROUTINE znot_t_v6(uref, znott)
695 use machine , only : kind_phys
696 IMPLICIT NONE
697
698
699
700 REAL(kind=kind_phys), INTENT(IN) :: uref
701 REAL(kind=kind_phys), INTENT(OUT):: znott
702 real(kind=kind_phys), parameter :: p00 = 1.100000000000000e-04,
703 & p15 = -9.144581627678278e-10, p14 = 7.020346616456421e-08,
704 & p13 = -2.155602086883837e-06, p12 = 3.333848806567684e-05,
705 & p11 = -2.628501274963990e-04, p10 = 8.634221567969181e-04,
706
707 & p25 = -8.654513012535990e-12, p24 = 1.232380050058077e-09,
708 & p23 = -6.837922749505057e-08, p22 = 1.871407733439947e-06,
709 & p21 = -2.552246987137160e-05, p20 = 1.428968311457630e-04,
710
711 & p35 = 3.207515102100162e-12, p34 = -2.945761895342535e-10,
712 & p33 = 8.788972147364181e-09, p32 = -3.814457439412957e-08,
713 & p31 = -2.448983648874671e-06, p30 = 3.436721779020359e-05,
714
715 & p45 = -3.530687797132211e-11, p44 = 3.939867958963747e-09,
716 & p43 = -1.227668406985956e-08, p42 = -1.367469811838390e-05,
717 & p41 = 5.988240863928883e-04, p40 = -7.746288511324971e-03,
718
719 & p56 = -1.187982453329086e-13, p55 = 4.801984186231693e-11,
720 & p54 = -8.049200462388188e-09, p53 = 7.169872601310186e-07,
721 & p52 = -3.581694433758150e-05, p51 = 9.503919224192534e-04,
722 & p50 = -1.036679430885215e-02,
723
724 & p60 = 4.751256171799112e-05
725
726 if (uref >= 0.0 .and. uref < 5.9 ) then
727 znott = p00
728 elseif (uref >= 5.9 .and. uref <= 15.4) then
729 znott = p10 + uref * (p11 + uref * (p12 + uref * (p13
730 & + uref * (p14 + uref * p15))))
731 elseif (uref > 15.4 .and. uref <= 21.6) then
732 znott = p20 + uref * (p21 + uref * (p22 + uref * (p23
733 & + uref * (p24 + uref * p25))))
734 elseif (uref > 21.6 .and. uref <= 42.2) then
735 znott = p30 + uref * (p31 + uref * (p32 + uref * (p33
736 & + uref * (p34 + uref * p35))))
737 elseif ( uref > 42.2 .and. uref <= 53.3) then
738 znott = p40 + uref * (p41 + uref * (p42 + uref * (p43
739 & + uref * (p44 + uref * p45))))
740 elseif ( uref > 53.3 .and. uref <= 80.0) then
741 znott = p50 + uref * (p51 + uref * (p52 + uref * (p53
742 & + uref * (p54 + uref * (p55 + uref * p56)))))
743 elseif ( uref > 80.0) then
744 znott = p60
745 else
746 print*, 'Wrong input uref value:',uref
747 endif
748
749 END SUBROUTINE znot_t_v6
750
751
756!
757!! Bin Liu, NOAA/NCEP/EMC 2018
758!
759!! uref(m/s) : wind speed at 10-m height
760!! znotm(meter): areodynamical roughness scale over water
761 SUBROUTINE znot_m_v7(uref, znotm)
762 use machine , only : kind_phys
763 IMPLICIT NONE
764
765
766
767 REAL(kind=kind_phys), INTENT(IN) :: uref
768 REAL(kind=kind_phys), INTENT(OUT):: znotm
769
770 real(kind=kind_phys), parameter :: p13 = -1.296521881682694e-02,
771 & p12 = 2.855780863283819e-01, p11 = -1.597898515251717e+00,
772 & p10 = -8.396975715683501e+00,
773
774 & p25 = 3.790846746036765e-10, p24 = 3.281964357650687e-09,
775 & p23 = 1.962282433562894e-07, p22 = -1.240239171056262e-06,
776 & p21 = 1.739759082358234e-07, p20 = 2.147264020369413e-05,
777
778 & p35 = 1.897534489606422e-07, p34 = -3.019495980684978e-05,
779 & p33 = 1.931392924987349e-03, p32 = -6.797293095862357e-02,
780 & p31 = 1.346757797103756e+00, p30 = -1.707846930193362e+01,
781
782 & p40 = 3.371427455376717e-04
783
784 if (uref >= 0.0 .and. uref <= 6.5 ) then
785 znotm = exp( p10 + uref * (p11 + uref * (p12 + uref * p13)))
786 elseif (uref > 6.5 .and. uref <= 15.7) then
787 znotm = p20 + uref * (p21 + uref * (p22 + uref * (p23
788 & + uref * (p24 + uref * p25))))
789 elseif (uref > 15.7 .and. uref <= 53.0) then
790 znotm = exp( p30 + uref * (p31 + uref * (p32 + uref * (p33
791 & + uref * (p34 + uref * p35)))))
792 elseif ( uref > 53.0) then
793 znotm = p40
794 else
795 print*, 'Wrong input uref value:',uref
796 endif
797
798 END SUBROUTINE znot_m_v7
799
809 SUBROUTINE znot_t_v7(uref, znott)
810 use machine , only : kind_phys
811 IMPLICIT NONE
812
813!
814
815 REAL(kind=kind_phys), INTENT(IN) :: uref
816 REAL(kind=kind_phys), INTENT(OUT):: znott
817
818 real(kind=kind_phys), parameter :: p00 = 1.100000000000000e-04,
819
820 & p15 = -9.193764479895316e-10, p14 = 7.052217518653943e-08,
821 & p13 = -2.163419217747114e-06, p12 = 3.342963077911962e-05,
822 & p11 = -2.633566691328004e-04, p10 = 8.644979973037803e-04,
823
824 & p25 = -9.402722450219142e-12, p24 = 1.325396583616614e-09,
825 & p23 = -7.299148051141852e-08, p22 = 1.982901461144764e-06,
826 & p21 = -2.680293455916390e-05, p20 = 1.484341646128200e-04,
827
828 & p35 = 7.921446674311864e-12, p34 = -1.019028029546602e-09,
829 & p33 = 5.251986927351103e-08, p32 = -1.337841892062716e-06,
830 & p31 = 1.659454106237737e-05, p30 = -7.558911792344770e-05,
831
832 & p45 = -2.694370426850801e-10, p44 = 5.817362913967911e-08,
833 & p43 = -5.000813324746342e-06, p42 = 2.143803523428029e-04,
834 & p41 = -4.588070983722060e-03, p40 = 3.924356617245624e-02,
835
836 & p56 = -1.663918773476178e-13, p55 = 6.724854483077447e-11,
837 & p54 = -1.127030176632823e-08, p53 = 1.003683177025925e-06,
838 & p52 = -5.012618091180904e-05, p51 = 1.329762020689302e-03,
839 & p50 = -1.450062148367566e-02, p60 = 6.840803042788488e-05
840
841 if (uref >= 0.0 .and. uref < 5.9 ) then
842 znott = p00
843 elseif (uref >= 5.9 .and. uref <= 15.4) then
844 znott = p10 + uref * (p11 + uref * (p12 + uref * (p13
845 & + uref * (p14 + uref * p15))))
846 elseif (uref > 15.4 .and. uref <= 21.6) then
847 znott = p20 + uref * (p21 + uref * (p22 + uref * (p23
848 & + uref * (p24 + uref * p25))))
849 elseif (uref > 21.6 .and. uref <= 42.6) then
850 znott = p30 + uref * (p31 + uref * (p32 + uref * (p33
851 & + uref * (p34 + uref * p35))))
852 elseif ( uref > 42.6 .and. uref <= 53.0) then
853 znott = p40 + uref * (p41 + uref * (p42 + uref * (p43
854 & + uref * (p44 + uref * p45))))
855 elseif ( uref > 53.0 .and. uref <= 80.0) then
856 znott = p50 + uref * (p51 + uref * (p52 + uref * (p53
857 & + uref * (p54 + uref * (p55 + uref * p56)))))
858 elseif ( uref > 80.0) then
859 znott = p60
860 else
861 print*, 'Wrong input uref value:',uref
862 endif
863
864 END SUBROUTINE znot_t_v7
866 end module sfc_diff
subroutine znot_m_v7(uref, znotm)
Calculate areodynamical roughness over water with input 10-m wind For low-to-moderate winds,...
Definition sfc_diff.f:762
subroutine znot_t_v7(uref, znott)
Calculate scalar roughness over water with input 10-m wind For low-to-moderate winds,...
Definition sfc_diff.f:810
subroutine, public sfc_diff_run(im, rvrdm1, eps, epsm1, grav, ps, t1, q1, z1, garea, wind, prsl1, prslki, prsik1, prslk1, sigmaf, vegtype, shdmax, ivegsrc, z0pert, ztpert, flag_iter, redrag, flag_lakefreeze, lakefrac, fice, u10m, v10m, sfc_z0_type, u1, v1, usfco, vsfco, use_oceanuv, wet, dry, icy, thsfc_loc, tskin_wat, tskin_lnd, tskin_ice, tsurf_wat, tsurf_lnd, tsurf_ice, z0rl_wat, z0rl_lnd, z0rl_ice, z0rl_wav, ustar_wat, ustar_lnd, ustar_ice, cm_wat, cm_lnd, cm_ice, ch_wat, ch_lnd, ch_ice, rb_wat, rb_lnd, rb_ice, stress_wat, stress_lnd, stress_ice, fm_wat, fm_lnd, fm_ice, fh_wat, fh_lnd, fh_ice, fm10_wat, fm10_lnd, fm10_ice, fh2_wat, fh2_lnd, fh2_ice, ztmax_wat, ztmax_lnd, ztmax_ice, zvfun, errmsg, errflg)
Definition sfc_diff.f:84
subroutine znot_m_v6(uref, znotm)
add fitted z0,zt curves for hurricane application (used in HWRF/HMON) Weiguo Wang,...
Definition sfc_diff.f:641
subroutine znot_t_v6(uref, znott)
Calculate scalar roughness over water with input 10-m wind For low-to-moderate winds,...
Definition sfc_diff.f:695
subroutine, public stability(z1, zvfun, gdx, tv1, thv1, wind, z0max, ztmax, tvs, grav, thsfc_loc, rb, fm, fh, fm10, fh2, cm, ch, stress, ustar)
Definition sfc_diff.f:469
This module contains the CCPP-compliant GFS surface layer scheme.
Definition sfc_diff.f:7