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)
 Main routine for land submodel. More...
 

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

Main routine for land submodel.

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

290  use scale_grid_index
292  use scale_history, only: &
293  hist_in
294  use scale_atmos_saturation, only: &
295  qsat => atmos_saturation_pres2qsat_all
296  use scale_bulkflux, only: &
297  bulkflux
298  implicit none
299 
300  ! parameter
301  logical, parameter :: LSOLAR = .false. ! [true=both, false=SSG only]
302 
303  real(RP), parameter :: Uabs_min = 0.1_rp
304 
305  ! arguments
306  real(RP), intent(out) :: TR_URB_t (IA,JA)
307  real(RP), intent(out) :: TB_URB_t (IA,JA)
308  real(RP), intent(out) :: TG_URB_t (IA,JA)
309  real(RP), intent(out) :: TC_URB_t (IA,JA)
310  real(RP), intent(out) :: QC_URB_t (IA,JA)
311  real(RP), intent(out) :: UC_URB_t (IA,JA)
312  real(RP), intent(out) :: TRL_URB_t (UKS:UKE,IA,JA)
313  real(RP), intent(out) :: TBL_URB_t (UKS:UKE,IA,JA)
314  real(RP), intent(out) :: TGL_URB_t (UKS:UKE,IA,JA)
315  real(RP), intent(out) :: RAINR_URB_t(IA,JA)
316  real(RP), intent(out) :: RAINB_URB_t(IA,JA)
317  real(RP), intent(out) :: RAING_URB_t(IA,JA)
318  real(RP), intent(out) :: ROFF_URB_t (IA,JA)
319 
320  real(RP), intent(out) :: SFC_TEMP(IA,JA)
321  real(RP), intent(out) :: ALBD_LW (IA,JA)
322  real(RP), intent(out) :: ALBD_SW (IA,JA)
323  real(RP), intent(out) :: MWFLX (IA,JA)
324  real(RP), intent(out) :: MUFLX (IA,JA)
325  real(RP), intent(out) :: MVFLX (IA,JA)
326  real(RP), intent(out) :: SHFLX (IA,JA)
327  real(RP), intent(out) :: LHFLX (IA,JA)
328  real(RP), intent(out) :: GHFLX (IA,JA)
329  real(RP), intent(out) :: Z0M (IA,JA)
330  real(RP), intent(out) :: Z0H (IA,JA)
331  real(RP), intent(out) :: Z0E (IA,JA)
332  real(RP), intent(out) :: U10 (IA,JA)
333  real(RP), intent(out) :: V10 (IA,JA)
334  real(RP), intent(out) :: T2 (IA,JA)
335  real(RP), intent(out) :: Q2 (IA,JA)
336 
337  real(RP), intent(in) :: TMPA (IA,JA)
338  real(RP), intent(in) :: PRSA (IA,JA)
339  real(RP), intent(in) :: W1 (IA,JA)
340  real(RP), intent(in) :: U1 (IA,JA)
341  real(RP), intent(in) :: V1 (IA,JA)
342  real(RP), intent(in) :: DENS (IA,JA)
343  real(RP), intent(in) :: QA (IA,JA)
344  real(RP), intent(in) :: Z1 (IA,JA)
345  real(RP), intent(in) :: PBL (IA,JA)
346  real(RP), intent(in) :: PRSS (IA,JA)
347  real(RP), intent(in) :: LWD (IA,JA,2)
348  real(RP), intent(in) :: SWD (IA,JA,2)
349  real(RP), intent(in) :: PREC (IA,JA)
350 
351  real(RP), intent(in) :: TR_URB (IA,JA)
352  real(RP), intent(in) :: TB_URB (IA,JA)
353  real(RP), intent(in) :: TG_URB (IA,JA)
354  real(RP), intent(in) :: TC_URB (IA,JA)
355  real(RP), intent(in) :: QC_URB (IA,JA)
356  real(RP), intent(in) :: UC_URB (IA,JA)
357  real(RP), intent(in) :: TRL_URB (UKS:UKE,IA,JA)
358  real(RP), intent(in) :: TBL_URB (UKS:UKE,IA,JA)
359  real(RP), intent(in) :: TGL_URB (UKS:UKE,IA,JA)
360  real(RP), intent(in) :: RAINR_URB(IA,JA)
361  real(RP), intent(in) :: RAINB_URB(IA,JA)
362  real(RP), intent(in) :: RAING_URB(IA,JA)
363  real(RP), intent(in) :: ROFF_URB (IA,JA)
364 
365  real(RP), intent(in) :: LON
366  real(RP), intent(in) :: LAT
367  integer, intent(in) :: NOWDATE(6)
368  real(DP), intent(in) :: dt
369 
370  ! work
371  real(RP) :: TR
372  real(RP) :: TB
373  real(RP) :: TG
374  real(RP) :: TC
375  real(RP) :: QC
376  real(RP) :: UC
377  real(RP) :: TRL(UKS:UKE)
378  real(RP) :: TBL(UKS:UKE)
379  real(RP) :: TGL(UKS:UKE)
380  real(RP) :: RAINR
381  real(RP) :: RAINB
382  real(RP) :: RAING
383  real(RP) :: ROFF
384 
385  real(RP) :: SHR (IA,JA)
386  real(RP) :: SHB (IA,JA)
387  real(RP) :: SHG (IA,JA)
388  real(RP) :: LHR (IA,JA)
389  real(RP) :: LHB (IA,JA)
390  real(RP) :: LHG (IA,JA)
391  real(RP) :: GHR (IA,JA)
392  real(RP) :: GHB (IA,JA)
393  real(RP) :: GHG (IA,JA)
394  real(RP) :: RNR (IA,JA)
395  real(RP) :: RNB (IA,JA)
396  real(RP) :: RNG (IA,JA)
397  real(RP) :: RNgrd(IA,JA)
398 
399  real(RP) :: Ustar ! friction velocity [m]
400  real(RP) :: Tstar ! friction temperature [K]
401  real(RP) :: Qstar ! friction mixing rate [kg/kg]
402  real(RP) :: Uabs ! modified absolute velocity [m/s]
403 
404  real(RP) :: QVsat ! saturation water vapor mixing ratio at surface [kg/kg]
405 
406  real(RP) :: FracU10 ! calculation parameter for U10 [-]
407  real(RP) :: FracT2 ! calculation parameter for T2 [-]
408  real(RP) :: FracQ2 ! calculation parameter for Q2 [-]
409 
410  integer :: k, i, j
411  !---------------------------------------------------------------------------
412 
413  if( io_l ) write(io_fid_log,*) '*** Urban surface physics step: Single Layer Canopy'
414 
415  do j = js, je
416  do i = is, ie
417 
418  if( is_urb(i,j) ) then
419 
420  uabs = max( sqrt( u1(i,j)**2 + v1(i,j)**2 + w1(i,j)**2 ), uabs_min )
421 
422  ! save
423  tr = tr_urb(i,j)
424  tb = tb_urb(i,j)
425  tg = tg_urb(i,j)
426  tc = tc_urb(i,j)
427  qc = qc_urb(i,j)
428  uc = uc_urb(i,j)
429 
430  do k = uks, uke
431  trl(k) = trl_urb(k,i,j)
432  tbl(k) = tbl_urb(k,i,j)
433  tgl(k) = tgl_urb(k,i,j)
434  end do
435 
436  rainr = rainr_urb(i,j)
437  rainb = rainb_urb(i,j)
438  raing = raing_urb(i,j)
439  roff = roff_urb(i,j)
440 
441  call slc_main( tr, & ! [INOUT]
442  tb, & ! [INOUT]
443  tg, & ! [INOUT]
444  tc, & ! [INOUT]
445  qc, & ! [INOUT]
446  uc, & ! [INOUT]
447  trl(:), & ! [INOUT]
448  tbl(:), & ! [INOUT]
449  tgl(:), & ! [INOUT]
450  rainr, & ! [INOUT]
451  rainb, & ! [INOUT]
452  raing, & ! [INOUT]
453  roff, & ! [INOUT]
454  albd_lw(i,j), & ! [OUT]
455  albd_sw(i,j), & ! [OUT]
456  shr(i,j), & ! [OUT]
457  shb(i,j), & ! [OUT]
458  shg(i,j), & ! [OUT]
459  lhr(i,j), & ! [OUT]
460  lhb(i,j), & ! [OUT]
461  lhg(i,j), & ! [OUT]
462  ghr(i,j), & ! [OUT]
463  ghb(i,j), & ! [OUT]
464  ghg(i,j), & ! [OUT]
465  rnr(i,j), & ! [OUT]
466  rnb(i,j), & ! [OUT]
467  rng(i,j), & ! [OUT]
468  sfc_temp(i,j), & ! [OUT]
469  rngrd(i,j), & ! [OUT]
470  shflx(i,j), & ! [OUT]
471  lhflx(i,j), & ! [OUT]
472  ghflx(i,j), & ! [OUT]
473  u10(i,j), & ! [OUT]
474  v10(i,j), & ! [OUT]
475  t2(i,j), & ! [OUT]
476  q2(i,j), & ! [OUT]
477  lsolar, & ! [IN]
478  prsa(i,j), & ! [IN]
479  prss(i,j), & ! [IN]
480  tmpa(i,j), & ! [IN]
481  qa(i,j), & ! [IN]
482  uabs, & ! [IN]
483  u1(i,j), & ! [IN]
484  v1(i,j), & ! [IN]
485  z1(i,j), & ! [IN]
486  swd(i,j,:), & ! [IN]
487  lwd(i,j,:), & ! [IN]
488  prec(i,j), & ! [IN]
489  dens(i,j), & ! [IN]
490  lon, & ! [IN]
491  lat, & ! [IN]
492  nowdate(:), & ! [IN]
493  dt, i, j ) ! [IN]
494 
495  ! calculate tendency
496  tr_urb_t(i,j) = ( tr - tr_urb(i,j) ) / dt
497  tb_urb_t(i,j) = ( tb - tb_urb(i,j) ) / dt
498  tg_urb_t(i,j) = ( tg - tg_urb(i,j) ) / dt
499  tc_urb_t(i,j) = ( tc - tc_urb(i,j) ) / dt
500  qc_urb_t(i,j) = ( qc - qc_urb(i,j) ) / dt
501  uc_urb_t(i,j) = ( uc - uc_urb(i,j) ) / dt
502 
503  do k = uks, uke
504  trl_urb_t(k,i,j) = ( trl(k) - trl_urb(k,i,j) ) / dt
505  tbl_urb_t(k,i,j) = ( tbl(k) - tbl_urb(k,i,j) ) / dt
506  tgl_urb_t(k,i,j) = ( tgl(k) - tgl_urb(k,i,j) ) / dt
507  end do
508 
509  rainr_urb_t(i,j) = ( rainr - rainr_urb(i,j) ) / dt
510  rainb_urb_t(i,j) = ( rainb - rainb_urb(i,j) ) / dt
511  raing_urb_t(i,j) = ( raing - raing_urb(i,j) ) / dt
512  roff_urb_t(i,j) = ( roff - roff_urb(i,j) ) / dt
513 
514  ! saturation at the surface
515  call qsat( qvsat, sfc_temp(i,j), prss(i,j) )
516 
517  call bulkflux( ustar, & ! [OUT]
518  tstar, & ! [OUT]
519  qstar, & ! [OUT]
520  uabs, & ! [OUT]
521  fracu10, & ! [OUT]
522  fract2, & ! [OUT]
523  fracq2, & ! [OUT]
524  tmpa(i,j), & ! [IN]
525  sfc_temp(i,j), & ! [IN]
526  prsa(i,j), & ! [IN]
527  prss(i,j), & ! [IN]
528  qa(i,j), & ! [IN]
529  qvsat, & ! [IN]
530  u1(i,j), & ! [IN]
531  v1(i,j), & ! [IN]
532  z1(i,j), & ! [IN]
533  pbl(i,j), & ! [IN]
534  z0c, & ! [IN]
535  z0hc, & ! [IN]
536  z0hc ) ! [IN]
537 
538  mwflx(i,j) = -dens(i,j) * ustar**2 / uabs * w1(i,j)
539  muflx(i,j) = -dens(i,j) * ustar**2 / uabs * u1(i,j)
540  mvflx(i,j) = -dens(i,j) * ustar**2 / uabs * v1(i,j)
541 
542  z0m(i,j) = z0c
543  z0h(i,j) = z0hc
544  z0e(i,j) = z0hc
545 
546  else
547  ! not calculate urban module
548  tr_urb_t(i,j) = 0.0_rp
549  tb_urb_t(i,j) = 0.0_rp
550  tg_urb_t(i,j) = 0.0_rp
551  tc_urb_t(i,j) = 0.0_rp
552  qc_urb_t(i,j) = 0.0_rp
553  uc_urb_t(i,j) = 0.0_rp
554  trl_urb_t(:,i,j) = 0.0_rp
555  tbl_urb_t(:,i,j) = 0.0_rp
556  tgl_urb_t(:,i,j) = 0.0_rp
557  rainr_urb_t(i,j) = 0.0_rp
558  rainb_urb_t(i,j) = 0.0_rp
559  raing_urb_t(i,j) = 0.0_rp
560  roff_urb_t(i,j) = 0.0_rp
561  sfc_temp(i,j) = 300.0_rp ! constant value
562  albd_lw(i,j) = 0.0_rp
563  albd_sw(i,j) = 0.0_rp
564  mwflx(i,j) = 0.0_rp
565  muflx(i,j) = 0.0_rp
566  mvflx(i,j) = 0.0_rp
567  shflx(i,j) = 0.0_rp
568  lhflx(i,j) = 0.0_rp
569  ghflx(i,j) = 0.0_rp
570  z0m(i,j) = 0.0_rp
571  z0h(i,j) = 0.0_rp
572  z0e(i,j) = 0.0_rp
573  u10(i,j) = 0.0_rp
574  v10(i,j) = 0.0_rp
575  t2(i,j) = 0.0_rp
576  q2(i,j) = 0.0_rp
577  !
578  shr(i,j) = 0.0_rp
579  shb(i,j) = 0.0_rp
580  shg(i,j) = 0.0_rp
581  lhr(i,j) = 0.0_rp
582  lhb(i,j) = 0.0_rp
583  lhg(i,j) = 0.0_rp
584  ghr(i,j) = 0.0_rp
585  ghb(i,j) = 0.0_rp
586  ghg(i,j) = 0.0_rp
587  rnr(i,j) = 0.0_rp
588  rnb(i,j) = 0.0_rp
589  rng(i,j) = 0.0_rp
590  rngrd(i,j) = 0.0_rp
591  endif
592 
593  end do
594  end do
595 
596  call hist_in( shr(:,:), 'URBAN_SHR', 'urban sensible heat flux on roof', 'W/m2' )
597  call hist_in( shb(:,:), 'URBAN_SHB', 'urban sensible heat flux on wall', 'W/m2' )
598  call hist_in( shg(:,:), 'URBAN_SHG', 'urban sensible heat flux on road', 'W/m2' )
599  call hist_in( lhr(:,:), 'URBAN_LHR', 'urban latent heat flux on roof', 'W/m2' )
600  call hist_in( lhb(:,:), 'URBAN_LHB', 'urban latent heat flux on wall', 'W/m2' )
601  call hist_in( lhg(:,:), 'URBAN_LHG', 'urban latent heat flux on road', 'W/m2' )
602  call hist_in( ghr(:,:), 'URBAN_GHR', 'urban ground heat flux on roof', 'W/m2' )
603  call hist_in( ghb(:,:), 'URBAN_GHB', 'urban ground heat flux on wall', 'W/m2' )
604  call hist_in( ghg(:,:), 'URBAN_GHG', 'urban ground heat flux on road', 'W/m2' )
605  call hist_in( rnr(:,:), 'URBAN_RNR', 'urban net radiation on roof', 'W/m2' )
606  call hist_in( rnb(:,:), 'URBAN_RNB', 'urban net radiation on wall', 'W/m2' )
607  call hist_in( rng(:,:), 'URBAN_RNG', 'urban net radiation on road', 'W/m2' )
608  call hist_in( rngrd(:,:), 'URBAN_RNgrd', 'urban grid average of net radiation', 'W/m2' )
609 
610  return
module ATMOSPHERE / Saturation adjustment
module grid index
module Surface bulk flux
module HISTORY
module urban grid index
Here is the call graph for this function:
Here is the caller graph for this function: