SCALE-RM
Functions/Subroutines
scale_urban_phy_slc Module Reference

module URBAN / Surface fluxes with Single-layer Canpoy Model More...

Functions/Subroutines

subroutine, public urban_phy_slc_setup (URBAN_TYPE)
 Setup. More...
 
subroutine, public urban_phy_slc (TR_URB_t, TB_URB_t, TG_URB_t, TC_URB_t, QC_URB_t, UC_URB_t, TRL_URB_t, TBL_URB_t, TGL_URB_t, RAINR_URB_t, RAINB_URB_t, RAING_URB_t, ROFF_URB_t, SFC_TEMP, ALBD_LW, ALBD_SW, MWFLX, MUFLX, MVFLX, SHFLX, LHFLX, GHFLX, Z0M, Z0H, Z0E, U10, V10, T2, Q2, TMPA, PRSA, W1, U1, V1, DENS, QA, Z1, PBL, PRSS, LWD, SWD, PREC, TR_URB, TB_URB, TG_URB, TC_URB, QC_URB, UC_URB, TRL_URB, TBL_URB, TGL_URB, RAINR_URB, RAINB_URB, RAING_URB, ROFF_URB, LON, LAT, NOWDATE, dt)
 

Detailed Description

module URBAN / Surface fluxes with Single-layer Canpoy Model

Description
Surface fluxes between atmosphere and urban based on Single-layer Urban Canopy Model (Kusaka et al. 2001, BLM)
Author
Team SCALE
History
NAMELIST
  • PARAM_URBAN_PHY_SLC
    nametypedefault valuecomment
    DTS_MAX real(RP) 0.1_RP maximum dT during one minute [K/sec]
    ZR real(RP) 10.0_RP roof level ( building height) [m]
    ROOF_WIDTH real(RP) 9.0_RP roof level ( building height) [m]
    ROAD_WIDTH real(RP) 11.0_RP roof level ( building height) [m]
    SIGMA_ZED real(RP) 1.0_RP Standard deviation of roof height [m]
    AH real(RP) 17.5_RP Sensible Anthropogenic heat [W/m^2]
    ALH real(RP) 0.0_RP Latent Anthropogenic heat [W/m^2]
    BETR real(RP) 0.0_RP Evaporation efficiency of roof [-]
    BETB real(RP) 0.0_RP of building [-]
    BETG real(RP) 0.0_RP of ground [-]
    STRGR real(RP) 0.0_RP rain strage on roof [-]
    STRGB real(RP) 0.0_RP on wall [-]
    STRGG real(RP) 0.0_RP on ground [-]
    CAPR real(RP) 1.2E6_RP heat capacity of roof [J m-3 K]
    CAPB real(RP) 1.2E6_RP of wall [J m-3 K]
    CAPG real(RP) 1.2E6_RP of ground [J m-3 K]
    AKSR real(RP) 2.28_RP thermal conductivity of roof [W m-1 K]
    AKSB real(RP) 2.28_RP of wall [W m-1 K]
    AKSG real(RP) 2.28_RP of ground [W m-1 K]
    ALBR real(RP) 0.2_RP surface albedo of roof
    ALBB real(RP) 0.2_RP surface albedo of wall
    ALBG real(RP) 0.2_RP surface albedo of ground
    EPSR real(RP) 0.90_RP Surface emissivity of roof
    EPSB real(RP) 0.90_RP Surface emissivity of wall
    EPSG real(RP) 0.90_RP Surface emissivity of ground
    Z0R real(RP) 0.01_RP roughness length for momentum of building roof
    Z0B real(RP) 0.0001_RP roughness length for momentum of building wall
    Z0G real(RP) 0.01_RP roughness length for momentum of ground
    TRLEND real(RP) 293.00_RP lower boundary condition of roof temperature [K]
    TBLEND real(RP) 293.00_RP lower boundary condition of wall temperature [K]
    TGLEND real(RP) 293.00_RP lower boundary condition of ground temperature [K]
    BOUND integer

History Output
namedescriptionunitvariable
URBAN_GHB urban ground heat flux on wall W/m2 GHB
URBAN_GHG urban ground heat flux on road W/m2 GHG
URBAN_GHR urban ground heat flux on roof W/m2 GHR
URBAN_LHB urban latent heat flux on wall W/m2 LHB
URBAN_LHG urban latent heat flux on road W/m2 LHG
URBAN_LHR urban latent heat flux on roof W/m2 LHR
URBAN_RNB urban net radiation on wall W/m2 RNB
URBAN_RNG urban net radiation on road W/m2 RNG
URBAN_RNR urban net radiation on roof W/m2 RNR
URBAN_RNgrd urban grid average of net radiation W/m2 RNgrd
URBAN_SHB urban sensible heat flux on wall W/m2 SHB
URBAN_SHG urban sensible heat flux on road W/m2 SHG
URBAN_SHR urban sensible heat flux on roof W/m2 SHR

Function/Subroutine Documentation

◆ urban_phy_slc_setup()

subroutine, public scale_urban_phy_slc::urban_phy_slc_setup ( character(len=*), intent(in)  URBAN_TYPE)

Setup.

Definition at line 115 of file scale_urban_phy_slc.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_urban, scale_process::prc_mpistop(), scale_urban_grid_index::uke, and scale_urban_grid_index::uks.

Referenced by scale_urban_phy::urban_phy_setup().

115  use scale_process, only: &
117  use scale_landuse, only: &
119  implicit none
120 
121  character(len=*), intent(in) :: urban_type
122  integer :: i, j
123 
124  namelist / param_urban_phy_slc / &
125  dts_max, &
126  zr, &
127  roof_width, &
128  road_width, &
129  sigma_zed, &
130  ah, &
131  alh, &
132  betr, &
133  betb, &
134  betg, &
135  strgr, &
136  strgb, &
137  strgg, &
138  capr, &
139  capb, &
140  capg, &
141  aksr, &
142  aksb, &
143  aksg, &
144  albr, &
145  albb, &
146  albg, &
147  epsr, &
148  epsb, &
149  epsg, &
150  z0r, &
151  z0b, &
152  z0g, &
153  trlend, &
154  tblend, &
155  tglend, &
156  bound
157 
158  integer :: ierr
159  !---------------------------------------------------------------------------
160 
161  if( io_l ) write(io_fid_log,*)
162  if( io_l ) write(io_fid_log,*) '++++++ Module[SLC] / Categ[URBAN PHY] / Origin[SCALElib]'
163 
164  allocate( dzr(uks:uke) )
165  allocate( dzb(uks:uke) )
166  allocate( dzg(uks:uke) )
167  dzr(uks:uke) = (/0.01_rp,0.01_rp,0.03_rp,0.05_rp,0.10_rp/)
168  dzb(uks:uke) = (/0.01_rp,0.01_rp,0.03_rp,0.05_rp,0.10_rp/)
169  dzg(uks:uke) = (/0.01_rp,0.01_rp,0.03_rp,0.05_rp,0.10_rp/)
170 
171  !--- read namelist
172  rewind(io_fid_conf)
173  read(io_fid_conf,nml=param_urban_phy_slc,iostat=ierr)
174  if( ierr < 0 ) then !--- missing
175  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
176  elseif( ierr > 0 ) then !--- fatal error
177  write(*,*) 'xxx Not appropriate names in namelist PARAM_URBAN_PHY_SLC. Check!'
178  call prc_mpistop
179  endif
180  if( io_lnml ) write(io_fid_log,nml=param_urban_phy_slc)
181 
182  ahdiurnal(:) = (/ 0.356, 0.274, 0.232, 0.251, 0.375, 0.647, 0.919, 1.135, 1.249, 1.328, &
183  1.365, 1.363, 1.375, 1.404, 1.457, 1.526, 1.557, 1.521, 1.372, 1.206, &
184  1.017, 0.876, 0.684, 0.512 /)
185 
186  ! set other urban parameters
187  call urban_param_setup
188 
189  ! judge to run slab land model
190  allocate( is_urb(ia,ja) )
191 
192  do j = js, je
193  do i = is, ie
194  if( landuse_fact_urban(i,j) > 0.0_rp )then
195  is_urb(i,j) = .true.
196  else
197  is_urb(i,j) = .false.
198  endif
199  enddo
200  enddo
201 
202  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.
real(rp), dimension(:,:), allocatable, public landuse_fact_urban
urban factor
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
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:

◆ urban_phy_slc()

subroutine, public scale_urban_phy_slc::urban_phy_slc ( real(rp), dimension (ia,ja), intent(out)  TR_URB_t,
real(rp), dimension (ia,ja), intent(out)  TB_URB_t,
real(rp), dimension (ia,ja), intent(out)  TG_URB_t,
real(rp), dimension (ia,ja), intent(out)  TC_URB_t,
real(rp), dimension (ia,ja), intent(out)  QC_URB_t,
real(rp), dimension (ia,ja), intent(out)  UC_URB_t,
real(rp), dimension (uks:uke,ia,ja), intent(out)  TRL_URB_t,
real(rp), dimension (uks:uke,ia,ja), intent(out)  TBL_URB_t,
real(rp), dimension (uks:uke,ia,ja), intent(out)  TGL_URB_t,
real(rp), dimension(ia,ja), intent(out)  RAINR_URB_t,
real(rp), dimension(ia,ja), intent(out)  RAINB_URB_t,
real(rp), dimension(ia,ja), intent(out)  RAING_URB_t,
real(rp), dimension (ia,ja), intent(out)  ROFF_URB_t,
real(rp), dimension(ia,ja), intent(out)  SFC_TEMP,
real(rp), dimension (ia,ja), intent(out)  ALBD_LW,
real(rp), dimension (ia,ja), intent(out)  ALBD_SW,
real(rp), dimension (ia,ja), intent(out)  MWFLX,
real(rp), dimension (ia,ja), intent(out)  MUFLX,
real(rp), dimension (ia,ja), intent(out)  MVFLX,
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)  Z0M,
real(rp), dimension (ia,ja), intent(out)  Z0H,
real(rp), dimension (ia,ja), intent(out)  Z0E,
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)  W1,
real(rp), dimension (ia,ja), intent(in)  U1,
real(rp), dimension (ia,ja), intent(in)  V1,
real(rp), dimension (ia,ja), intent(in)  DENS,
real(rp), dimension (ia,ja), intent(in)  QA,
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,2), intent(in)  LWD,
real(rp), dimension (ia,ja,2), intent(in)  SWD,
real(rp), dimension (ia,ja), intent(in)  PREC,
real(rp), dimension (ia,ja), intent(in)  TR_URB,
real(rp), dimension (ia,ja), intent(in)  TB_URB,
real(rp), dimension (ia,ja), intent(in)  TG_URB,
real(rp), dimension (ia,ja), intent(in)  TC_URB,
real(rp), dimension (ia,ja), intent(in)  QC_URB,
real(rp), dimension (ia,ja), intent(in)  UC_URB,
real(rp), dimension (uks:uke,ia,ja), intent(in)  TRL_URB,
real(rp), dimension (uks:uke,ia,ja), intent(in)  TBL_URB,
real(rp), dimension (uks:uke,ia,ja), intent(in)  TGL_URB,
real(rp), dimension(ia,ja), intent(in)  RAINR_URB,
real(rp), dimension(ia,ja), intent(in)  RAINB_URB,
real(rp), dimension(ia,ja), intent(in)  RAING_URB,
real(rp), dimension (ia,ja), intent(in)  ROFF_URB,
real(rp), intent(in)  LON,
real(rp), intent(in)  LAT,
integer, dimension(6), intent(in)  NOWDATE,
real(dp), intent(in)  dt 
)

Definition at line 265 of file scale_urban_phy_slc.F90.

References scale_bulkflux::bulkflux, scale_const::const_cpdry, scale_const::const_eps, scale_const::const_pi, 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(), scale_process::prc_myrank, scale_tracer::qa, scale_urban_grid_index::uke, and scale_urban_grid_index::uks.

Referenced by scale_urban_phy::urban_phy_setup().

265  use scale_grid_index
267  use scale_history, only: &
268  hist_in
269  use scale_atmos_saturation, only: &
270  qsat => atmos_saturation_pres2qsat_all
271  use scale_bulkflux, only: &
272  bulkflux
273  implicit none
274 
275  ! parameter
276  logical, parameter :: lsolar = .false. ! [true=both, false=SSG only]
277 
278  real(RP), parameter :: uabs_min = 0.1_rp
279 
280  ! arguments
281  real(RP), intent(out) :: tr_urb_t (ia,ja)
282  real(RP), intent(out) :: tb_urb_t (ia,ja)
283  real(RP), intent(out) :: tg_urb_t (ia,ja)
284  real(RP), intent(out) :: tc_urb_t (ia,ja)
285  real(RP), intent(out) :: qc_urb_t (ia,ja)
286  real(RP), intent(out) :: uc_urb_t (ia,ja)
287  real(RP), intent(out) :: trl_urb_t (uks:uke,ia,ja)
288  real(RP), intent(out) :: tbl_urb_t (uks:uke,ia,ja)
289  real(RP), intent(out) :: tgl_urb_t (uks:uke,ia,ja)
290  real(RP), intent(out) :: rainr_urb_t(ia,ja)
291  real(RP), intent(out) :: rainb_urb_t(ia,ja)
292  real(RP), intent(out) :: raing_urb_t(ia,ja)
293  real(RP), intent(out) :: roff_urb_t (ia,ja)
294 
295  real(RP), intent(out) :: sfc_temp(ia,ja)
296  real(RP), intent(out) :: albd_lw (ia,ja)
297  real(RP), intent(out) :: albd_sw (ia,ja)
298  real(RP), intent(out) :: mwflx (ia,ja)
299  real(RP), intent(out) :: muflx (ia,ja)
300  real(RP), intent(out) :: mvflx (ia,ja)
301  real(RP), intent(out) :: shflx (ia,ja)
302  real(RP), intent(out) :: lhflx (ia,ja)
303  real(RP), intent(out) :: ghflx (ia,ja)
304  real(RP), intent(out) :: z0m (ia,ja)
305  real(RP), intent(out) :: z0h (ia,ja)
306  real(RP), intent(out) :: z0e (ia,ja)
307  real(RP), intent(out) :: u10 (ia,ja)
308  real(RP), intent(out) :: v10 (ia,ja)
309  real(RP), intent(out) :: t2 (ia,ja)
310  real(RP), intent(out) :: q2 (ia,ja)
311 
312  real(RP), intent(in) :: tmpa (ia,ja)
313  real(RP), intent(in) :: prsa (ia,ja)
314  real(RP), intent(in) :: w1 (ia,ja)
315  real(RP), intent(in) :: u1 (ia,ja)
316  real(RP), intent(in) :: v1 (ia,ja)
317  real(RP), intent(in) :: dens (ia,ja)
318  real(RP), intent(in) :: qa (ia,ja)
319  real(RP), intent(in) :: z1 (ia,ja)
320  real(RP), intent(in) :: pbl (ia,ja)
321  real(RP), intent(in) :: prss (ia,ja)
322  real(RP), intent(in) :: lwd (ia,ja,2)
323  real(RP), intent(in) :: swd (ia,ja,2)
324  real(RP), intent(in) :: prec (ia,ja)
325 
326  real(RP), intent(in) :: tr_urb (ia,ja)
327  real(RP), intent(in) :: tb_urb (ia,ja)
328  real(RP), intent(in) :: tg_urb (ia,ja)
329  real(RP), intent(in) :: tc_urb (ia,ja)
330  real(RP), intent(in) :: qc_urb (ia,ja)
331  real(RP), intent(in) :: uc_urb (ia,ja)
332  real(RP), intent(in) :: trl_urb (uks:uke,ia,ja)
333  real(RP), intent(in) :: tbl_urb (uks:uke,ia,ja)
334  real(RP), intent(in) :: tgl_urb (uks:uke,ia,ja)
335  real(RP), intent(in) :: rainr_urb(ia,ja)
336  real(RP), intent(in) :: rainb_urb(ia,ja)
337  real(RP), intent(in) :: raing_urb(ia,ja)
338  real(RP), intent(in) :: roff_urb (ia,ja)
339 
340  real(RP), intent(in) :: lon
341  real(RP), intent(in) :: lat
342  integer, intent(in) :: nowdate(6)
343  real(DP), intent(in) :: dt
344 
345  ! work
346  real(RP) :: tr
347  real(RP) :: tb
348  real(RP) :: tg
349  real(RP) :: tc
350  real(RP) :: qc
351  real(RP) :: uc
352  real(RP) :: trl(uks:uke)
353  real(RP) :: tbl(uks:uke)
354  real(RP) :: tgl(uks:uke)
355  real(RP) :: rainr
356  real(RP) :: rainb
357  real(RP) :: raing
358  real(RP) :: roff
359 
360  real(RP) :: shr (ia,ja)
361  real(RP) :: shb (ia,ja)
362  real(RP) :: shg (ia,ja)
363  real(RP) :: lhr (ia,ja)
364  real(RP) :: lhb (ia,ja)
365  real(RP) :: lhg (ia,ja)
366  real(RP) :: ghr (ia,ja)
367  real(RP) :: ghb (ia,ja)
368  real(RP) :: ghg (ia,ja)
369  real(RP) :: rnr (ia,ja)
370  real(RP) :: rnb (ia,ja)
371  real(RP) :: rng (ia,ja)
372  real(RP) :: rngrd(ia,ja)
373 
374  real(RP) :: ustar ! friction velocity [m]
375  real(RP) :: tstar ! friction temperature [K]
376  real(RP) :: qstar ! friction mixing rate [kg/kg]
377  real(RP) :: uabs ! modified absolute velocity [m/s]
378 
379  real(RP) :: qvs ! saturation water vapor mixing ratio at surface [kg/kg]
380 
381  integer :: k, i, j
382  !---------------------------------------------------------------------------
383 
384  if( io_l ) write(io_fid_log,*) '*** Urban step: Single Layer Canopy'
385 
386 
387  do j = js, je
388  do i = is, ie
389 
390  if( is_urb(i,j) ) then
391 
392  uabs = max( sqrt( u1(i,j)**2 + v1(i,j)**2 + w1(i,j)**2 ), uabs_min )
393 
394  ! save
395  tr = tr_urb(i,j)
396  tb = tb_urb(i,j)
397  tg = tg_urb(i,j)
398  tc = tc_urb(i,j)
399  qc = qc_urb(i,j)
400  uc = uc_urb(i,j)
401 
402  do k = uks, uke
403  trl(k) = trl_urb(k,i,j)
404  tbl(k) = tbl_urb(k,i,j)
405  tgl(k) = tgl_urb(k,i,j)
406  end do
407 
408  rainr = rainr_urb(i,j)
409  rainb = rainb_urb(i,j)
410  raing = raing_urb(i,j)
411  roff = roff_urb(i,j)
412 
413  call slc_main( tr, & ! [INOUT]
414  tb, & ! [INOUT]
415  tg, & ! [INOUT]
416  tc, & ! [INOUT]
417  qc, & ! [INOUT]
418  uc, & ! [INOUT]
419  trl(:), & ! [INOUT]
420  tbl(:), & ! [INOUT]
421  tgl(:), & ! [INOUT]
422  rainr, & ! [INOUT]
423  rainb, & ! [INOUT]
424  raing, & ! [INOUT]
425  roff, & ! [INOUT]
426  albd_lw(i,j), & ! [OUT]
427  albd_sw(i,j), & ! [OUT]
428  shr(i,j), & ! [OUT]
429  shb(i,j), & ! [OUT]
430  shg(i,j), & ! [OUT]
431  lhr(i,j), & ! [OUT]
432  lhb(i,j), & ! [OUT]
433  lhg(i,j), & ! [OUT]
434  ghr(i,j), & ! [OUT]
435  ghb(i,j), & ! [OUT]
436  ghg(i,j), & ! [OUT]
437  rnr(i,j), & ! [OUT]
438  rnb(i,j), & ! [OUT]
439  rng(i,j), & ! [OUT]
440  sfc_temp(i,j), & ! [OUT]
441  rngrd(i,j), & ! [OUT]
442  shflx(i,j), & ! [OUT]
443  lhflx(i,j), & ! [OUT]
444  ghflx(i,j), & ! [OUT]
445  u10(i,j), & ! [OUT]
446  v10(i,j), & ! [OUT]
447  t2(i,j), & ! [OUT]
448  q2(i,j), & ! [OUT]
449  lsolar, & ! [IN]
450  prsa(i,j), & ! [IN]
451  prss(i,j), & ! [IN]
452  tmpa(i,j), & ! [IN]
453  qa(i,j), & ! [IN]
454  uabs, & ! [IN]
455  u1(i,j), & ! [IN]
456  v1(i,j), & ! [IN]
457  z1(i,j), & ! [IN]
458  swd(i,j,:), & ! [IN]
459  lwd(i,j,:), & ! [IN]
460  prec(i,j), & ! [IN]
461  dens(i,j), & ! [IN]
462  lon, & ! [IN]
463  lat, & ! [IN]
464  nowdate(:), & ! [IN]
465  dt, i, j ) ! [IN]
466 
467  ! calculate tendency
468  tr_urb_t(i,j) = ( tr - tr_urb(i,j) ) / dt
469  tb_urb_t(i,j) = ( tb - tb_urb(i,j) ) / dt
470  tg_urb_t(i,j) = ( tg - tg_urb(i,j) ) / dt
471  tc_urb_t(i,j) = ( tc - tc_urb(i,j) ) / dt
472  qc_urb_t(i,j) = ( qc - qc_urb(i,j) ) / dt
473  uc_urb_t(i,j) = ( uc - uc_urb(i,j) ) / dt
474 
475  do k = uks, uke
476  trl_urb_t(k,i,j) = ( trl(k) - trl_urb(k,i,j) ) / dt
477  tbl_urb_t(k,i,j) = ( tbl(k) - tbl_urb(k,i,j) ) / dt
478  tgl_urb_t(k,i,j) = ( tgl(k) - tgl_urb(k,i,j) ) / dt
479  end do
480 
481  rainr_urb_t(i,j) = ( rainr - rainr_urb(i,j) ) / dt
482  rainb_urb_t(i,j) = ( rainb - rainb_urb(i,j) ) / dt
483  raing_urb_t(i,j) = ( raing - raing_urb(i,j) ) / dt
484  roff_urb_t(i,j) = ( roff - roff_urb(i,j) ) / dt
485 
486  ! saturation at the surface
487  call qsat( qvs, sfc_temp(i,j), prss(i,j) )
488 
489  call bulkflux( ustar, & ! [OUT]
490  tstar, & ! [OUT]
491  qstar, & ! [OUT]
492  uabs, & ! [OUT]
493  tmpa(i,j), & ! [IN]
494  sfc_temp(i,j), & ! [IN]
495  prsa(i,j), & ! [IN]
496  prss(i,j), & ! [IN]
497  qa(i,j), & ! [IN]
498  qvs, & ! [IN]
499  u1(i,j), & ! [IN]
500  v1(i,j), & ! [IN]
501  z1(i,j), & ! [IN]
502  pbl(i,j), & ! [IN]
503  z0c, & ! [IN]
504  z0hc, & ! [IN]
505  z0hc ) ! [IN]
506 
507  mwflx(i,j) = -dens(i,j) * ustar**2 / uabs * w1(i,j)
508  muflx(i,j) = -dens(i,j) * ustar**2 / uabs * u1(i,j)
509  mvflx(i,j) = -dens(i,j) * ustar**2 / uabs * v1(i,j)
510 
511  z0m(i,j) = z0c
512  z0h(i,j) = z0hc
513  z0e(i,j) = z0hc
514 
515  else
516  ! not calculate urban module
517  tr_urb_t(i,j) = 0.0_rp
518  tb_urb_t(i,j) = 0.0_rp
519  tg_urb_t(i,j) = 0.0_rp
520  tc_urb_t(i,j) = 0.0_rp
521  qc_urb_t(i,j) = 0.0_rp
522  uc_urb_t(i,j) = 0.0_rp
523  trl_urb_t(:,i,j) = 0.0_rp
524  tbl_urb_t(:,i,j) = 0.0_rp
525  tgl_urb_t(:,i,j) = 0.0_rp
526  rainr_urb_t(i,j) = 0.0_rp
527  rainb_urb_t(i,j) = 0.0_rp
528  raing_urb_t(i,j) = 0.0_rp
529  roff_urb_t(i,j) = 0.0_rp
530  sfc_temp(i,j) = 300.0_rp ! constant value
531  albd_lw(i,j) = 0.0_rp
532  albd_sw(i,j) = 0.0_rp
533  mwflx(i,j) = 0.0_rp
534  muflx(i,j) = 0.0_rp
535  mvflx(i,j) = 0.0_rp
536  shflx(i,j) = 0.0_rp
537  lhflx(i,j) = 0.0_rp
538  ghflx(i,j) = 0.0_rp
539  z0m(i,j) = 0.0_rp
540  z0h(i,j) = 0.0_rp
541  z0e(i,j) = 0.0_rp
542  u10(i,j) = 0.0_rp
543  v10(i,j) = 0.0_rp
544  t2(i,j) = 0.0_rp
545  q2(i,j) = 0.0_rp
546  !
547  shr(i,j) = 0.0_rp
548  shb(i,j) = 0.0_rp
549  shg(i,j) = 0.0_rp
550  lhr(i,j) = 0.0_rp
551  lhb(i,j) = 0.0_rp
552  lhg(i,j) = 0.0_rp
553  ghr(i,j) = 0.0_rp
554  ghb(i,j) = 0.0_rp
555  ghg(i,j) = 0.0_rp
556  rnr(i,j) = 0.0_rp
557  rnb(i,j) = 0.0_rp
558  rng(i,j) = 0.0_rp
559  rngrd(i,j) = 0.0_rp
560  endif
561 
562  end do
563  end do
564 
565  call hist_in( shr(:,:), 'URBAN_SHR', 'urban sensible heat flux on roof', 'W/m2' )
566  call hist_in( shb(:,:), 'URBAN_SHB', 'urban sensible heat flux on wall', 'W/m2' )
567  call hist_in( shg(:,:), 'URBAN_SHG', 'urban sensible heat flux on road', 'W/m2' )
568  call hist_in( lhr(:,:), 'URBAN_LHR', 'urban latent heat flux on roof', 'W/m2' )
569  call hist_in( lhb(:,:), 'URBAN_LHB', 'urban latent heat flux on wall', 'W/m2' )
570  call hist_in( lhg(:,:), 'URBAN_LHG', 'urban latent heat flux on road', 'W/m2' )
571  call hist_in( ghr(:,:), 'URBAN_GHR', 'urban ground heat flux on roof', 'W/m2' )
572  call hist_in( ghb(:,:), 'URBAN_GHB', 'urban ground heat flux on wall', 'W/m2' )
573  call hist_in( ghg(:,:), 'URBAN_GHG', 'urban ground heat flux on road', 'W/m2' )
574  call hist_in( rnr(:,:), 'URBAN_RNR', 'urban net radiation on roof', 'W/m2' )
575  call hist_in( rnb(:,:), 'URBAN_RNB', 'urban net radiation on wall', 'W/m2' )
576  call hist_in( rng(:,:), 'URBAN_RNG', 'urban net radiation on road', 'W/m2' )
577  call hist_in( rngrd(:,:), 'URBAN_RNgrd', 'urban grid average of net radiation', 'W/m2' )
578 
579  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
module ATMOSPHERE / Saturation adjustment
integer, public qa
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
integer, public ie
end point of inner domain: x, local
module Surface bulk flux
module HISTORY
module urban grid index
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: