SCALE-RM
scale_land_sfc_thin-slab.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_stdio
19  !-----------------------------------------------------------------------------
20  implicit none
21  private
22  !-----------------------------------------------------------------------------
23  !
24  !++ Public procedure
25  !
26  public :: land_sfc_thin_slab_setup
27  public :: land_sfc_thin_slab
28 
29  !-----------------------------------------------------------------------------
30  !
31  !++ Public parameters & variables
32  !
33  !-----------------------------------------------------------------------------
34  !
35  !++ Private procedure
36  !
37  !-----------------------------------------------------------------------------
38  !
39  !++ Private parameters & variables
40  !
41  !-----------------------------------------------------------------------------
42  integer, private :: LAND_SFC_THIN_SLAB_itr_max = 100 ! maximum iteration number
43 
44  real(RP), private :: LAND_SFC_THIN_SLAB_dTS_max = 5.0e-2_rp ! maximum delta surface temperature [K/s]
45  real(RP), private :: LAND_SFC_THIN_SLAB_res_min = 1.0e+0_rp ! minimum value of residual
46  real(RP), private :: LAND_SFC_THIN_SLAB_err_min = 1.0e-2_rp ! minimum value of error
47  real(RP), private :: LAND_SFC_THIN_SLAB_dreslim = 1.0e+2_rp ! limiter of d(residual)
48 
49 contains
50  !-----------------------------------------------------------------------------
52  subroutine land_sfc_thin_slab_setup( LAND_TYPE )
53  use scale_process, only: &
55  implicit none
56 
57  character(len=*), intent(in) :: land_type
58 
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
65 
66  integer :: ierr
67  !---------------------------------------------------------------------------
68 
69  if( io_l ) write(io_fid_log,*)
70  if( io_l ) write(io_fid_log,*) '++++++ Module[THIN-SLAB] / Categ[LAND SFC] / Origin[SCALElib]'
71 
72  !--- read namelist
73  rewind(io_fid_conf)
74  read(io_fid_conf,nml=param_land_sfc_thin_slab,iostat=ierr)
75  if( ierr < 0 ) then !--- missing
76  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
77  elseif( ierr > 0 ) then !--- fatal error
78  write(*,*) 'xxx Not appropriate names in namelist PARAM_LAND_SFC_THIN_SLAB. Check!'
79  call prc_mpistop
80  endif
81  if( io_nml ) write(io_fid_nml,nml=param_land_sfc_thin_slab)
82 
83  return
84  end subroutine land_sfc_thin_slab_setup
85 
86  !-----------------------------------------------------------------------------
87  subroutine land_sfc_thin_slab( &
88  LST_t, &
89  ZMFLX, &
90  XMFLX, &
91  YMFLX, &
92  SHFLX, &
93  LHFLX, &
94  GHFLX, &
95  U10, &
96  V10, &
97  T2, &
98  Q2, &
99  TMPA, &
100  PRSA, &
101  WA, &
102  UA, &
103  VA, &
104  RHOA, &
105  QVA, &
106  Z1, &
107  PBL, &
108  PRSS, &
109  LWD, &
110  SWD, &
111  TG, &
112  LST, &
113  QVEF, &
114  ALB_LW, &
115  ALB_SW, &
116  DZG, &
117  TCS, &
118  Z0M, &
119  Z0H, &
120  Z0E, &
121  dt )
122  use scale_process, only: &
123  prc_myrank, &
125  use scale_const, only: &
126  pre00 => const_pre00, &
127  rdry => const_rdry, &
128  cpdry => const_cpdry, &
129  stb => const_stb
130  use scale_landuse, only: &
132  use scale_atmos_hydrometeor, only: &
133  hydrometeor_lhv => atmos_hydrometeor_lhv
134  use scale_atmos_saturation, only: &
135  qsat => atmos_saturation_pres2qsat_all
136  use scale_bulkflux, only: &
137  bulkflux
138  implicit none
139 
140  ! parameters
141  real(RP), parameter :: dts0 = 1.0e-4_rp ! delta surface temp.
142 
143  real(RP), parameter :: redf_min = 1.0e-2_rp ! minimum reduced factor
144  real(RP), parameter :: redf_max = 1.0e+0_rp ! maximum reduced factor
145  real(RP), parameter :: tfa = 0.5e+0_rp ! factor a in Tomita (2009)
146  real(RP), parameter :: tfb = 1.1e+0_rp ! factor b in Tomita (2009)
147 
148  ! arguments
149  real(RP), intent(out) :: lst_t(ia,ja) ! tendency of LST
150  real(RP), intent(out) :: zmflx(ia,ja) ! z-momentum flux at the surface [kg/m/s2]
151  real(RP), intent(out) :: xmflx(ia,ja) ! x-momentum flux at the surface [kg/m/s2]
152  real(RP), intent(out) :: ymflx(ia,ja) ! y-momentum flux at the surface [kg/m/s2]
153  real(RP), intent(out) :: shflx(ia,ja) ! sensible heat flux at the surface [J/m2/s]
154  real(RP), intent(out) :: lhflx(ia,ja) ! latent heat flux at the surface [J/m2/s]
155  real(RP), intent(out) :: ghflx(ia,ja) ! ground heat flux at the surface [J/m2/s]
156  real(RP), intent(out) :: u10 (ia,ja) ! velocity u at 10m [m/s]
157  real(RP), intent(out) :: v10 (ia,ja) ! velocity v at 10m [m/s]
158  real(RP), intent(out) :: t2 (ia,ja) ! temperature at 2m [K]
159  real(RP), intent(out) :: q2 (ia,ja) ! water vapor at 2m [kg/kg]
160 
161  real(RP), intent(in) :: tmpa(ia,ja) ! temperature at the lowest atmospheric layer [K]
162  real(RP), intent(in) :: prsa(ia,ja) ! pressure at the lowest atmospheric layer [Pa]
163  real(RP), intent(in) :: wa (ia,ja) ! velocity w at the lowest atmospheric layer [m/s]
164  real(RP), intent(in) :: ua (ia,ja) ! velocity u at the lowest atmospheric layer [m/s]
165  real(RP), intent(in) :: va (ia,ja) ! velocity v at the lowest atmospheric layer [m/s]
166  real(RP), intent(in) :: rhoa(ia,ja) ! density at the lowest atmospheric layer [kg/m3]
167  real(RP), intent(in) :: qva (ia,ja) ! ratio of water vapor mass to total mass at the lowest atmospheric layer [kg/kg]
168  real(RP), intent(in) :: z1 (ia,ja) ! cell center height at the lowest atmospheric layer [m]
169  real(RP), intent(in) :: pbl (ia,ja) ! the top of atmospheric mixing layer [m]
170  real(RP), intent(in) :: prss(ia,ja) ! pressure at the surface [Pa]
171  real(RP), intent(in) :: lwd (ia,ja) ! downward long-wave radiation flux at the surface [J/m2/s]
172  real(RP), intent(in) :: swd (ia,ja) ! downward short-wave radiation flux at the surface [J/m2/s]
173 
174  real(RP), intent(in) :: tg (ia,ja) ! soil temperature [K]
175  real(RP), intent(in) :: lst (ia,ja) ! land surface temperature [K]
176  real(RP), intent(in) :: qvef (ia,ja) ! efficiency of evaporation (0-1)
177  real(RP), intent(in) :: alb_lw(ia,ja) ! surface albedo for LW (0-1)
178  real(RP), intent(in) :: alb_sw(ia,ja) ! surface albedo for SW (0-1)
179  real(RP), intent(in) :: dzg (ia,ja) ! soil depth [m]
180  real(RP), intent(in) :: tcs (ia,ja) ! thermal conductivity for soil [J/m/K/s]
181  real(RP), intent(in) :: z0m (ia,ja) ! roughness length for momemtum [m]
182  real(RP), intent(in) :: z0h (ia,ja) ! roughness length for heat [m]
183  real(RP), intent(in) :: z0e (ia,ja) ! roughness length for vapor [m]
184  real(DP), intent(in) :: dt ! delta time
185 
186  ! works
187  real(RP) :: lst1(ia,ja)
188 
189  real(RP) :: res ! residual
190  real(RP) :: dres ! d(residual)/dLST
191  real(RP) :: oldres ! residual in previous step
192  real(RP) :: redf ! reduced factor
193 
194  real(RP) :: ustar, dustar ! friction velocity [m]
195  real(RP) :: tstar, dtstar ! friction potential temperature [K]
196  real(RP) :: qstar, dqstar ! friction water vapor mass ratio [kg/kg]
197  real(RP) :: uabs, duabs ! modified absolute velocity [m/s]
198  real(RP) :: qvsat, dqvsat ! saturation water vapor mixing ratio at surface [kg/kg]
199  real(RP) :: qvs, dqvs ! water vapor mixing ratio at surface [kg/kg]
200 
201  real(RP) :: fracu10 ! calculation parameter for U10 [-]
202  real(RP) :: fract2 ! calculation parameter for T2 [-]
203  real(RP) :: fracq2 ! calculation parameter for Q2 [-]
204 
205  real(RP) :: lhv(ia,ja) ! latent heat of vaporization [J/kg]
206 
207  integer :: i, j, n
208  !---------------------------------------------------------------------------
209 
210  if( io_l ) write(io_fid_log,*) '*** Land surface step: Thin-Slab'
211 
212  call hydrometeor_lhv( lhv(:,:), tmpa(:,:) )
213 
214  ! copy land surfce temperature for iteration
215  do j = js, je
216  do i = is, ie
217  lst1(i,j) = lst(i,j)
218  end do
219  end do
220 
221  ! update surface temperature
222  do j = js, je
223  do i = is, ie
224 
225  if( landuse_fact_land(i,j) > 0.0_rp ) then
226 
227  redf = 1.0_rp
228  oldres = huge(0.0_rp)
229 
230  ! modified Newton-Raphson method (Tomita 2009)
231  do n = 1, land_sfc_thin_slab_itr_max
232 
233  call qsat( qvsat, & ! [OUT]
234  lst1(i,j), & ! [IN]
235  prss(i,j) ) ! [IN]
236  call qsat( dqvsat, & ! [OUT]
237  lst1(i,j)+dts0, & ! [IN]
238  prss(i,j) ) ! [IN]
239 
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
242 
243  call bulkflux( &
244  ustar, & ! [OUT]
245  tstar, & ! [OUT]
246  qstar, & ! [OUT]
247  uabs, & ! [OUT]
248  fracu10, & ! [OUT] ! not used
249  fract2, & ! [OUT] ! not used
250  fracq2, & ! [OUT] ! not used
251  tmpa(i,j), & ! [IN]
252  lst1(i,j), & ! [IN]
253  prsa(i,j), & ! [IN]
254  prss(i,j), & ! [IN]
255  qva(i,j), & ! [IN]
256  qvs, & ! [IN]
257  ua(i,j), & ! [IN]
258  va(i,j), & ! [IN]
259  z1(i,j), & ! [IN]
260  pbl(i,j), & ! [IN]
261  z0m(i,j), & ! [IN]
262  z0h(i,j), & ! [IN]
263  z0e(i,j) ) ! [IN]
264 
265  call bulkflux( &
266  dustar, & ! [OUT]
267  dtstar, & ! [OUT]
268  dqstar, & ! [OUT]
269  duabs, & ! [OUT]
270  fracu10, & ! [OUT] ! not used
271  fract2, & ! [OUT] ! not used
272  fracq2, & ! [OUT] ! not used
273  tmpa(i,j), & ! [IN]
274  lst1(i,j)+dts0, & ! [IN]
275  prsa(i,j), & ! [IN]
276  prss(i,j), & ! [IN]
277  qva(i,j), & ! [IN]
278  dqvs, & ! [IN]
279  ua(i,j), & ! [IN]
280  va(i,j), & ! [IN]
281  z1(i,j), & ! [IN]
282  pbl(i,j), & ! [IN]
283  z0m(i,j), & ! [IN]
284  z0h(i,j), & ! [IN]
285  z0e(i,j) ) ! [IN]
286 
287  ! calculation for residual
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)
293 
294  ! calculation for d(residual)/dLST
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)
299 
300  ! convergence test with residual and error levels
301  if( abs( res ) < land_sfc_thin_slab_res_min .or. &
302  abs( res/dres ) < land_sfc_thin_slab_err_min ) then
303  exit
304  end if
305 
306  ! stop iteration to prevent numerical error
307  if( abs(dres) * land_sfc_thin_slab_dreslim < abs(res) ) then
308  exit
309  end if
310 
311  ! calculate reduced factor
312  if( dres < 0.0_rp ) then
313  if( abs(res) > abs(oldres) ) then
314  redf = max( tfa*abs(redf), redf_min )
315  else
316  redf = min( tfb*abs(redf), redf_max )
317  end if
318  else
319  redf = -1.0_rp
320  end if
321 
322  ! estimate next surface temperature
323  lst1(i,j) = lst1(i,j) - redf * res / dres
324 
325  ! save residual in this step
326  oldres = res
327 
328  end do
329 
330  ! update land surface temperature with limitation
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 ) )
334 
335  if( n > land_sfc_thin_slab_itr_max ) then
336  ! land surface temperature was not converged
337  if( io_l ) write(io_fid_log,'(A)' ) 'Warning: land surface tempearture was not converged.'
338  if( io_l ) write(io_fid_log,'(A)' ) ''
339  if( io_l ) write(io_fid_log,'(A,I32)' ) 'DEBUG --- PRC_myrank [no unit] :', prc_myrank
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
342  if( io_l ) write(io_fid_log,'(A)' ) ''
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
346  if( io_l ) write(io_fid_log,'(A,F32.16)') ''
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)
359  if( io_l ) write(io_fid_log,'(A)' ) ''
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)
370  if( io_l ) write(io_fid_log,'(A)' ) ''
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)
380 
381  ! check NaN
382  if ( .NOT. ( res > -1.0_rp .OR. res < 1.0_rp ) ) then ! must be NaN
383  write(*,*) 'xxx NaN is detected for land surface temperature.'
384  write(*,*) ''
385  write(*,*) 'DEBUG --- PRC_myrank :', prc_myrank
386  write(*,*) 'DEBUG --- number of i :', i
387  write(*,*) 'DEBUG --- number of j :', j
388  write(*,*) ''
389  write(*,*) 'DEBUG --- Residual [J/m2/s] :', res
390  write(*,*) 'DEBUG --- delta Residual [J/m2/s] :', dres
391  write(*,*) ''
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)
404  write(*,*) ''
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)
415  write(*,*) ''
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)
425 
426  call prc_mpistop
427  endif
428 
429  end if
430 
431  end if
432 
433  ! calculate tendency
434  lst_t(i,j) = ( lst1(i,j) - lst(i,j) ) / dt
435 
436  end do
437  end do
438 
439  ! calculate surface flux
440  do j = js, je
441  do i = is, ie
442 
443  if( landuse_fact_land(i,j) > 0.0_rp ) then
444 
445  call qsat( qvsat, & ! [OUT]
446  lst1(i,j), & ! [IN]
447  prss(i,j) ) ! [IN]
448 
449  qvs = ( 1.0_rp - qvef(i,j) ) * qva(i,j) + qvef(i,j) * qvsat
450 
451  call bulkflux( &
452  ustar, & ! [OUT]
453  tstar, & ! [OUT]
454  qstar, & ! [OUT]
455  uabs, & ! [OUT]
456  fracu10, & ! [OUT]
457  fract2, & ! [OUT]
458  fracq2, & ! [OUT]
459  tmpa(i,j), & ! [IN]
460  lst1(i,j), & ! [IN]
461  prsa(i,j), & ! [IN]
462  prss(i,j), & ! [IN]
463  qva(i,j), & ! [IN]
464  qvs, & ! [IN]
465  ua(i,j), & ! [IN]
466  va(i,j), & ! [IN]
467  z1(i,j), & ! [IN]
468  pbl(i,j), & ! [IN]
469  z0m(i,j), & ! [IN]
470  z0h(i,j), & ! [IN]
471  z0e(i,j) ) ! [IN]
472 
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)
476 
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)
480 
481  ghflx(i,j) = -2.0_rp * tcs(i,j) * ( lst1(i,j) - tg(i,j) ) / dzg(i,j)
482 
483  ! calculation for residual
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)
487 
488  ! put residual in ground heat flux
489  ghflx(i,j) = ghflx(i,j) - res
490 
491  ! diagnostic variables considering unstable/stable state
492  !U10(i,j) = FracU10 * UA(i,j)
493  !V10(i,j) = FracU10 * VA(i,j)
494  !T2 (i,j) = ( 1.0_RP - FracT2 ) * LST1(i,j) + FracT2 * TMPA(i,j)
495  !Q2 (i,j) = ( 1.0_RP - FracQ2 ) * QVS + FracQ2 * QVA (i,j)
496 
497  ! diagnostic variables for neutral state
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) ) )
504 
505  else
506 
507  ! not calculate surface flux
508  zmflx(i,j) = 0.0_rp
509  xmflx(i,j) = 0.0_rp
510  ymflx(i,j) = 0.0_rp
511  shflx(i,j) = 0.0_rp
512  lhflx(i,j) = 0.0_rp
513  ghflx(i,j) = 0.0_rp
514  u10(i,j) = 0.0_rp
515  v10(i,j) = 0.0_rp
516  t2(i,j) = 0.0_rp
517  q2(i,j) = 0.0_rp
518 
519  end if
520 
521  end do
522  end do
523 
524  return
525  end subroutine land_sfc_thin_slab
526 
527 end module scale_land_sfc_thin_slab
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]
Definition: scale_const.F90:58
module ATMOSPHERE / Saturation adjustment
subroutine, public prc_mpistop
Abort MPI.
real(rp), parameter, public const_stb
Stefan-Boltzman constant [W/m2/K4].
Definition: scale_const.F90:51
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:61
module STDIO
Definition: scale_stdio.F90:12
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
module grid index
logical, public io_nml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:62
integer, public ia
of whole cells: x, local, with HALO
module LANDUSE
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:90
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)
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
integer, public prc_myrank
process num in local communicator
integer, public ie
end point of inner domain: x, local
module Surface bulk flux
module PRECISION
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
real(rp), dimension(:,:), allocatable, public landuse_fact_land
land factor
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, public io_fid_nml
Log file ID (only for output namelist)
Definition: scale_stdio.F90:57
integer, public ja
of whole cells: y, local, with HALO