SCALE-RM
Functions/Subroutines
scale_land_sfc_slab Module Reference

module LAND / Surface fluxes with slab land model More...

Functions/Subroutines

subroutine, public land_sfc_slab_setup (LAND_TYPE)
 Setup. More...
 
subroutine, public land_sfc_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)
 

Detailed Description

module LAND / Surface fluxes with slab land model

Description
Surface flux with slab land model
Author
Team SCALE
NAMELIST
  • PARAM_LAND_SFC_SLAB
    nametypedefault valuecomment
    LAND_SFC_SLAB_ITR_MAX integer 100 maximum iteration number
    LAND_SFC_SLAB_DTS_MAX real(RP) 5.0E-2_RP maximum delta surface temperature [K/s]
    LAND_SFC_SLAB_RES_MIN real(RP) 1.0E+0_RP minimum value of residual
    LAND_SFC_SLAB_ERR_MIN real(RP) 1.0E-2_RP minimum value of error
    LAND_SFC_SLAB_DRESLIM real(RP) 1.0E+2_RP limiter of d(residual)

History Output
No history output

Function/Subroutine Documentation

◆ land_sfc_slab_setup()

subroutine, public scale_land_sfc_slab::land_sfc_slab_setup ( character(len=*), intent(in)  LAND_TYPE)

Setup.

Definition at line 57 of file scale_land_sfc_slab.F90.

References scale_grid_index::ia, scale_grid_index::ie, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, scale_grid_index::is, scale_grid_index::ja, scale_grid_index::je, scale_grid_index::js, scale_landuse::landuse_fact_land, and scale_process::prc_mpistop().

Referenced by scale_land_sfc::land_sfc_setup().

57  use scale_process, only: &
59  use scale_landuse, only: &
61  implicit none
62 
63  character(len=*), intent(in) :: land_type
64 
65  namelist / param_land_sfc_slab / &
66  land_sfc_slab_itr_max, &
67  land_sfc_slab_dts_max, &
68  land_sfc_slab_res_min, &
69  land_sfc_slab_err_min, &
70  land_sfc_slab_dreslim
71 
72  integer :: i, j
73  integer :: ierr
74  !---------------------------------------------------------------------------
75 
76  if( io_l ) write(io_fid_log,*)
77  if( io_l ) write(io_fid_log,*) '++++++ Module[SLAB] / Categ[LAND SFC] / Origin[SCALElib]'
78 
79  !--- read namelist
80  rewind(io_fid_conf)
81  read(io_fid_conf,nml=param_land_sfc_slab,iostat=ierr)
82  if( ierr < 0 ) then !--- missing
83  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
84  elseif( ierr > 0 ) then !--- fatal error
85  write(*,*) 'xxx Not appropriate names in namelist PARAM_LAND_SFC_SLAB. Check!'
86  call prc_mpistop
87  endif
88  if( io_lnml ) write(io_fid_log,nml=param_land_sfc_slab)
89 
90  if( land_type == 'CONST' ) then
91  lst_update = .false.
92  else if( land_type == 'SLAB' ) then
93  lst_update = .true.
94  else
95  write(*,*) 'xxx wrong LAND_TYPE. Check!'
96  call prc_mpistop
97  end if
98 
99  ! judge to run slab land model
100  allocate( is_lnd(ia,ja) )
101 
102  do j = js, je
103  do i = is, ie
104  if( landuse_fact_land(i,j) > 0.0_rp ) then
105  is_lnd(i,j) = .true.
106  else
107  is_lnd(i,j) = .false.
108  end if
109  end do
110  end do
111 
112  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
subroutine, public prc_mpistop
Abort MPI.
integer, public ia
of x whole cells (local, with HALO)
module LANDUSE
integer, public js
start point of inner domain: y, local
module PROCESS
integer, public ie
end point of inner domain: x, local
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 ja
of y whole cells (local, with HALO)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ land_sfc_slab()

