SCALE-RM
Functions/Subroutines
scale_urban_dyn_kusaka01 Module Reference

module urban / dynamics / Kusaka01 More...

Functions/Subroutines

subroutine, public urban_dyn_kusaka01_setup (UIA, UIS, UIE, UJA, UJS, UJE, Z0M, Z0H, Z0E)
 Setup. More...
 
subroutine, public urban_dyn_kusaka01 (UKA, UKS, UKE, UIA, UIS, UIE, UJA, UJS, UJE, TMPA, PRSA, W1, U1, V1, DENS, QA, LHV, Z1, PBL, RHOS, PRSS, LWD, SWD, RAIN, SNOW, CDZ, fact_urban, tloc, dsec, dt, TRL_URB, TBL_URB, TGL_URB, TR_URB, TB_URB, TG_URB, TC_URB, QC_URB, UC_URB, RAINR_URB, RAINB_URB, RAING_URB, ROFF_URB, SFC_TEMP, ALBEDO, MWFLX, MUFLX, MVFLX, SHFLX, LHFLX, GHFLX, Z0M, Z0H, Z0E, U10, V10, T2, Q2)
 Main routine for land submodel. More...
 

Detailed Description

module urban / dynamics / Kusaka01

Description
Single-layer Urban Canopy Model (Kusaka et al. 2001, BLM)
Author
Team SCALE
NAMELIST
  • PARAM_URBAN_DYN_KUSAKA01
    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 URBAN_GHB
URBAN_GHG urban ground heat flux on road W/m2 URBAN_GHG
URBAN_GHR urban ground heat flux on roof W/m2 URBAN_GHR
URBAN_LHB urban latent heat flux on wall W/m2 URBAN_LHB
URBAN_LHG urban latent heat flux on road W/m2 URBAN_LHG
URBAN_LHR urban latent heat flux on roof W/m2 URBAN_LHR
URBAN_RNB urban net radiation on wall W/m2 URBAN_RNB
URBAN_RNG urban net radiation on road W/m2 URBAN_RNG
URBAN_RNR urban net radiation on roof W/m2 URBAN_RNR
URBAN_RNgrd urban grid average of net radiation W/m2 URBAN_RNgrd
URBAN_SHB urban sensible heat flux on wall W/m2 URBAN_SHB
URBAN_SHG urban sensible heat flux on road W/m2 URBAN_SHG
URBAN_SHR urban sensible heat flux on roof W/m2 URBAN_SHR

Function/Subroutine Documentation

◆ urban_dyn_kusaka01_setup()

subroutine, public scale_urban_dyn_kusaka01::urban_dyn_kusaka01_setup ( integer, intent(in)  UIA,
integer, intent(in)  UIS,
integer, intent(in)  UIE,
integer, intent(in)  UJA,
integer, intent(in)  UJS,
integer, intent(in)  UJE,
real(rp), dimension(uia,uja), intent(out)  Z0M,
real(rp), dimension(uia,uja), intent(out)  Z0H,
real(rp), dimension(uia,uja), intent(out)  Z0E 
)

Setup.

Definition at line 118 of file scale_urban_dyn_kusaka01.F90.

References scale_file_history::file_history_reg(), scale_io::io_fid_conf, and scale_prc::prc_abort().

Referenced by mod_urban_driver::urban_driver_setup().

