42 integer,
private :: LAND_SFC_THIN_SLAB_itr_max = 100
44 real(RP),
private :: LAND_SFC_THIN_SLAB_dTS_max = 5.0e-2_rp
45 real(RP),
private :: LAND_SFC_THIN_SLAB_res_min = 1.0e+0_rp
46 real(RP),
private :: LAND_SFC_THIN_SLAB_err_min = 1.0e-2_rp
47 real(RP),
private :: LAND_SFC_THIN_SLAB_dreslim = 1.0e+2_rp
57 character(len=*),
intent(in) :: land_type
59 namelist / param_land_sfc_thin_slab / &
60 land_sfc_thin_slab_itr_max, &
61 land_sfc_thin_slab_dts_max, &
62 land_sfc_thin_slab_res_min, &
63 land_sfc_thin_slab_err_min, &
64 land_sfc_thin_slab_dreslim
70 if(
io_l )
write(
io_fid_log,*)
'++++++ Module[THIN-SLAB] / Categ[LAND SFC] / Origin[SCALElib]' 74 read(
io_fid_conf,nml=param_land_sfc_thin_slab,iostat=ierr)
76 if(
io_l )
write(
io_fid_log,*)
'*** Not found namelist. Default used.' 77 elseif( ierr > 0 )
then 78 write(*,*)
'xxx Not appropriate names in namelist PARAM_LAND_SFC_THIN_SLAB. Check!' 133 hydrometeor_lhv => atmos_hydrometeor_lhv
135 qsat => atmos_saturation_pres2qsat_all
141 real(RP),
parameter :: dts0 = 1.0e-4_rp
143 real(RP),
parameter :: redf_min = 1.0e-2_rp
144 real(RP),
parameter :: redf_max = 1.0e+0_rp
145 real(RP),
parameter :: tfa = 0.5e+0_rp
146 real(RP),
parameter :: tfb = 1.1e+0_rp
149 real(RP),
intent(out) :: lst_t(
ia,
ja)
150 real(RP),
intent(out) :: zmflx(
ia,
ja)
151 real(RP),
intent(out) :: xmflx(
ia,
ja)
152 real(RP),
intent(out) :: ymflx(
ia,
ja)
153 real(RP),
intent(out) :: shflx(
ia,
ja)
154 real(RP),
intent(out) :: lhflx(
ia,
ja)
155 real(RP),
intent(out) :: ghflx(
ia,
ja)
156 real(RP),
intent(out) :: u10 (
ia,
ja)
157 real(RP),
intent(out) :: v10 (
ia,
ja)
158 real(RP),
intent(out) :: t2 (
ia,
ja)
159 real(RP),
intent(out) :: q2 (
ia,
ja)
161 real(RP),
intent(in) :: tmpa(
ia,
ja)
162 real(RP),
intent(in) :: prsa(
ia,
ja)
163 real(RP),
intent(in) :: wa (
ia,
ja)
164 real(RP),
intent(in) :: ua (
ia,
ja)
165 real(RP),
intent(in) :: va (
ia,
ja)
166 real(RP),
intent(in) :: rhoa(
ia,
ja)
167 real(RP),
intent(in) :: qva (
ia,
ja)
168 real(RP),
intent(in) :: z1 (
ia,
ja)
169 real(RP),
intent(in) :: pbl (
ia,
ja)
170 real(RP),
intent(in) :: prss(
ia,
ja)
171 real(RP),
intent(in) :: lwd (
ia,
ja)
172 real(RP),
intent(in) :: swd (
ia,
ja)
174 real(RP),
intent(in) :: tg (
ia,
ja)
175 real(RP),
intent(in) :: lst (
ia,
ja)
176 real(RP),
intent(in) :: qvef (
ia,
ja)
177 real(RP),
intent(in) :: alb_lw(
ia,
ja)
178 real(RP),
intent(in) :: alb_sw(
ia,
ja)
179 real(RP),
intent(in) :: dzg (
ia,
ja)
180 real(RP),
intent(in) :: tcs (
ia,
ja)
181 real(RP),
intent(in) :: z0m (
ia,
ja)
182 real(RP),
intent(in) :: z0h (
ia,
ja)
183 real(RP),
intent(in) :: z0e (
ia,
ja)
184 real(DP),
intent(in) :: dt
187 real(RP) :: lst1(
ia,
ja)
194 real(RP) :: ustar, dustar
195 real(RP) :: tstar, dtstar
196 real(RP) :: qstar, dqstar
197 real(RP) :: uabs, duabs
198 real(RP) :: qvsat, dqvsat
199 real(RP) :: qvs, dqvs
205 real(RP) :: lhv(
ia,
ja)
212 call hydrometeor_lhv( lhv(:,:), tmpa(:,:) )
228 oldres = huge(0.0_rp)
231 do n = 1, land_sfc_thin_slab_itr_max
240 qvs = ( 1.0_rp - qvef(i,j) ) * qva(i,j) + qvef(i,j) * qvsat
241 dqvs = ( 1.0_rp - qvef(i,j) ) * qva(i,j) + qvef(i,j) * dqvsat
288 res = ( 1.0_rp - alb_sw(i,j) ) * swd(i,j) &
289 + ( 1.0_rp - alb_lw(i,j) ) * ( lwd(i,j) - stb * lst1(i,j)**4 ) &
290 + cpdry * rhoa(i,j) * ustar * tstar &
291 + lhv(i,j) * rhoa(i,j) * ustar * qstar &
292 - 2.0_rp * tcs(i,j) * ( lst1(i,j) - tg(i,j) ) / dzg(i,j)
295 dres = -4.0_rp * ( 1.0_rp - alb_lw(i,j) ) * stb * lst1(i,j)**3 &
296 + cpdry * rhoa(i,j) * ( (dustar-ustar)/dts0 * tstar + ustar * (dtstar-tstar)/dts0 ) &
297 + lhv(i,j) * rhoa(i,j) * ( (dustar-ustar)/dts0 * qstar + ustar * (dqstar-qstar)/dts0 ) &
298 - 2.0_rp * tcs(i,j) / dzg(i,j)
301 if( abs( res ) < land_sfc_thin_slab_res_min .or. &
302 abs( res/dres ) < land_sfc_thin_slab_err_min )
then 307 if( abs(dres) * land_sfc_thin_slab_dreslim < abs(res) )
then 312 if( dres < 0.0_rp )
then 313 if( abs(res) > abs(oldres) )
then 314 redf = max( tfa*abs(redf), redf_min )
316 redf = min( tfb*abs(redf), redf_max )
323 lst1(i,j) = lst1(i,j) - redf * res / dres
331 lst1(i,j) = min( max( lst1(i,j), &
332 lst(i,j) - land_sfc_thin_slab_dts_max *
real( dt, kind=RP ) ), &
333 lst (i,j) + land_sfc_thin_slab_dts_max *
real( dt, kind=RP ) )
335 if( n > land_sfc_thin_slab_itr_max )
then 337 if(
io_l )
write(
io_fid_log,
'(A)' )
'Warning: land surface tempearture was not converged.' 340 if(
io_l )
write(
io_fid_log,
'(A,I32)' )
'DEBUG --- number of i [no unit] :', i
341 if(
io_l )
write(
io_fid_log,
'(A,I32)' )
'DEBUG --- number of j [no unit] :', j
343 if(
io_l )
write(
io_fid_log,
'(A,I32)' )
'DEBUG --- loop number [no unit] :', n
344 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- Residual [J/m2/s] :', res
345 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- delta Residual [J/m2/s] :', dres
347 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- temperature [K] :', tmpa(i,j)
348 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- pressure [Pa] :', prsa(i,j)
349 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- velocity w [m/s] :', wa(i,j)
350 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- velocity u [m/s] :', ua(i,j)
351 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- velocity v [m/s] :', va(i,j)
352 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- density [kg/m3] :', rhoa(i,j)
353 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- water vapor mass ratio [kg/kg] :', qva(i,j)
354 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- cell center height [m] :', z1(i,j)
355 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- atmospheric mixing layer height [m] :', pbl(i,j)
356 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- pressure at the surface [Pa] :', prss(i,j)
357 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- downward long-wave radiation [J/m2/s] :', lwd(i,j)
358 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- downward short-wave radiation [J/m2/s] :', swd(i,j)
360 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- soil temperature [K] :', tg(i,j)
361 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- land surface temperature [K] :', lst(i,j)
362 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- efficiency of evaporation [1] :', qvef(i,j)
363 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- surface albedo for LW [1] :', alb_lw(i,j)
364 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- surface albedo for SW [1] :', alb_sw(i,j)
365 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- soil depth [m] :', dzg(i,j)
366 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- thermal conductivity for soil [J/m/K/s] :', tcs(i,j)
367 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- roughness length for momemtum [m] :', z0m(i,j)
368 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- roughness length for heat [m] :', z0h(i,j)
369 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- roughness length for vapor [m] :', z0e(i,j)
371 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- latent heat [J/kg] :', lhv(i,j)
372 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- friction velocity [m] :', ustar
373 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- friction potential temperature [K] :', tstar
374 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- friction water vapor mass ratio [kg/kg] :', qstar
375 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- d(friction velocity) [m] :', dustar
376 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- d(friction potential temperature) [K] :', dtstar
377 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- d(friction water vapor mass ratio) [kg/kg] :', dqstar
378 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- modified absolute velocity [m/s] :', uabs
379 if(
io_l )
write(
io_fid_log,
'(A,F32.16)')
'DEBUG --- next land surface temperature [K] :', lst1(i,j)
382 if ( .NOT. ( res > -1.0_rp .OR. res < 1.0_rp ) )
then 383 write(*,*)
'xxx NaN is detected for land surface temperature.' 385 write(*,*)
'DEBUG --- PRC_myrank :',
prc_myrank 386 write(*,*)
'DEBUG --- number of i :', i
387 write(*,*)
'DEBUG --- number of j :', j
389 write(*,*)
'DEBUG --- Residual [J/m2/s] :', res
390 write(*,*)
'DEBUG --- delta Residual [J/m2/s] :', dres
392 write(*,*)
'DEBUG --- temperature [K] :', tmpa(i,j)
393 write(*,*)
'DEBUG --- pressure [Pa] :', prsa(i,j)
394 write(*,*)
'DEBUG --- velocity w [m/s] :', wa(i,j)
395 write(*,*)
'DEBUG --- velocity u [m/s] :', ua(i,j)
396 write(*,*)
'DEBUG --- velocity v [m/s] :', va(i,j)
397 write(*,*)
'DEBUG --- density [kg/m3] :', rhoa(i,j)
398 write(*,*)
'DEBUG --- water vapor mass ratio [kg/kg] :', qva(i,j)
399 write(*,*)
'DEBUG --- cell center height [m] :', z1(i,j)
400 write(*,*)
'DEBUG --- atmospheric mixing layer height [m] :', pbl(i,j)
401 write(*,*)
'DEBUG --- pressure at the surface [Pa] :', prss(i,j)
402 write(*,*)
'DEBUG --- downward long-wave radiation [J/m2/s] :', lwd(i,j)
403 write(*,*)
'DEBUG --- downward short-wave radiation [J/m2/s] :', swd(i,j)
405 write(*,*)
'DEBUG --- soil temperature [K] :', tg(i,j)
406 write(*,*)
'DEBUG --- land surface temperature [K] :', lst(i,j)
407 write(*,*)
'DEBUG --- efficiency of evaporation [1] :', qvef(i,j)
408 write(*,*)
'DEBUG --- surface albedo for LW [1] :', alb_lw(i,j)
409 write(*,*)
'DEBUG --- surface albedo for SW [1] :', alb_sw(i,j)
410 write(*,*)
'DEBUG --- soil depth [m] :', dzg(i,j)
411 write(*,*)
'DEBUG --- thermal conductivity for soil [J/m/K/s] :', tcs(i,j)
412 write(*,*)
'DEBUG --- roughness length for momemtum [m] :', z0m(i,j)
413 write(*,*)
'DEBUG --- roughness length for heat [m] :', z0h(i,j)
414 write(*,*)
'DEBUG --- roughness length for vapor [m] :', z0e(i,j)
416 write(*,*)
'DEBUG --- latent heat [J/kg] :', lhv(i,j)
417 write(*,*)
'DEBUG --- friction velocity [m] :', ustar
418 write(*,*)
'DEBUG --- friction potential temperature [K] :', tstar
419 write(*,*)
'DEBUG --- friction water vapor mass ratio [kg/kg] :', qstar
420 write(*,*)
'DEBUG --- d(friction velocity) [m] :', dustar
421 write(*,*)
'DEBUG --- d(friction potential temperature) [K] :', dtstar
422 write(*,*)
'DEBUG --- d(friction water vapor mass ratio) [kg/kg] :', dqstar
423 write(*,*)
'DEBUG --- modified absolute velocity [m/s] :', uabs
424 write(*,*)
'DEBUG --- next land surface temperature [K] :', lst1(i,j)
434 lst_t(i,j) = ( lst1(i,j) - lst(i,j) ) / dt
449 qvs = ( 1.0_rp - qvef(i,j) ) * qva(i,j) + qvef(i,j) * qvsat
473 zmflx(i,j) = -rhoa(i,j) * ustar**2 / uabs * wa(i,j)
474 xmflx(i,j) = -rhoa(i,j) * ustar**2 / uabs * ua(i,j)
475 ymflx(i,j) = -rhoa(i,j) * ustar**2 / uabs * va(i,j)
477 shflx(i,j) = - rhoa(i,j) * ustar * tstar &
478 * cpdry * ( prss(i,j) / pre00 )**( rdry/cpdry )
479 lhflx(i,j) = - rhoa(i,j) * ustar * qstar * lhv(i,j)
481 ghflx(i,j) = -2.0_rp * tcs(i,j) * ( lst1(i,j) - tg(i,j) ) / dzg(i,j)
484 res = ( 1.0_rp - alb_sw(i,j) ) * swd(i,j) &
485 + ( 1.0_rp - alb_lw(i,j) ) * ( lwd(i,j) - stb * lst1(i,j)**4 ) &
486 - shflx(i,j) - lhflx(i,j) + ghflx(i,j)
489 ghflx(i,j) = ghflx(i,j) - res
498 u10(i,j) = ua(i,j) * log( 10.0_rp / z0m(i,j) ) / log( z1(i,j) / z0m(i,j) )
499 v10(i,j) = va(i,j) * log( 10.0_rp / z0m(i,j) ) / log( z1(i,j) / z0m(i,j) )
500 t2(i,j) = lst1(i,j) + ( tmpa(i,j) - lst1(i,j) ) * ( log( 2.0_rp / z0m(i,j) ) * log( 2.0_rp / z0h(i,j) ) ) &
501 / ( log( z1(i,j) / z0m(i,j) ) * log( z1(i,j) / z0h(i,j) ) )
502 q2(i,j) = qvs + ( qva(i,j) - qvs ) * ( log( 2.0_rp / z0m(i,j) ) * log( 2.0_rp / z0e(i,j) ) ) &
503 / ( log( z1(i,j) / z0m(i,j) ) * log( z1(i,j) / z0e(i,j) ) )
integer, public is
start point of inner domain: x, local
module LAND / Surface fluxes with thin-slab land model
integer, public je
end point of inner domain: y, local
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
module ATMOSPHERE / Saturation adjustment
subroutine, public prc_mpistop
Abort MPI.
real(rp), parameter, public const_stb
Stefan-Boltzman constant [W/m2/K4].
logical, public io_l
output log or not? (this process)
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
logical, public io_nml
output log or not? (for namelist, this process)
integer, public ia
of whole cells: x, local, with HALO
real(rp), public const_pre00
pressure reference [Pa]
procedure(bc), pointer, public bulkflux
integer, public js
start point of inner domain: y, local
subroutine, public land_sfc_thin_slab_setup(LAND_TYPE)
Setup.
subroutine, public land_sfc_thin_slab(LST_t, ZMFLX, XMFLX, YMFLX, SHFLX, LHFLX, GHFLX, U10, V10, T2, Q2, TMPA, PRSA, WA, UA, VA, RHOA, QVA, Z1, PBL, PRSS, LWD, SWD, TG, LST, QVEF, ALB_LW, ALB_SW, DZG, TCS, Z0M, Z0H, Z0E, dt)
integer, public prc_myrank
process num in local communicator
integer, public ie
end point of inner domain: x, local
integer, public io_fid_conf
Config file ID.
real(rp), dimension(:,:), allocatable, public landuse_fact_land
land factor
integer, public io_fid_log
Log file ID.
integer, public io_fid_nml
Log file ID (only for output namelist)
integer, public ja
of whole cells: y, local, with HALO