subroutine, public scale_land_sfc_slab::land_sfc_slab ( real(rp), dimension(ia,ja), intent(out)  LST_t,
real(rp), dimension(ia,ja), intent(out)  ZMFLX,
real(rp), dimension(ia,ja), intent(out)  XMFLX,
real(rp), dimension(ia,ja), intent(out)  YMFLX,
real(rp), dimension(ia,ja), intent(out)  SHFLX,
real(rp), dimension(ia,ja), intent(out)  LHFLX,
real(rp), dimension(ia,ja), intent(out)  GHFLX,
real(rp), dimension (ia,ja), intent(out)  U10,
real(rp), dimension (ia,ja), intent(out)  V10,
real(rp), dimension (ia,ja), intent(out)  T2,
real(rp), dimension (ia,ja), intent(out)  Q2,
real(rp), dimension(ia,ja), intent(in)  TMPA,
real(rp), dimension(ia,ja), intent(in)  PRSA,
real(rp), dimension (ia,ja), intent(in)  WA,
real(rp), dimension (ia,ja), intent(in)  UA,
real(rp), dimension (ia,ja), intent(in)  VA,
real(rp), dimension(ia,ja), intent(in)  RHOA,
real(rp), dimension (ia,ja), intent(in)  QVA,
real(rp), dimension (ia,ja), intent(in)  Z1,
real(rp), dimension (ia,ja), intent(in)  PBL,
real(rp), dimension(ia,ja), intent(in)  PRSS,
real(rp), dimension (ia,ja), intent(in)  LWD,
real(rp), dimension (ia,ja), intent(in)  SWD,
real(rp), dimension (ia,ja), intent(in)  TG,
real(rp), dimension (ia,ja), intent(in)  LST,
real(rp), dimension (ia,ja), intent(in)  QVEF,
real(rp), dimension(ia,ja), intent(in)  ALB_LW,
real(rp), dimension(ia,ja), intent(in)  ALB_SW,
real(rp), dimension (ia,ja), intent(in)  DZG,
real(rp), dimension (ia,ja), intent(in)  TCS,
real(rp), dimension (ia,ja), intent(in)  Z0M,
real(rp), dimension (ia,ja), intent(in)  Z0H,
real(rp), dimension (ia,ja), intent(in)  Z0E,
real(dp), intent(in)  dt 
)

Definition at line 151 of file scale_land_sfc_slab.F90.

References scale_bulkflux::bulkflux, scale_const::const_cpdry, scale_const::const_eps, scale_const::const_rdry, scale_const::const_stb, scale_grid_index::ie, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, dc_log::log(), scale_process::prc_mpistop(), and scale_process::prc_myrank.

Referenced by scale_land_sfc::land_sfc_setup().