118  use scale_prc, only: &
119  prc_abort
120  use scale_file_history, only: &
122  implicit none
123 
124  integer, intent(in) :: uia, uis, uie
125  integer, intent(in) :: uja, ujs, uje
126  real(RP), intent(out) :: z0m(uia,uja)
127  real(RP), intent(out) :: z0h(uia,uja)
128  real(RP), intent(out) :: z0e(uia,uja)
129 
130  namelist / param_urban_dyn_kusaka01 / &
131  dts_max, &
132  zr, &
133  roof_width, &
134  road_width, &
135  sigma_zed, &
136  ah, &
137  alh, &
138  betr, &
139  betb, &
140  betg, &
141  strgr, &
142  strgb, &
143  strgg, &
144  capr, &
145  capb, &
146  capg, &
147  aksr, &
148  aksb, &
149  aksg, &
150  albr, &
151  albb, &
152  albg, &
153  epsr, &
154  epsb, &
155  epsg, &
156  z0r, &
157  z0b, &
158  z0g, &
159  trlend, &
160  tblend, &
161  tglend, &
162  bound
163 
164  integer :: i, j
165  integer :: ierr
166  !---------------------------------------------------------------------------
167 
168  log_newline
169  log_info("URBAN_DYN_kusaka01_setup",*) 'Setup'
170 
171  !--- read namelist
172  rewind(io_fid_conf)
173  read(io_fid_conf,nml=param_urban_dyn_kusaka01,iostat=ierr)
174  if( ierr < 0 ) then !--- missing
175  log_info("URBAN_DYN_kusaka01_setup",*) 'Not found namelist. Default used.'
176  elseif( ierr > 0 ) then !--- fatal error
177  log_error("URBAN_DYN_kusaka01_setup",*) 'Not appropriate names in namelist PARAM_URBAN_DYN_KUSAKA01. Check!'
178  call prc_abort
179  endif
180  log_nml(param_urban_dyn_kusaka01)
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  do j = ujs, uje
191  do i = uis, uie
192  z0m(i,j) = z0c
193  z0h(i,j) = z0hc
194  z0e(i,j) = z0hc
195  enddo
196  enddo
197 
198  call file_history_reg( 'URBAN_SHR', 'urban sensible heat flux on roof', 'W/m2', i_shr , ndims=2 )
199  call file_history_reg( 'URBAN_SHB', 'urban sensible heat flux on wall', 'W/m2', i_shb , ndims=2 )
200  call file_history_reg( 'URBAN_SHG', 'urban sensible heat flux on road', 'W/m2', i_shg , ndims=2 )
201  call file_history_reg( 'URBAN_LHR', 'urban latent heat flux on roof', 'W/m2', i_lhr , ndims=2 )
202  call file_history_reg( 'URBAN_LHB', 'urban latent heat flux on wall', 'W/m2', i_lhb , ndims=2 )
203  call file_history_reg( 'URBAN_LHG', 'urban latent heat flux on road', 'W/m2', i_lhg , ndims=2 )
204  call file_history_reg( 'URBAN_GHR', 'urban ground heat flux on roof', 'W/m2', i_ghr , ndims=2 )
205  call file_history_reg( 'URBAN_GHB', 'urban ground heat flux on wall', 'W/m2', i_ghb , ndims=2 )
206  call file_history_reg( 'URBAN_GHG', 'urban ground heat flux on road', 'W/m2', i_ghg , ndims=2 )
207  call file_history_reg( 'URBAN_RNR', 'urban net radiation on roof', 'W/m2', i_rnr , ndims=2 )
208  call file_history_reg( 'URBAN_RNB', 'urban net radiation on wall', 'W/m2', i_rnb , ndims=2 )
209  call file_history_reg( 'URBAN_RNG', 'urban net radiation on road', 'W/m2', i_rng , ndims=2 )
210  call file_history_reg( 'URBAN_RNgrd', 'urban grid average of net radiation', 'W/m2', i_rngrd, ndims=2 )
211 
212  return
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
subroutine, public file_history_reg(name, desc, unit, itemid, standard_name, ndims, dim_type, cell_measures, fill_halo)
Register/Append variable to history file.
module PROCESS
Definition: scale_prc.F90:11
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module file_history
Here is the call graph for this function:
Here is the caller graph for this function:

◆ urban_dyn_kusaka01()

subroutine, public scale_urban_dyn_kusaka01::urban_dyn_kusaka01 ( integer, intent(in)  UKA,
integer, intent(in)  UKS,
integer, intent(in)  UKE,
integer, intent(in)  UIA,
integer, intent(in)  UIS,
integer, intent(in)  UIE,
integer, intent(in)  UJA,
integer, intent(in)  UJS,
integer, intent(in)  UJE,
real(rp), dimension(uia,uja), intent(in)  TMPA,
real(rp), dimension(uia,uja), intent(in)  PRSA,
real(rp), dimension (uia,uja), intent(in)  W1,
real(rp), dimension (uia,uja), intent(in)  U1,
real(rp), dimension (uia,uja), intent(in)  V1,
real(rp), dimension(uia,uja), intent(in)  DENS,
real(rp), dimension (uia,uja), intent(in)  QA,
real(rp), dimension (uia,uja), intent(in)  LHV,
real(rp), dimension (uia,uja), intent(in)  Z1,
real(rp), dimension (uia,uja), intent(in)  PBL,
real(rp), dimension(uia,uja), intent(in)  RHOS,
real(rp), dimension(uia,uja), intent(in)  PRSS,
real(rp), dimension (uia,uja,2), intent(in)  LWD,
real(rp), dimension (uia,uja,2), intent(in)  SWD,
real(rp), dimension(uia,uja), intent(in)  RAIN,
real(rp), dimension(uia,uja), intent(in)  SNOW,
real(rp), dimension(uka), intent(in)  CDZ,
real(rp), dimension(uia,uja), intent(in)  fact_urban,
integer, intent(in)  tloc,
real(rp), intent(in)  dsec,
real(dp), intent(in)  dt,
real(rp), dimension (uks:uke,uia,uja), intent(inout)  TRL_URB,
real(rp), dimension (uks:uke,uia,uja), intent(inout)  TBL_URB,
real(rp), dimension (uks:uke,uia,uja), intent(inout)  TGL_URB,
real(rp), dimension (uia,uja), intent(inout)  TR_URB,
real(rp), dimension (uia,uja), intent(inout)  TB_URB,
real(rp), dimension (uia,uja), intent(inout)  TG_URB,
real(rp), dimension (uia,uja), intent(inout)  TC_URB,
real(rp), dimension (uia,uja), intent(inout)  QC_URB,
real(rp), dimension (uia,uja), intent(inout)  UC_URB,
real(rp), dimension(uia,uja), intent(inout)  RAINR_URB,
real(rp), dimension(uia,uja), intent(inout)  RAINB_URB,
real(rp), dimension(uia,uja), intent(inout)  RAING_URB,
real(rp), dimension (uia,uja), intent(inout)  ROFF_URB,
real(rp), dimension(uia,uja), intent(out)  SFC_TEMP,
real(rp), dimension (uia,uja,n_rad_dir,n_rad_rgn), intent(out)  ALBEDO,
real(rp), dimension (uia,uja), intent(out)  MWFLX,
real(rp), dimension (uia,uja), intent(out)  MUFLX,
real(rp), dimension (uia,uja), intent(out)  MVFLX,
real(rp), dimension (uia,uja), intent(out)  SHFLX,
real(rp), dimension (uia,uja), intent(out)  LHFLX,
real(rp), dimension (uia,uja), intent(out)  GHFLX,
real(rp), dimension (uia,uja), intent(out)  Z0M,
real(rp), dimension (uia,uja), intent(out)  Z0H,
real(rp), dimension (uia,uja), intent(out)  Z0E,
real(rp), dimension (uia,uja), intent(out)  U10,
real(rp), dimension (uia,uja), intent(out)  V10,
real(rp), dimension (uia,uja), intent(out)  T2,
real(rp), dimension (uia,uja), intent(out)  Q2 
)

Main routine for land submodel.

Definition at line 241 of file scale_urban_dyn_kusaka01.F90.

References scale_bulkflux::bulkflux, scale_const::const_cpdry, scale_const::const_eps, scale_const::const_pi, scale_const::const_rdry, scale_const::const_rvap, scale_cpl_sfc_index::i_r_diffuse, scale_cpl_sfc_index::i_r_direct, scale_cpl_sfc_index::i_r_ir, scale_cpl_sfc_index::i_r_nir, scale_cpl_sfc_index::i_r_vis, scale_prc::prc_abort(), and scale_prc::prc_myrank.

Referenced by mod_urban_driver::urban_driver_calc_tendency().

241  use scale_const, only: &
242  rdry => const_rdry, &
243  rvap => const_rvap
244  use scale_atmos_saturation, only: &
245  qsat => atmos_saturation_pres2qsat_all
246  use scale_bulkflux, only: &
247  bulkflux
248  implicit none
249  integer, intent(in) :: uka, uks, uke
250  integer, intent(in) :: uia, uis, uie
251  integer, intent(in) :: uja, ujs, uje
252 
253  real(RP), intent(in) :: tmpa(uia,uja)
254  real(RP), intent(in) :: prsa(uia,uja)
255  real(RP), intent(in) :: w1 (uia,uja)
256  real(RP), intent(in) :: u1 (uia,uja)
257  real(RP), intent(in) :: v1 (uia,uja)
258  real(RP), intent(in) :: dens(uia,uja)
259  real(RP), intent(in) :: qa (uia,uja)
260  real(RP), intent(in) :: lhv (uia,uja)
261  real(RP), intent(in) :: z1 (uia,uja)
262  real(RP), intent(in) :: pbl (uia,uja)
263  real(RP), intent(in) :: rhos(uia,uja) ! density at the surface [kg/m3]
264  real(RP), intent(in) :: prss(uia,uja)
265  real(RP), intent(in) :: lwd (uia,uja,2)
266  real(RP), intent(in) :: swd (uia,uja,2)
267  real(RP), intent(in) :: rain(uia,uja)
268  real(RP), intent(in) :: snow(uia,uja)
269  real(RP), intent(in) :: cdz(uka)
270  real(RP), intent(in) :: fact_urban(uia,uja)
271  integer, intent(in) :: tloc
272  real(RP), intent(in) :: dsec
273  real(DP), intent(in) :: dt
274 
275  real(RP), intent(inout) :: tr_urb (uia,uja)
276  real(RP), intent(inout) :: tb_urb (uia,uja)
277  real(RP), intent(inout) :: tg_urb (uia,uja)
278  real(RP), intent(inout) :: tc_urb (uia,uja)
279  real(RP), intent(inout) :: qc_urb (uia,uja)
280  real(RP), intent(inout) :: uc_urb (uia,uja)
281  real(RP), intent(inout) :: trl_urb (uks:uke,uia,uja)
282  real(RP), intent(inout) :: tbl_urb (uks:uke,uia,uja)
283  real(RP), intent(inout) :: tgl_urb (uks:uke,uia,uja)
284  real(RP), intent(inout) :: rainr_urb(uia,uja)
285  real(RP), intent(inout) :: rainb_urb(uia,uja)
286  real(RP), intent(inout) :: raing_urb(uia,uja)
287  real(RP), intent(inout) :: roff_urb (uia,uja)
288 
289  real(RP), intent(out) :: sfc_temp(uia,uja)
290  real(RP), intent(out) :: albedo (uia,uja,n_rad_dir,n_rad_rgn)
291  real(RP), intent(out) :: mwflx (uia,uja)
292  real(RP), intent(out) :: muflx (uia,uja)
293  real(RP), intent(out) :: mvflx (uia,uja)
294  real(RP), intent(out) :: shflx (uia,uja)
295  real(RP), intent(out) :: lhflx (uia,uja)
296  real(RP), intent(out) :: ghflx (uia,uja)
297  real(RP), intent(out) :: z0m (uia,uja)
298  real(RP), intent(out) :: z0h (uia,uja)
299  real(RP), intent(out) :: z0e (uia,uja)
300  real(RP), intent(out) :: u10 (uia,uja)
301  real(RP), intent(out) :: v10 (uia,uja)
302  real(RP), intent(out) :: t2 (uia,uja)
303  real(RP), intent(out) :: q2 (uia,uja)
304 
305 
306  ! parameter
307  logical, parameter :: lsolar = .false. ! [true=both, false=SSG only]
308 
309  real(RP), parameter :: uabs_min = 0.1_rp
310 
311  ! work
312  real(RP) :: tr
313  real(RP) :: tb
314  real(RP) :: tg
315  real(RP) :: tc
316  real(RP) :: qc
317  real(RP) :: uc
318  real(RP) :: trl(uks:uke)
319  real(RP) :: tbl(uks:uke)
320  real(RP) :: tgl(uks:uke)
321  real(RP) :: rainr
322  real(RP) :: rainb
323  real(RP) :: raing
324  real(RP) :: roff
325  real(RP) :: albd_lw
326  real(RP) :: albd_sw
327 
328  real(RP) :: shr (uia,uja)
329  real(RP) :: shb (uia,uja)
330  real(RP) :: shg (uia,uja)
331  real(RP) :: lhr (uia,uja)
332  real(RP) :: lhb (uia,uja)
333  real(RP) :: lhg (uia,uja)
334  real(RP) :: ghr (uia,uja)
335  real(RP) :: ghb (uia,uja)
336  real(RP) :: ghg (uia,uja)
337  real(RP) :: rnr (uia,uja)
338  real(RP) :: rnb (uia,uja)
339  real(RP) :: rng (uia,uja)
340  real(RP) :: rngrd(uia,uja)
341 
342  real(RP) :: dzr(uka) ! thickness of each roof layer [m]
343  real(RP) :: dzb(uka) ! thickness of each building layer [m]
344  real(RP) :: dzg(uka) ! thickness of each road layer [m]
345 
346 
347  real(RP) :: ustar ! friction velocity [m]
348  real(RP) :: tstar ! friction temperature [K]
349  real(RP) :: qstar ! friction mixing rate [kg/kg]
350  real(RP) :: uabs ! modified absolute velocity [m/s]
351  real(RP) :: ra ! Aerodynamic resistance (=1/Ce) [1/s]
352 
353  real(RP) :: qvsat ! saturation water vapor mixing ratio at surface [kg/kg]
354  real(RP) :: rtot ! total gas constant
355  real(RP) :: qdry ! dry air mass ratio [kg/kg]
356 
357  real(RP) :: fracu10 ! calculation parameter for U10 [-]
358  real(RP) :: fract2 ! calculation parameter for T2 [-]
359  real(RP) :: fracq2 ! calculation parameter for Q2 [-]
360 
361  integer :: k, i, j
362  !---------------------------------------------------------------------------
363 
364  log_progress(*) 'urban / dynamics / Kusaka01'
365 
366  dzr(:) = cdz(:)
367  dzb(:) = cdz(:)
368  dzg(:) = cdz(:)
369 
370  do j = ujs, uje
371  do i = uis, uie
372 
373  if( fact_urban(i,j) > 0.0_rp ) then
374 
375  qdry = 1.0_rp - qa(i,j)
376  rtot = qdry * rdry + qa(i,j) * rvap
377 
378  uabs = max( sqrt( u1(i,j)**2 + v1(i,j)**2 + w1(i,j)**2 ), uabs_min )
379 
380  ! save
381  tr = tr_urb(i,j)
382  tb = tb_urb(i,j)
383  tg = tg_urb(i,j)
384  tc = tc_urb(i,j)
385  qc = qc_urb(i,j)
386  uc = uc_urb(i,j)
387 
388  do k = uks, uke
389  trl(k) = trl_urb(k,i,j)
390  tbl(k) = tbl_urb(k,i,j)
391  tgl(k) = tgl_urb(k,i,j)
392  end do
393 
394  rainr = rainr_urb(i,j)
395  rainb = rainb_urb(i,j)
396  raing = raing_urb(i,j)
397  roff = roff_urb(i,j)
398 
399  call slc_main( uka, uks, uke, uia, uis, uie, uja, ujs, uje, &
400  trl(:), & ! [INOUT]
401  tbl(:), & ! [INOUT]
402  tgl(:), & ! [INOUT]
403  tr, & ! [INOUT]
404  tb, & ! [INOUT]
405  tg, & ! [INOUT]
406  tc, & ! [INOUT]
407  qc, & ! [INOUT]
408  uc, & ! [INOUT]
409  rainr, & ! [INOUT]
410  rainb, & ! [INOUT]
411  raing, & ! [INOUT]
412  roff, & ! [INOUT]
413  albd_lw, & ! [OUT]
414  albd_sw, & ! [OUT]
415  shr(i,j), & ! [OUT]
416  shb(i,j), & ! [OUT]
417  shg(i,j), & ! [OUT]
418  lhr(i,j), & ! [OUT]
419  lhb(i,j), & ! [OUT]
420  lhg(i,j), & ! [OUT]
421  ghr(i,j), & ! [OUT]
422  ghb(i,j), & ! [OUT]
423  ghg(i,j), & ! [OUT]
424  rnr(i,j), & ! [OUT]
425  rnb(i,j), & ! [OUT]
426  rng(i,j), & ! [OUT]
427  sfc_temp(i,j), & ! [OUT]
428  rngrd(i,j), & ! [OUT]
429  shflx(i,j), & ! [OUT]
430  lhflx(i,j), & ! [OUT]
431  ghflx(i,j), & ! [OUT]
432  u10(i,j), & ! [OUT]
433  v10(i,j), & ! [OUT]
434  t2(i,j), & ! [OUT]
435  q2(i,j), & ! [OUT]
436  lsolar, & ! [IN]
437  prsa(i,j), & ! [IN]
438  prss(i,j), & ! [IN]
439  tmpa(i,j), & ! [IN]
440  qa(i,j), & ! [IN]
441  uabs, & ! [IN]
442  u1(i,j), & ! [IN]
443  v1(i,j), & ! [IN]
444  lhv(i,j), & ! [IN]
445  z1(i,j), & ! [IN]
446  swd(i,j,:), & ! [IN]
447  lwd(i,j,:), & ! [IN]
448  rain(i,j), & ! [IN]
449  snow(i,j), & ! [IN]
450  dens(i,j), & ! [IN]
451  dzr(:), dzg(:), dzb(:), & ! [IN]
452  tloc, dsec, dt, & ! (in)
453  i, j ) ! [IN]
454 
455 
456  ! update
457  tr_urb(i,j) = tr
458  tb_urb(i,j) = tb
459  tg_urb(i,j) = tg
460  tc_urb(i,j) = tc
461  qc_urb(i,j) = qc
462  uc_urb(i,j) = uc
463  do k = uks, uke
464  trl_urb(k,i,j) = trl(k)
465  tbl_urb(k,i,j) = tbl(k)
466  tgl_urb(k,i,j) = tgl(k)
467  end do
468  rainr_urb(i,j) = rainr
469  rainb_urb(i,j) = rainb
470  raing_urb(i,j) = raing
471  roff_urb(i,j) = roff
472 
473  albedo(i,j,i_r_direct ,i_r_ir ) = albd_lw
474  albedo(i,j,i_r_diffuse,i_r_ir ) = albd_lw
475  albedo(i,j,i_r_direct ,i_r_nir) = albd_sw
476  albedo(i,j,i_r_diffuse,i_r_nir) = albd_sw
477  albedo(i,j,i_r_direct ,i_r_vis) = albd_sw
478  albedo(i,j,i_r_diffuse,i_r_vis) = albd_sw
479 
480  ! saturation at the surface
481  call qsat( sfc_temp(i,j), prss(i,j), qdry, & ! [IN]
482  qvsat ) ! [OUT]
483 
484  call bulkflux( ustar, & ! [OUT]
485  tstar, & ! [OUT]
486  qstar, & ! [OUT]
487  uabs, & ! [OUT]
488  ra, & ! [OUT]
489  fracu10, & ! [OUT]
490  fract2, & ! [OUT]
491  fracq2, & ! [OUT]
492  tmpa(i,j), & ! [IN]
493  sfc_temp(i,j), & ! [IN]
494  prsa(i,j), & ! [IN]
495  prss(i,j), & ! [IN]
496  qa(i,j), & ! [IN]
497  qvsat, & ! [IN]
498  u1(i,j), & ! [IN]
499  v1(i,j), & ! [IN]
500  z1(i,j), & ! [IN]
501  pbl(i,j), & ! [IN]
502  z0c, & ! [IN]
503  z0hc, & ! [IN]
504  z0hc ) ! [IN]
505 
506  mwflx(i,j) = -rhos(i,j) * ustar**2 / uabs * w1(i,j)
507  muflx(i,j) = -rhos(i,j) * ustar**2 / uabs * u1(i,j)
508  mvflx(i,j) = -rhos(i,j) * ustar**2 / uabs * v1(i,j)
509 
510  z0m(i,j) = z0c
511  z0h(i,j) = z0hc
512  z0e(i,j) = z0hc
513 
514  else
515  sfc_temp(i,j) = 300.0_rp ! constant value
516  albedo(i,j,:,:) = 0.0_rp
517  mwflx(i,j) = 0.0_rp
518  muflx(i,j) = 0.0_rp
519  mvflx(i,j) = 0.0_rp
520  shflx(i,j) = 0.0_rp
521  lhflx(i,j) = 0.0_rp
522  ghflx(i,j) = 0.0_rp
523  z0m(i,j) = 0.0_rp
524  z0h(i,j) = 0.0_rp
525  z0e(i,j) = 0.0_rp
526  u10(i,j) = 0.0_rp
527  v10(i,j) = 0.0_rp
528  t2(i,j) = 0.0_rp
529  q2(i,j) = 0.0_rp
530  shr(i,j) = 0.0_rp
531  shb(i,j) = 0.0_rp
532  shg(i,j) = 0.0_rp
533  lhr(i,j) = 0.0_rp
534  lhb(i,j) = 0.0_rp
535  lhg(i,j) = 0.0_rp
536  ghr(i,j) = 0.0_rp
537  ghb(i,j) = 0.0_rp
538  ghg(i,j) = 0.0_rp
539  rnr(i,j) = 0.0_rp
540  rnb(i,j) = 0.0_rp
541  rng(i,j) = 0.0_rp
542  rngrd(i,j) = 0.0_rp
543  endif
544 
545  end do
546  end do
547 
548  call put_history( uia, uja, &
549  shr(:,:), shb(:,:), shg(:,:), &
550  lhr(:,:), lhb(:,:), lhg(:,:), &
551  ghr(:,:), ghb(:,:), ghg(:,:), &
552  rnr(:,:), rnb(:,:), rng(:,:), &
553  rngrd(:,:) )
554 
555  return
module atmosphere / saturation
integer, parameter, public i_r_vis
integer, public qa
integer, parameter, public n_rad_dir
integer, parameter, public n_rad_rgn
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:55
procedure(bc), pointer, public bulkflux
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:63
module CONSTANT
Definition: scale_const.F90:11
integer, parameter, public i_r_direct
module Surface bulk flux
integer, parameter, public i_r_nir
integer, parameter, public i_r_ir
integer, parameter, public i_r_diffuse
Here is the call graph for this function:
Here is the caller graph for this function: