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, Z0M, Z0H, Z0E)
 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,
real(rp), dimension(ia,ja), intent(out)  Z0M,
real(rp), dimension(ia,ja), intent(out)  Z0H,
real(rp), dimension(ia,ja), intent(out)  Z0E 
)

Setup.

Definition at line 119 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().

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

Referenced by scale_urban_phy::urban_phy_setup().

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