151  use scale_process, only: &
152  prc_myrank, &
154  use scale_const, only: &
155  eps => const_eps, &
156  rdry => const_rdry, &
157  cpdry => const_cpdry, &
158  stb => const_stb
159  use scale_grid_index
160  use scale_atmos_saturation, only: &
161  qsat => atmos_saturation_pres2qsat_all
162  use scale_atmos_thermodyn, only: &
163  atmos_thermodyn_templhv
164  use scale_bulkflux, only: &
165  bulkflux
166  implicit none
167 
168  ! parameters
169  real(RP), parameter :: dts0 = 1.0e-4_rp ! delta surface temp.
170 
171  real(RP), parameter :: redf_min = 1.0e-2_rp ! minimum reduced factor
172  real(RP), parameter :: redf_max = 1.0e+0_rp ! maximum reduced factor
173  real(RP), parameter :: tfa = 0.5e+0_rp ! factor a in Tomita (2009)
174  real(RP), parameter :: tfb = 1.1e+0_rp ! factor b in Tomita (2009)
175 
176  ! arguments
177  real(RP), intent(out) :: lst_t(ia,ja) ! tendency of LST
178  real(RP), intent(out) :: zmflx(ia,ja) ! z-momentum flux at the surface [kg/m2/s]
179  real(RP), intent(out) :: xmflx(ia,ja) ! x-momentum flux at the surface [kg/m2/s]
180  real(RP), intent(out) :: ymflx(ia,ja) ! y-momentum flux at the surface [kg/m2/s]
181  real(RP), intent(out) :: shflx(ia,ja) ! sensible heat flux at the surface [J/m2/s]
182  real(RP), intent(out) :: lhflx(ia,ja) ! latent heat flux at the surface [J/m2/s]
183  real(RP), intent(out) :: ghflx(ia,ja) ! ground heat flux at the surface [J/m2/s]
184  real(RP), intent(out) :: u10 (ia,ja) ! velocity u at 10m [m/s]
185  real(RP), intent(out) :: v10 (ia,ja) ! velocity v at 10m [m/s]
186  real(RP), intent(out) :: t2 (ia,ja) ! temperature at 2m [K]
187  real(RP), intent(out) :: q2 (ia,ja) ! water vapor at 2m [kg/kg]
188 
189  real(RP), intent(in) :: tmpa(ia,ja) ! temperature at the lowest atmospheric layer [K]
190  real(RP), intent(in) :: prsa(ia,ja) ! pressure at the lowest atmospheric layer [Pa]
191  real(RP), intent(in) :: wa (ia,ja) ! velocity w at the lowest atmospheric layer [m/s]
192  real(RP), intent(in) :: ua (ia,ja) ! velocity u at the lowest atmospheric layer [m/s]
193  real(RP), intent(in) :: va (ia,ja) ! velocity v at the lowest atmospheric layer [m/s]
194  real(RP), intent(in) :: rhoa(ia,ja) ! density at the lowest atmospheric layer [kg/m3]
195  real(RP), intent(in) :: qva (ia,ja) ! ratio of water vapor mass to total mass at the lowest atmospheric layer [kg/kg]
196  real(RP), intent(in) :: z1 (ia,ja) ! cell center height at the lowest atmospheric layer [m]
197  real(RP), intent(in) :: pbl (ia,ja) ! the top of atmospheric mixing layer [m]
198  real(RP), intent(in) :: prss(ia,ja) ! pressure at the surface [Pa]
199  real(RP), intent(in) :: lwd (ia,ja) ! downward long-wave radiation flux at the surface [J/m2/s]
200  real(RP), intent(in) :: swd (ia,ja) ! downward short-wave radiation flux at the surface [J/m2/s]
201 
202  real(RP), intent(in) :: tg (ia,ja) ! soil temperature [K]
203  real(RP), intent(in) :: lst (ia,ja) ! land surface temperature [K]
204  real(RP), intent(in) :: qvef (ia,ja) ! efficiency of evaporation [0-1]
205  real(RP), intent(in) :: alb_lw(ia,ja) ! surface albedo for LW [0-1]
206  real(RP), intent(in) :: alb_sw(ia,ja) ! surface albedo for SW [0-1]
207  real(RP), intent(in) :: dzg (ia,ja) ! soil depth [m]
208  real(RP), intent(in) :: tcs (ia,ja) ! thermal conductivity for soil [J/m/K/s]
209  real(RP), intent(in) :: z0m (ia,ja) ! roughness length for momemtum [m]
210  real(RP), intent(in) :: z0h (ia,ja) ! roughness length for heat [m]
211  real(RP), intent(in) :: z0e (ia,ja) ! roughness length for vapor [m]
212  real(DP), intent(in) :: dt ! delta time
213 
214  ! works
215  real(RP) :: lst1(ia,ja)
216 
217  real(RP) :: res ! residual
218  real(RP) :: dres ! d(residual)/dLST
219  real(RP) :: oldres ! residual in previous step
220  real(RP) :: redf ! reduced factor
221 
222  real(RP) :: ustar, dustar ! friction velocity [m]
223  real(RP) :: tstar, dtstar ! friction potential temperature [K]
224  real(RP) :: qstar, dqstar ! friction water vapor mass ratio [kg/kg]
225  real(RP) :: uabs, duabs ! modified absolute velocity [m/s]
226  real(RP) :: sqv, dsqv ! saturation water vapor mixing ratio at surface [kg/kg]
227 
228  real(RP) :: lhv(ia,ja) ! latent heat for vaporization depending on temperature [J/kg]
229 
230  integer :: i, j, n
231  !---------------------------------------------------------------------------
232 
233  ! copy land surfce temperature for iteration
234  do j = js, je
235  do i = is, ie
236  lst1(i,j) = lst(i,j)
237  end do
238  end do
239 
240  call atmos_thermodyn_templhv( lhv, tmpa )
241 
242  ! update surface temperature
243  if( lst_update ) then
244 
245  do j = js, je
246  do i = is, ie
247 
248  if( is_lnd(i,j) ) then
249 
250  redf = 1.0_rp
251  oldres = huge(0.0_rp)
252 
253  ! modified Newton-Raphson method (Tomita 2009)
254  do n = 1, land_sfc_slab_itr_max
255 
256  call qsat( sqv, & ! [OUT]
257  lst1(i,j), & ! [IN]
258  prss(i,j) ) ! [IN]
259  call qsat( dsqv, & ! [OUT]
260  lst1(i,j)+dts0, & ! [IN]
261  prss(i,j) ) ! [IN]
262 
263  call bulkflux( &
264  ustar, & ! [OUT]
265  tstar, & ! [OUT]
266  qstar, & ! [OUT]
267  uabs, & ! [OUT]
268  tmpa(i,j), & ! [IN]
269  lst1(i,j), & ! [IN]
270  prsa(i,j), & ! [IN]
271  prss(i,j), & ! [IN]
272  qva(i,j), & ! [IN]
273  sqv, & ! [IN]
274  ua(i,j), & ! [IN]
275  va(i,j), & ! [IN]
276  z1(i,j), & ! [IN]
277  pbl(i,j), & ! [IN]
278  z0m(i,j), & ! [IN]
279  z0h(i,j), & ! [IN]
280  z0e(i,j) ) ! [IN]
281 
282  call bulkflux( &
283  dustar, & ! [OUT]
284  dtstar, & ! [OUT]
285  dqstar, & ! [OUT]
286  duabs, & ! [OUT]
287  tmpa(i,j), & ! [IN]
288  lst1(i,j)+dts0, & ! [IN]
289  prsa(i,j), & ! [IN]
290  prss(i,j), & ! [IN]
291  qva(i,j), & ! [IN]
292  dsqv, & ! [IN]
293  ua(i,j), & ! [IN]
294  va(i,j), & ! [IN]
295  z1(i,j), & ! [IN]
296  pbl(i,j), & ! [IN]
297  z0m(i,j), & ! [IN]
298  z0h(i,j), & ! [IN]
299  z0e(i,j) ) ! [IN]
300 
301  ! calculation for residual
302  res = ( 1.0_rp - alb_sw(i,j) ) * swd(i,j) &
303  + ( 1.0_rp - alb_lw(i,j) ) * ( lwd(i,j) - stb * lst1(i,j)**4 ) &
304  + cpdry * rhoa(i,j) * ustar * tstar &
305  + lhv(i,j) * rhoa(i,j) * ustar * qstar * qvef(i,j) &
306  - 2.0_rp * tcs(i,j) * ( lst1(i,j) - tg(i,j) ) / dzg(i,j)
307 
308  ! calculation for d(residual)/dLST
309  dres = -4.0_rp * ( 1.0_rp - alb_lw(i,j) ) * stb * lst1(i,j)**3 &
310  + cpdry * rhoa(i,j) * ( (dustar-ustar)/dts0 * tstar + ustar * (dtstar-tstar)/dts0 ) &
311  + lhv(i,j) * rhoa(i,j) * ( (dustar-ustar)/dts0 * qstar + ustar * (dqstar-qstar)/dts0 ) * qvef(i,j) &
312  - 2.0_rp * tcs(i,j) / dzg(i,j)
313 
314  ! convergence test with residual and error levels
315  if( abs( res ) < land_sfc_slab_res_min .OR. &
316  abs( res/dres ) < land_sfc_slab_err_min ) then
317  exit
318  end if
319 
320  ! stop iteration to prevent numerical error
321  if( abs(dres) * land_sfc_slab_dreslim < abs(res) ) then
322  exit
323  end if
324 
325  ! calculate reduced factor
326  if( dres < 0.0_rp ) then
327  if( abs(res) > abs(oldres) ) then
328  redf = max( tfa*abs(redf), redf_min )
329  else
330  redf = min( tfb*abs(redf), redf_max )
331  end if
332  else
333  redf = -1.0_rp
334  end if
335 
336  ! estimate next surface temperature
337  lst1(i,j) = lst1(i,j) - redf * res / dres
338 
339  ! save residual in this step
340  oldres = res
341 
342  end do
343 
344  ! update land surface temperature with limitation
345  lst1(i,j) = min( max( lst1(i,j), &
346  lst(i,j) - land_sfc_slab_dts_max * dt ), &
347  lst(i,j) + land_sfc_slab_dts_max * dt )
348 
349  if( n > land_sfc_slab_itr_max ) then
350  ! land surface temperature was not converged
351  if( io_l ) write(io_fid_log,'(A)' ) 'Warning: land surface tempearture was not converged.'
352  if( io_l ) write(io_fid_log,'(A)' ) ''
353  if( io_l ) write(io_fid_log,'(A,I32)' ) 'DEBUG --- PRC_myrank [no unit] :', prc_myrank
354  if( io_l ) write(io_fid_log,'(A,I32)' ) 'DEBUG --- number of i [no unit] :', i
355  if( io_l ) write(io_fid_log,'(A,I32)' ) 'DEBUG --- number of j [no unit] :', j
356  if( io_l ) write(io_fid_log,'(A)' ) ''
357  if( io_l ) write(io_fid_log,'(A,I32)' ) 'DEBUG --- loop number [no unit] :', n
358  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- Residual [J/m2/s] :', res
359  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- delta Residual [J/m2/s] :', dres
360  if( io_l ) write(io_fid_log,'(A,F32.16)') ''
361  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- temperature [K] :', tmpa(i,j)
362  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- pressure [Pa] :', prsa(i,j)
363  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- velocity w [m/s] :', wa(i,j)
364  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- velocity u [m/s] :', ua(i,j)
365  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- velocity v [m/s] :', va(i,j)
366  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- density [kg/m3] :', rhoa(i,j)
367  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- water vapor mass ratio [kg/kg] :', qva(i,j)
368  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- cell center height [m] :', z1(i,j)
369  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- atmospheric mixing layer height [m] :', pbl(i,j)
370  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- pressure at the surface [Pa] :', prss(i,j)
371  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- downward long-wave radiation [J/m2/s] :', lwd(i,j)
372  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- downward short-wave radiation [J/m2/s] :', swd(i,j)
373  if( io_l ) write(io_fid_log,'(A)' ) ''
374  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- soil temperature [K] :', tg(i,j)
375  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- land surface temperature [K] :', lst(i,j)
376  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- efficiency of evaporation [0-1] :', qvef(i,j)
377  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- surface albedo for LW [0-1] :', alb_lw(i,j)
378  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- surface albedo for SW [0-1] :', alb_sw(i,j)
379  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- soil depth [m] :', dzg(i,j)
380  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- thermal conductivity for soil [J/m/K/s] :', tcs(i,j)
381  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- roughness length for momemtum [m] :', z0m(i,j)
382  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- roughness length for heat [m] :', z0h(i,j)
383  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- roughness length for vapor [m] :', z0e(i,j)
384  if( io_l ) write(io_fid_log,'(A)' ) ''
385  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- latent heat [J/kg] :', lhv(i,j)
386  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- friction velocity [m] :', ustar
387  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- friction potential temperature [K] :', tstar
388  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- friction water vapor mass ratio [kg/kg] :', qstar
389  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- d(friction velocity) [m] :', dustar
390  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- d(friction potential temperature) [K] :', dtstar
391  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- d(friction water vapor mass ratio) [kg/kg] :', dqstar
392  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- modified absolute velocity [m/s] :', uabs
393  if( io_l ) write(io_fid_log,'(A,F32.16)') 'DEBUG --- next land surface temperature [K] :', lst1(i,j)
394 
395  ! check NaN
396  if ( .NOT. ( res > -1.0_rp .OR. res < 1.0_rp ) ) then ! must be NaN
397  write(*,*) 'xxx NaN is detected for land surface temperature.'
398  write(*,*) ''
399  write(*,*) 'DEBUG --- PRC_myrank :', prc_myrank
400  write(*,*) 'DEBUG --- number of i :', i
401  write(*,*) 'DEBUG --- number of j :', j
402  write(*,*) ''
403  write(*,*) 'DEBUG --- Residual [J/m2/s] :', res
404  write(*,*) 'DEBUG --- delta Residual [J/m2/s] :', dres
405  write(*,*) ''
406  write(*,*) 'DEBUG --- temperature [K] :', tmpa(i,j)
407  write(*,*) 'DEBUG --- pressure [Pa] :', prsa(i,j)
408  write(*,*) 'DEBUG --- velocity w [m/s] :', wa(i,j)
409  write(*,*) 'DEBUG --- velocity u [m/s] :', ua(i,j)
410  write(*,*) 'DEBUG --- velocity v [m/s] :', va(i,j)
411  write(*,*) 'DEBUG --- density [kg/m3] :', rhoa(i,j)
412  write(*,*) 'DEBUG --- water vapor mass ratio [kg/kg] :', qva(i,j)
413  write(*,*) 'DEBUG --- cell center height [m] :', z1(i,j)
414  write(*,*) 'DEBUG --- atmospheric mixing layer height [m] :', pbl(i,j)
415  write(*,*) 'DEBUG --- pressure at the surface [Pa] :', prss(i,j)
416  write(*,*) 'DEBUG --- downward long-wave radiation [J/m2/s] :', lwd(i,j)
417  write(*,*) 'DEBUG --- downward short-wave radiation [J/m2/s] :', swd(i,j)
418  write(*,*) ''
419  write(*,*) 'DEBUG --- soil temperature [K] :', tg(i,j)
420  write(*,*) 'DEBUG --- land surface temperature [K] :', lst(i,j)
421  write(*,*) 'DEBUG --- efficiency of evaporation [0-1] :', qvef(i,j)
422  write(*,*) 'DEBUG --- surface albedo for LW [0-1] :', alb_lw(i,j)
423  write(*,*) 'DEBUG --- surface albedo for SW [0-1] :', alb_sw(i,j)
424  write(*,*) 'DEBUG --- soil depth [m] :', dzg(i,j)
425  write(*,*) 'DEBUG --- thermal conductivity for soil [J/m/K/s] :', tcs(i,j)
426  write(*,*) 'DEBUG --- roughness length for momemtum [m] :', z0m(i,j)
427  write(*,*) 'DEBUG --- roughness length for heat [m] :', z0h(i,j)
428  write(*,*) 'DEBUG --- roughness length for vapor [m] :', z0e(i,j)
429  write(*,*) ''
430  write(*,*) 'DEBUG --- latent heat [J/kg] :', lhv(i,j)
431  write(*,*) 'DEBUG --- friction velocity [m] :', ustar
432  write(*,*) 'DEBUG --- friction potential temperature [K] :', tstar
433  write(*,*) 'DEBUG --- friction water vapor mass ratio [kg/kg] :', qstar
434  write(*,*) 'DEBUG --- d(friction velocity) [m] :', dustar
435  write(*,*) 'DEBUG --- d(friction potential temperature) [K] :', dtstar
436  write(*,*) 'DEBUG --- d(friction water vapor mass ratio) [kg/kg] :', dqstar
437  write(*,*) 'DEBUG --- modified absolute velocity [m/s] :', uabs
438  write(*,*) 'DEBUG --- next land surface temperature [K] :', lst1(i,j)
439 
440  call prc_mpistop
441  endif
442 
443  end if
444 
445  end if
446 
447  ! calculate tendency
448  lst_t(i,j) = ( lst1(i,j) - lst(i,j) ) / dt
449 
450  end do
451  end do
452 
453  ! not update temperature
454  else
455 
456  do j = js, je
457  do i = is, ie
458  lst_t(i,j) = 0.0_rp
459  end do
460  end do
461 
462  end if
463 
464 
465  ! calculate surface flux
466  do j = js, je
467  do i = is, ie
468 
469  if( is_lnd(i,j) ) then
470 
471  call qsat( sqv, & ! [OUT]
472  lst1(i,j), & ! [IN]
473  prss(i,j) ) ! [IN]
474 
475  call bulkflux( &
476  ustar, & ! [OUT]
477  tstar, & ! [OUT]
478  qstar, & ! [OUT]
479  uabs, & ! [OUT]
480  tmpa(i,j), & ! [IN]
481  lst1(i,j), & ! [IN]
482  prsa(i,j), & ! [IN]
483  prss(i,j), & ! [IN]
484  qva(i,j), & ! [IN]
485  sqv, & ! [IN]
486  ua(i,j), & ! [IN]
487  va(i,j), & ! [IN]
488  z1(i,j), & ! [IN]
489  pbl(i,j), & ! [IN]
490  z0m(i,j), & ! [IN]
491  z0h(i,j), & ! [IN]
492  z0e(i,j) ) ! [IN]
493 
494  zmflx(i,j) = -rhoa(i,j) * ustar**2 / uabs * wa(i,j)
495  xmflx(i,j) = -rhoa(i,j) * ustar**2 / uabs * ua(i,j)
496  ymflx(i,j) = -rhoa(i,j) * ustar**2 / uabs * va(i,j)
497  shflx(i,j) = -cpdry * rhoa(i,j) * ustar * tstar
498  lhflx(i,j) = -lhv(i,j) * rhoa(i,j) * ustar * qstar * qvef(i,j)
499  ghflx(i,j) = -2.0_rp * tcs(i,j) * ( lst1(i,j) - tg(i,j) ) / dzg(i,j)
500 
501  ! calculation for residual
502  res = ( 1.0_rp - alb_sw(i,j) ) * swd(i,j) &
503  + ( 1.0_rp - alb_lw(i,j) ) * ( lwd(i,j) - stb * lst1(i,j)**4 ) &
504  - shflx(i,j) - lhflx(i,j) + ghflx(i,j)
505 
506  ! put residual in ground heat flux
507  ghflx(i,j) = ghflx(i,j) - res
508 
509  ! diagnostic variables
510  u10(i,j) = ua(i,j) * log( 10.0_rp / z0m(i,j) ) / log( z1(i,j) / z0m(i,j) )
511  v10(i,j) = va(i,j) * log( 10.0_rp / z0m(i,j) ) / log( z1(i,j) / z0m(i,j) )
512  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) ) ) &
513  / ( log( z1(i,j) / z0m(i,j) ) * log( z1(i,j) / z0h(i,j) ) )
514  q2(i,j) = sqv + ( qva(i,j) - sqv ) * ( log( 2.0_rp / z0m(i,j) ) * log( 2.0_rp / z0e(i,j) ) ) &
515  / ( log( z1(i,j) / z0m(i,j) ) * log( z1(i,j) / z0e(i,j) ) )
516  else
517 
518  ! not calculate surface flux
519  zmflx(i,j) = 0.0_rp
520  xmflx(i,j) = 0.0_rp
521  ymflx(i,j) = 0.0_rp
522  shflx(i,j) = 0.0_rp
523  lhflx(i,j) = 0.0_rp
524  ghflx(i,j) = 0.0_rp
525  u10(i,j) = 0.0_rp
526  v10(i,j) = 0.0_rp
527  t2(i,j) = 0.0_rp
528  q2(i,j) = 0.0_rp
529 
530  end if
531 
532  end do
533  end do
534 
535  return
integer, public is
start point of inner domain: x, local
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.
integer, public va
Definition: scale_index.F90:38
real(rp), parameter, public const_stb
Stefan-Boltzman constant [W/m2/K4].
Definition: scale_const.F90:51
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
module grid index
integer, public ia
of x whole cells (local, with HALO)
procedure(bc), pointer, public bulkflux
integer, public js
start point of inner domain: y, local
module PROCESS
subroutine, public log(type, message)
Definition: dc_log.f90:133
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
real(rp), public const_eps
small number
Definition: scale_const.F90:36
module ATMOSPHERE / Thermodynamics
integer, public ja
of y whole cells (local, with HALO)
Here is the call graph for this function:
Here is the caller graph for this function: