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, fact_urban, Z0M, Z0H, Z0E, ZD, AH_URB, AHL_URB, AH_TOFFSET)
 Setup. More...
 
subroutine, public urban_dyn_kusaka01_finalize
 Finalize. More...
 
subroutine, public urban_dyn_kusaka01 (UKA, UKS, UKE, UIA, UIS, UIE, UJA, UJS, UJE, TMPA, PRSA, U1, V1, DENS, QA, LHV, Z1, RHOS, PRSS, LWD, SWD, RAIN, EFLX, Z0M, Z0H, Z0E, ZD, CDZ, TanSL_X, TanSL_Y, fact_urban, 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, Ustar, Tstar, Qstar, Wstar, RLmo, U10, V10, T2, Q2)
 Main routine for land submodel. More...
 
subroutine put_history (UIA, UJA, SHR, SHB, SHG, LHR, LHB, LHG, GHR, GHB, GHG, RNR, RNB, RNG, RNgrd)
 

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 step [K/step]
    BOUND integer
    DEBUG logical .false.
    URBAN_DYN_KUSAKA01_PARAM_IN_FILENAME character(len=H_LONG) '' urban parameter table
    URBAN_DYN_KUSAKA01_GRIDDED_Z0M_IN_FILENAME character(len=H_LONG) '' gridded data of Z0M
    URBAN_DYN_KUSAKA01_GRIDDED_Z0H_IN_FILENAME character(len=H_LONG) '' gridded data of Z0H

  • PARAM_URBAN_DATA
    nametypedefault valuecomment
    ZR real(RP) 10.0_RP roof level (building height) [m]
    ROOF_WIDTH real(RP) 9.0_RP roof width [m]
    ROAD_WIDTH real(RP) 11.0_RP road width [m]
    SIGMA_ZED real(RP) 1.0_RP Standard deviation of roof height [m]
    AH_TBL real(RP) 17.5_RP Sensible Anthropogenic heat from urban subgrid [W/m^2]
    AHL_TBL real(RP) 0.0_RP Latent Anthropogenic heat from urban subgrid [W/m^2]
    BETR real(RP)
    BETB real(RP)
    BETG real(RP)
    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]

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(in)  fact_urban,
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)  ZD,
real(rp), dimension (uia,uja,1:24), intent(out)  AH_URB,
real(rp), dimension (uia,uja,1:24), intent(out)  AHL_URB,
real(rp), intent(out)  AH_TOFFSET 
)

Setup.

Definition at line 139 of file scale_urban_dyn_kusaka01.F90.

139  use scale_prc, only: &
140  prc_myrank, &
141  prc_abort
142  use scale_const, only: &
143  undef => const_undef
144  use scale_file_history, only: &
146  implicit none
147 
148  integer, intent(in) :: UIA, UIS, UIE
149  integer, intent(in) :: UJA, UJS, UJE
150  real(RP), intent(in) :: fact_urban(UIA,UJA)
151  real(RP), intent(out) :: Z0M(UIA,UJA)
152  real(RP), intent(out) :: Z0H(UIA,UJA)
153  real(RP), intent(out) :: Z0E(UIA,UJA)
154  real(RP), intent(out) :: ZD (UIA,UJA)
155  real(RP), intent(out) :: AH_URB (UIA,UJA,1:24)
156  real(RP), intent(out) :: AHL_URB (UIA,UJA,1:24)
157  real(RP), intent(out) :: AH_TOFFSET
158 
159  character(len=H_LONG) :: URBAN_DYN_KUSAKA01_PARAM_IN_FILENAME = ''
160  character(len=H_LONG) :: URBAN_DYN_KUSAKA01_GRIDDED_Z0M_IN_FILENAME = ''
161  character(len=H_LONG) :: URBAN_DYN_KUSAKA01_GRIDDED_Z0M_IN_VARNAME = 'URBAN_Z0M'
162  character(len=H_LONG) :: URBAN_DYN_KUSAKA01_GRIDDED_Z0H_IN_FILENAME = ''
163  character(len=H_LONG) :: URBAN_DYN_KUSAKA01_GRIDDED_Z0H_IN_VARNAME = 'URBAN_Z0H'
164  !character(len=H_LONG) :: URBAN_DYN_KUSAKA01_GRIDDED_ZD_IN_FILENAME = '' !< gridded data of ZD
165  !character(len=H_LONG) :: URBAN_DYN_KUSAKA01_GRIDDED_ZD_IN_VARNAME = 'URBAN_ZD' !< var name of gridded data for Zd
166  character(len=H_LONG) :: URBAN_DYN_KUSAKA01_GRIDDED_AH_IN_FILENAME = ''
167  character(len=H_LONG) :: URBAN_DYN_KUSAKA01_GRIDDED_AH_IN_VARNAME = 'URBAN_AH'
168  character(len=H_LONG) :: URBAN_DYN_KUSAKA01_GRIDDED_AHL_IN_FILENAME = ''
169  character(len=H_LONG) :: URBAN_DYN_KUSAKA01_GRIDDED_AHL_IN_VARNAME = 'URBAN_AHL'
170 
171  namelist / param_urban_dyn_kusaka01 / &
172  dts_max, &
173  bound, &
174  debug, &
175  urban_dyn_kusaka01_param_in_filename, &
176  urban_dyn_kusaka01_gridded_z0m_in_filename, &
177  urban_dyn_kusaka01_gridded_z0h_in_filename, &
178  !URBAN_DYN_KUSAKA01_GRIDDED_ZD_IN_FILENAME, &
179  urban_dyn_kusaka01_gridded_ah_in_filename, &
180  urban_dyn_kusaka01_gridded_ahl_in_filename
181 
182  real(RP) :: udata(UIA,UJA)
183  real(RP) :: udata2(UIA,UJA,24)
184  real(RP) :: rtime
185  integer :: itime
186 
187  real(RP) :: ahdiurnal(1:24) ! AH diurnal profile
188 
189  integer :: i, j, k
190  integer :: ierr
191  !---------------------------------------------------------------------------
192 
193  log_newline
194  log_info("URBAN_DYN_kusaka01_setup",*) 'Setup'
195 
196  !$acc data copyin(fact_urban) &
197  !$acc copyout(Z0M, Z0H, Z0E, ZD, AH_URB, AHL_URB)
198 
199  !--- read namelist
200  rewind(io_fid_conf)
201  read(io_fid_conf,nml=param_urban_dyn_kusaka01,iostat=ierr)
202  if( ierr < 0 ) then !--- missing
203  log_info("URBAN_DYN_kusaka01_setup",*) 'Not found namelist. Default used.'
204  elseif( ierr > 0 ) then !--- fatal error
205  log_error("URBAN_DYN_kusaka01_setup",*) 'Not appropriate names in namelist PARAM_URBAN_DYN_KUSAKA01. Check!'
206  call prc_abort
207  endif
208  log_nml(param_urban_dyn_kusaka01)
209 
210  !-- read urban parameter from file
211  if( urban_dyn_kusaka01_param_in_filename /= '' ) then
212  call read_urban_param_table( trim(urban_dyn_kusaka01_param_in_filename) )
213  endif
214 
215  ! set other urban parameters
216  call urban_param_setup
217 
218  ! Local time: 1 to 24
219  ahdiurnal(:) = (/ 0.356, 0.274, 0.232, 0.251, 0.375, 0.647, 0.919, 1.135, 1.249, 1.328, &
220  1.365, 1.363, 1.375, 1.404, 1.457, 1.526, 1.557, 1.521, 1.372, 1.206, &
221  1.017, 0.876, 0.684, 0.512 /)
222 
223  ! Shift to UTC based on local solar timezone of domain center
224  rtime = modulo(base_lon, 360.0_rp) / 15.0_rp
225  itime = nint(rtime)
226  ah_toffset = rtime - itime ! currently not used: difference of time from domain center
227  ahdiurnal(:) = cshift(ahdiurnal(:), -1*itime) ! convert from LT to UTC
228 
229  do j = ujs, uje
230  do i = uis, uie
231  z0m(i,j) = z0c_tbl
232  z0h(i,j) = z0hc_tbl
233  z0e(i,j) = z0hc_tbl
234  zd(i,j) = zdc_tbl
235  do k = 1, 24 ! UTC
236  ah_urb(i,j,k) = ah_tbl * ahdiurnal(k) * fact_urban(i,j)
237  ahl_urb(i,j,k) = ahl_tbl * ahdiurnal(k) * fact_urban(i,j)
238  enddo
239  enddo
240  enddo
241 
242  !-- replace gridded Z0M data if there is a file
243  if( urban_dyn_kusaka01_gridded_z0m_in_filename /= '' ) then
244  udata(:,:) = 0.0_rp
245  call read_urban_gridded_data_2d( &
246  uia, uja, &
247  urban_dyn_kusaka01_gridded_z0m_in_filename, &
248  urban_dyn_kusaka01_gridded_z0m_in_varname, &
249  udata(:,:) )
250 
251  ! replace to gridded data
252  do j = ujs, uje
253  do i = uis, uie
254  if( udata(i,j) /= undef )then
255  if ( udata(i,j) > 0.0_rp ) then
256  z0m(i,j) = udata(i,j)
257  else if ( udata(i,j) < 0.0_rp ) then
258  log_error("URBAN_DYN_kusaka01_setup",*) 'Gridded Z0M data includes data less than 0. Please check data!',i,j
259  call prc_abort
260  else ! Z0M = 0[m]
261  log_warn("URBAN_DYN_kusaka01_setup",*) 'Gridded Z0M data includes 0; default or table value is used to avoid zero division',prc_myrank,i,j
262  endif
263  endif
264  enddo
265  enddo
266  endif
267 
268  !-- read gridded Z0H & Z0E data from a file
269  if( urban_dyn_kusaka01_gridded_z0h_in_filename /= '' ) then
270  udata = 0.0_rp
271  call read_urban_gridded_data_2d( &
272  uia, uja, &
273  urban_dyn_kusaka01_gridded_z0h_in_filename, &
274  urban_dyn_kusaka01_gridded_z0h_in_varname, &
275  udata )
276 
277  ! replace to gridded data
278  do j = ujs, uje
279  do i = uis, uie
280  if( udata(i,j) /= undef )then
281  if ( udata(i,j) > 0.0_rp ) then
282  z0h(i,j) = udata(i,j)
283  z0e(i,j) = udata(i,j)
284  else if ( udata(i,j) < 0.0_rp ) then
285  log_error("URBAN_DYN_kusaka01_setup",*) 'Gridded Z0H data includes data less than 0. Please check data!',i,j
286  call prc_abort
287  else ! Z0H = 0[m]
288  log_warn("URBAN_DYN_kusaka01_setup",*) 'Gridded Z0H data includes 0; default or table value is used to avoid zero division',prc_myrank,i,j
289  endif
290  endif
291  enddo
292  enddo
293  endif
294 
295  ! currently NOT USED
296  !-- read gridded ZD data from a file
297  !if( URBAN_DYN_KUSAKA01_GRIDDED_ZD_IN_FILENAME /= '' ) then
298  ! udata = 0.0_RP
299  ! call read_urban_gridded_data_2D( &
300  ! UIA, UJA, &
301  ! URBAN_DYN_KUSAKA01_GRIDDED_ZD_IN_FILENAME, &
302  ! URBAN_DYN_KUSAKA01_GRIDDED_ZD_IN_VARNAME, &
303  ! udata )
304  !
305  ! ! replace to gridded data
306  ! do j = UJS, UJE
307  ! do i = UIS, UIE
308  ! if( udata(i,j) /= UNDEF )then
309  ! if ( udata(i,j) >= 0.0_RP ) then
310  ! ZD(i,j) = udata(i,j)
311  ! else
312  ! LOG_ERROR("URBAN_DYN_kusaka01_setup",*) 'Gridded ZD data includes data less than 0. Please check data!',i,j
313  ! call PRC_abort
314  ! endif
315  ! endif
316  ! enddo
317  ! enddo
318  !endif
319 
320  !-- read gridded AH data from a file
321  if( urban_dyn_kusaka01_gridded_ah_in_filename /= '' ) then
322  udata2 = 0.0_rp
323  call read_urban_gridded_data_3d( &
324  uia, uja, &
325  urban_dyn_kusaka01_gridded_ah_in_filename, &
326  urban_dyn_kusaka01_gridded_ah_in_varname, &
327  udata2 )
328 
329  ! replace to gridded data
330  do k = 1, 24
331  do j = ujs, uje
332  do i = uis, uie
333  if( udata2(i,j,k) /= undef )then
334  ah_urb(i,j,k) = udata2(i,j,k)
335  endif
336  enddo
337  enddo
338  enddo
339  endif
340 
341  !-- read gridded AHL data from a file
342  if( urban_dyn_kusaka01_gridded_ahl_in_filename /= '' ) then
343  udata2 = 0.0_rp
344  call read_urban_gridded_data_3d( &
345  uia, uja, &
346  urban_dyn_kusaka01_gridded_ahl_in_filename, &
347  urban_dyn_kusaka01_gridded_ahl_in_varname, &
348  udata2 )
349 
350  ! replace to gridded data
351  do k = 1, 24
352  do j = ujs, uje
353  do i = uis, uie
354  if( udata2(i,j,k) /= undef )then
355  ahl_urb(i,j,k) = udata2(i,j,k)
356  endif
357  enddo
358  enddo
359  enddo
360  endif
361 
362  call file_history_reg( 'URBAN_SHR', 'urban sensible heat flux on roof', 'W/m2', i_shr , ndims=2 )
363  call file_history_reg( 'URBAN_SHB', 'urban sensible heat flux on wall', 'W/m2', i_shb , ndims=2 )
364  call file_history_reg( 'URBAN_SHG', 'urban sensible heat flux on road', 'W/m2', i_shg , ndims=2 )
365  call file_history_reg( 'URBAN_LHR', 'urban latent heat flux on roof', 'W/m2', i_lhr , ndims=2 )
366  call file_history_reg( 'URBAN_LHB', 'urban latent heat flux on wall', 'W/m2', i_lhb , ndims=2 )
367  call file_history_reg( 'URBAN_LHG', 'urban latent heat flux on road', 'W/m2', i_lhg , ndims=2 )
368  call file_history_reg( 'URBAN_GHR', 'urban ground heat flux on roof', 'W/m2', i_ghr , ndims=2 )
369  call file_history_reg( 'URBAN_GHB', 'urban ground heat flux on wall', 'W/m2', i_ghb , ndims=2 )
370  call file_history_reg( 'URBAN_GHG', 'urban ground heat flux on road', 'W/m2', i_ghg , ndims=2 )
371  call file_history_reg( 'URBAN_RNR', 'urban net radiation on roof', 'W/m2', i_rnr , ndims=2 )
372  call file_history_reg( 'URBAN_RNB', 'urban net radiation on wall', 'W/m2', i_rnb , ndims=2 )
373  call file_history_reg( 'URBAN_RNG', 'urban net radiation on road', 'W/m2', i_rng , ndims=2 )
374  call file_history_reg( 'URBAN_RNgrd', 'urban grid average of net radiation', 'W/m2', i_rngrd, ndims=2 )
375 
376  !$acc update device(Z0M, Z0H, Z0E, ZD, AH_URB, AHL_URB)
377  !$acc end data
378 
379  !$acc update device(DTS_MAX,BOUND,debug)
380  !$acc update device(ZR,BETR_CONST,BETB_CONST,BETG_CONST,STRGR,STRGB,STRGG,CAPR,CAPB,CAPG,AKSR,AKSB,AKSG,ALBR,ALBB,ALBG,EPSR,EPSB,EPSG,Z0R,TRLEND,TBLEND,TGLEND)
381  !$acc update device(R,RW,HGT,Z0HR,Z0HB,Z0HG,SVF)
382 
383 #ifdef DEBUG_KUSAKA01
384  cnt_num1 = 0
385  cnt_num2 = 0
386  cnt_itr1 = 0
387  cnt_itr2 = 0
388  max_itr1 = 0
389  max_itr2 = 0
390  !$acc update device(cnt_num1,cnt_num2,cnt_itr1,cnt_itr2,max_itr1,max_itr2)
391 #endif
392 
393  return

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

Referenced by mod_urban_driver::urban_driver_setup().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ urban_dyn_kusaka01_finalize()

subroutine, public scale_urban_dyn_kusaka01::urban_dyn_kusaka01_finalize

Finalize.

Definition at line 399 of file scale_urban_dyn_kusaka01.F90.

399  use mpi
400  use scale_prc, only: &
403  implicit none
404 
405 #ifdef DEBUG_KUSAKA01
406  integer :: iwork1(4), iwork2(2), ierr
407  !$acc update host(cnt_num1,cnt_num2,cnt_itr1,cnt_itr2,max_itr1,max_itr2)
408  iwork1(1) = cnt_num1
409  iwork1(2) = cnt_num2
410  iwork1(3) = cnt_itr1
411  iwork1(4) = cnt_itr2
412  call mpi_reduce( mpi_in_place, iwork1, 4, mpi_integer, mpi_sum, 0, prc_local_comm_world, ierr)
413  iwork2(1) = max_itr1
414  iwork2(2) = max_itr2
415  call mpi_reduce( mpi_in_place, iwork2, 2, mpi_integer, mpi_max, 0, prc_local_comm_world, ierr)
416 
417  if ( prc_ismaster ) then
418  log_info("URBAN_DYN_kusaka01_finalize",*) "Averaged iteration count"
419  log_info_cont(*) "TR: ", real(iwork1(3),dp) / iwork1(1), ", (max ", iwork2(1), ")"
420  log_info_cont(*) "TB/TG: ", real(iwork1(4),dp) / iwork1(2), ", (max ", iwork2(2), ")"
421  end if
422 #endif
423 

References scale_precision::dp, scale_prc::prc_ismaster, and scale_prc::prc_local_comm_world.

Referenced by mod_urban_driver::urban_driver_finalize().

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)  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)  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)  EFLX,
real(rp), dimension (uia,uja), intent(in)  Z0M,
real(rp), dimension (uia,uja), intent(in)  Z0H,
real(rp), dimension (uia,uja), intent(in)  Z0E,
real(rp), dimension (uia,uja), intent(in)  ZD,
real(rp), dimension(uka), intent(in)  CDZ,
real(rp), dimension(uia,uja), intent(in)  TanSL_X,
real(rp), dimension(uia,uja), intent(in)  TanSL_Y,
real(rp), dimension(uia,uja), intent(in)  fact_urban,
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(out)  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)  Ustar,
real(rp), dimension (uia,uja), intent(out)  Tstar,
real(rp), dimension (uia,uja), intent(out)  Qstar,
real(rp), dimension (uia,uja), intent(out)  Wstar,
real(rp), dimension (uia,uja), intent(out)  RLmo,
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 456 of file scale_urban_dyn_kusaka01.F90.

456  use scale_prc, only: &
457  prc_myrank, &
458  prc_abort
459  use scale_const, only: &
460  eps => const_eps, &
461  undef => const_undef, &
462  rdry => const_rdry, &
463  rvap => const_rvap
464  use scale_bulkflux, only: &
465  bulkflux
466  implicit none
467  integer, intent(in) :: UKA, UKS, UKE
468  integer, intent(in) :: UIA, UIS, UIE
469  integer, intent(in) :: UJA, UJS, UJE
470 
471  real(RP), intent(in) :: TMPA(UIA,UJA)
472  real(RP), intent(in) :: PRSA(UIA,UJA)
473  real(RP), intent(in) :: U1 (UIA,UJA)
474  real(RP), intent(in) :: V1 (UIA,UJA)
475  real(RP), intent(in) :: DENS(UIA,UJA)
476  real(RP), intent(in) :: QA (UIA,UJA)
477  real(RP), intent(in) :: LHV (UIA,UJA)
478  real(RP), intent(in) :: Z1 (UIA,UJA)
479  real(RP), intent(in) :: RHOS(UIA,UJA) ! density at the surface [kg/m3]
480  real(RP), intent(in) :: PRSS(UIA,UJA)
481  real(RP), intent(in) :: LWD (UIA,UJA,2)
482  real(RP), intent(in) :: SWD (UIA,UJA,2)
483  real(RP), intent(in) :: RAIN(UIA,UJA)
484  real(RP), intent(in) :: EFLX(UIA,UJA)
485  real(RP), intent(in) :: Z0M (UIA,UJA)
486  real(RP), intent(in) :: Z0H (UIA,UJA)
487  real(RP), intent(in) :: Z0E (UIA,UJA)
488  real(RP), intent(in) :: ZD (UIA,UJA)
489  real(RP), intent(in) :: CDZ(UKA)
490  real(RP), intent(in) :: TanSL_X(UIA,UJA)
491  real(RP), intent(in) :: TanSL_Y(UIA,UJA)
492  real(RP), intent(in) :: fact_urban(UIA,UJA)
493  real(DP), intent(in) :: dt
494 
495  real(RP), intent(inout) :: TR_URB (UIA,UJA)
496  real(RP), intent(inout) :: TB_URB (UIA,UJA)
497  real(RP), intent(inout) :: TG_URB (UIA,UJA)
498  real(RP), intent(inout) :: TC_URB (UIA,UJA)
499  real(RP), intent(inout) :: QC_URB (UIA,UJA)
500  real(RP), intent(inout) :: UC_URB (UIA,UJA)
501  real(RP), intent(inout) :: TRL_URB (UKS:UKE,UIA,UJA)
502  real(RP), intent(inout) :: TBL_URB (UKS:UKE,UIA,UJA)
503  real(RP), intent(inout) :: TGL_URB (UKS:UKE,UIA,UJA)
504  real(RP), intent(inout) :: RAINR_URB(UIA,UJA)
505  real(RP), intent(inout) :: RAINB_URB(UIA,UJA)
506  real(RP), intent(inout) :: RAING_URB(UIA,UJA)
507  real(RP), intent(out) :: ROFF_URB (UIA,UJA)
508 
509  real(RP), intent(out) :: SFC_TEMP(UIA,UJA)
510  real(RP), intent(out) :: ALBEDO (UIA,UJA,N_RAD_DIR,N_RAD_RGN)
511  real(RP), intent(out) :: MWFLX (UIA,UJA)
512  real(RP), intent(out) :: MUFLX (UIA,UJA)
513  real(RP), intent(out) :: MVFLX (UIA,UJA)
514  real(RP), intent(out) :: SHFLX (UIA,UJA)
515  real(RP), intent(out) :: LHFLX (UIA,UJA)
516  real(RP), intent(out) :: GHFLX (UIA,UJA)
517  real(RP), intent(out) :: Ustar (UIA,UJA)
518  real(RP), intent(out) :: Tstar (UIA,UJA)
519  real(RP), intent(out) :: Qstar (UIA,UJA)
520  real(RP), intent(out) :: Wstar (UIA,UJA)
521  real(RP), intent(out) :: RLmo (UIA,UJA)
522  real(RP), intent(out) :: U10 (UIA,UJA)
523  real(RP), intent(out) :: V10 (UIA,UJA)
524  real(RP), intent(out) :: T2 (UIA,UJA)
525  real(RP), intent(out) :: Q2 (UIA,UJA)
526 
527  ! parameter
528  logical, parameter :: LSOLAR = .false. ! [true=both, false=SSG only]
529 
530  real(RP), parameter :: Uabs_min = 0.1_rp
531 
532  ! work
533  real(RP) :: TR
534  real(RP) :: TB
535  real(RP) :: TG
536  real(RP) :: TC
537  real(RP) :: QC
538  real(RP) :: UC
539  real(RP) :: TRL(UKS:UKE)
540  real(RP) :: TBL(UKS:UKE)
541  real(RP) :: TGL(UKS:UKE)
542  real(RP) :: RAINR
543  real(RP) :: RAINB
544  real(RP) :: RAING
545  real(RP) :: ALBD_LW
546  real(RP) :: ALBD_SW
547 
548  real(RP) :: SHR (UIA,UJA)
549  real(RP) :: SHB (UIA,UJA)
550  real(RP) :: SHG (UIA,UJA)
551  real(RP) :: LHR (UIA,UJA)
552  real(RP) :: LHB (UIA,UJA)
553  real(RP) :: LHG (UIA,UJA)
554  real(RP) :: GHR (UIA,UJA)
555  real(RP) :: GHB (UIA,UJA)
556  real(RP) :: GHG (UIA,UJA)
557  real(RP) :: RNR (UIA,UJA)
558  real(RP) :: RNB (UIA,UJA)
559  real(RP) :: RNG (UIA,UJA)
560  real(RP) :: RNgrd(UIA,UJA)
561 
562  real(RP) :: DZR(UKA) ! thickness of each roof layer [m]
563  real(RP) :: DZB(UKA) ! thickness of each building layer [m]
564  real(RP) :: DZG(UKA) ! thickness of each road layer [m]
565 
566  real(RP) :: SWDt(2)
567  real(RP) :: LWDt(2)
568 
569  real(RP) :: Uabs ! modified absolute velocity [m/s]
570 
571  real(RP) :: MFLX
572  real(RP) :: w
573 
574  ! work
575  real(RP) :: TRLP(UKA)
576  real(RP) :: TBLP(UKA)
577  real(RP) :: TGLP(UKA)
578  real(RP) :: A(UKA)
579  real(RP) :: B(UKA)
580  real(RP) :: C(UKA)
581  real(RP) :: D(UKA)
582  real(RP) :: P(UKA)
583  real(RP) :: Q(UKA)
584 
585  logical :: converged
586 
587 #ifdef DEBUG_KUSAKA01
588  integer :: cn1, cn2, ci1, ci2, mi1, mi2
589 #endif
590 
591  integer :: k, i, j
592  !---------------------------------------------------------------------------
593 
594  log_progress(*) 'urban / dynamics / Kusaka01'
595 
596  !$acc data copyin(TMPA,PRSA,U1,V1,DENS,QA,LHV,Z1,RHOS,PRSS,LWD,SWD,RAIN,EFLX,Z0M,Z0H,Z0E,ZD,CDZ,TanSL_X,TanSL_Y,fact_urban) &
597  !$acc copy(TR_URB,TB_URB,TG_URB,TC_URB,QC_URB,UC_URB,TRL_URB,TBL_URB,TGL_URB,RAINR_URB,RAINB_URB,RAING_URB) &
598  !$acc copyout(ROFF_URB,SFC_TEMP,ALBEDO,MWFLX,MUFLX,MVFLX,SHFLX,LHFLX,GHFLX,Ustar,Tstar,Qstar,Wstar,RLmo,U10,V10,T2,Q2) &
599  !$acc create(SHR,SHB,SHG,LHR,LHB,LHG,GHR,GHB,GHG,RNR,RNB,RNG,RNgrd,DZR,DZB,DZG)
600 
601 
602  !$acc kernels
603  dzr(:) = cdz(:)
604  dzb(:) = cdz(:)
605  dzg(:) = cdz(:)
606  !$acc end kernels
607 
608  converged = .true.
609 
610 #ifdef DEBUG_KUSAKA01
611  cn1 = 0; cn2 = 0; ci1 = 0; ci2 = 0; mi1 = 0; mi2 = 0
612 #endif
613 
614  !$omp parallel do schedule(dynamic) collapse(2) &
615 #ifdef DEBUG_KUSAKA01
616  !$omp reduction(+:cn1,cn2,ci1,ci2) reduction(max:mi1,mi2) &
617 #endif
618  !$omp private(w,Uabs,TR,TB,TG,TC,QC,UC,TRL,TBL,TGL,RAINR,RAINB,RAING,ALBD_LW,ALBD_SW,MFLX,SWDt,LWDt, &
619  !$omp TRLP,TBLP,TGLP,A,B,C,D,P,Q)
620  !$acc parallel
621  !$acc loop gang collapse(2) reduction(.and.: converged) independent &
622 #ifdef DEBUG_KUSAKA01
623  !$acc reduction(+:cn1,cn2,ci1,ci2) reduction(max:mi1,mi2) &
624 #endif
625  !$acc private(w,Uabs,TR,TB,TG,TC,QC,UC,TRL,TBL,TGL,RAINR,RAINB,RAING,ALBD_LW,ALBD_SW,MFLX,SWDt,LWDt, &
626  !$acc TRLP,TBLP,TGLP,A,B,C,D,P,Q)
627  do j = ujs, uje
628  do i = uis, uie
629 
630  if( fact_urban(i,j) > 0.0_rp ) then
631 
632  w = u1(i,j) * tansl_x(i,j) + v1(i,j) * tansl_y(i,j)
633  uabs = sqrt( u1(i,j)**2 + v1(i,j)**2 + w**2 )
634 
635  ! save
636  tr = tr_urb(i,j)
637  tb = tb_urb(i,j)
638  tg = tg_urb(i,j)
639  tc = tc_urb(i,j)
640  qc = qc_urb(i,j)
641  uc = uc_urb(i,j)
642 
643  do k = uks, uke
644  trl(k) = trl_urb(k,i,j)
645  tbl(k) = tbl_urb(k,i,j)
646  tgl(k) = tgl_urb(k,i,j)
647  end do
648 
649  rainr = rainr_urb(i,j)
650  rainb = rainb_urb(i,j)
651  raing = raing_urb(i,j)
652 
653  swdt(:) = swd(i,j,:)
654  lwdt(:) = lwd(i,j,:)
655 
656  call slc_main( uka, uks, uke, &
657  trl(:), & ! [INOUT]
658  tbl(:), & ! [INOUT]
659  tgl(:), & ! [INOUT]
660  tr, & ! [INOUT]
661  tb, & ! [INOUT]
662  tg, & ! [INOUT]
663  tc, & ! [INOUT]
664  qc, & ! [INOUT]
665  uc, & ! [INOUT]
666  rainr, & ! [INOUT]
667  rainb, & ! [INOUT]
668  raing, & ! [INOUT]
669  roff_urb(i,j), & ! [OUT]
670  albd_lw, & ! [OUT]
671  albd_sw, & ! [OUT]
672  shr(i,j), & ! [OUT]
673  shb(i,j), & ! [OUT]
674  shg(i,j), & ! [OUT]
675  lhr(i,j), & ! [OUT]
676  lhb(i,j), & ! [OUT]
677  lhg(i,j), & ! [OUT]
678  ghr(i,j), & ! [OUT]
679  ghb(i,j), & ! [OUT]
680  ghg(i,j), & ! [OUT]
681  rnr(i,j), & ! [OUT]
682  rnb(i,j), & ! [OUT]
683  rng(i,j), & ! [OUT]
684  sfc_temp(i,j), & ! [OUT]
685  rngrd(i,j), & ! [OUT]
686  shflx(i,j), & ! [OUT]
687  lhflx(i,j), & ! [OUT]
688  ghflx(i,j), & ! [OUT]
689  mflx, & ! [OUT]
690  ustar(i,j), & ! [OUT]
691  tstar(i,j), & ! [OUT]
692  qstar(i,j), & ! [OUT]
693  u10(i,j), & ! [OUT]
694  v10(i,j), & ! [OUT]
695  t2(i,j), & ! [OUT]
696  q2(i,j), & ! [OUT]
697  converged, & ! [OUT]
698  lsolar, & ! [IN]
699  prsa(i,j), & ! [IN]
700  prss(i,j), & ! [IN]
701  tmpa(i,j), & ! [IN]
702  rhos(i,j), & ! [IN]
703  qa(i,j), & ! [IN]
704  uabs, & ! [IN]
705  u1(i,j), & ! [IN]
706  v1(i,j), & ! [IN]
707  lhv(i,j), & ! [IN]
708  z1(i,j), & ! [IN]
709  swdt(:), & ! [IN]
710  lwdt(:), & ! [IN]
711  rain(i,j), & ! [IN]
712  eflx(i,j), & ! [IN]
713  dens(i,j), & ! [IN]
714  z0m(i,j), & ! [IN]
715  z0h(i,j), & ! [IN]
716  zd(i,j), & ! [IN]
717  dzr(:), & ! [IN]
718  dzg(:), dzb(:), & ! [IN]
719  dt, & ! [IN]
720  i, j, & ! [IN]
721 #ifdef DEBUG_KUSAKA01
722  cn1, cn2, ci1, ci2, mi1, mi2, & ! [INOUT]
723 #endif
724  trlp(:), tblp(:), tglp(:), & ! (work)
725  a(:), b(:), c(:), d(:), & ! (work)
726  p(:), q(:) ) ! (work)
727 
728  ! update
729  tr_urb(i,j) = tr
730  tb_urb(i,j) = tb
731  tg_urb(i,j) = tg
732  tc_urb(i,j) = tc
733  qc_urb(i,j) = qc
734  uc_urb(i,j) = uc
735  do k = uks, uke
736  trl_urb(k,i,j) = trl(k)
737  tbl_urb(k,i,j) = tbl(k)
738  tgl_urb(k,i,j) = tgl(k)
739  end do
740  rainr_urb(i,j) = rainr
741  rainb_urb(i,j) = rainb
742  raing_urb(i,j) = raing
743 
744  albedo(i,j,i_r_direct ,i_r_ir ) = albd_lw
745  albedo(i,j,i_r_diffuse,i_r_ir ) = albd_lw
746  albedo(i,j,i_r_direct ,i_r_nir) = albd_sw
747  albedo(i,j,i_r_diffuse,i_r_nir) = albd_sw
748  albedo(i,j,i_r_direct ,i_r_vis) = albd_sw
749  albedo(i,j,i_r_diffuse,i_r_vis) = albd_sw
750 
751  if ( uabs < eps ) then
752  mwflx(i,j) = 0.0_rp
753  muflx(i,j) = 0.0_rp
754  mvflx(i,j) = 0.0_rp
755  else
756  mwflx(i,j) = mflx * w / uabs
757  muflx(i,j) = mflx * u1(i,j) / uabs
758  mvflx(i,j) = mflx * v1(i,j) / uabs
759  end if
760 
761  else
762  sfc_temp(i,j) = undef
763  albedo(i,j,:,:) = undef
764  mwflx(i,j) = undef
765  muflx(i,j) = undef
766  mvflx(i,j) = undef
767  shflx(i,j) = undef
768  lhflx(i,j) = undef
769  ghflx(i,j) = undef
770  ustar(i,j) = undef
771  tstar(i,j) = undef
772  qstar(i,j) = undef
773  wstar(i,j) = undef
774  rlmo(i,j) = undef
775  u10(i,j) = undef
776  v10(i,j) = undef
777  t2(i,j) = undef
778  q2(i,j) = undef
779  shr(i,j) = undef
780  shb(i,j) = undef
781  shg(i,j) = undef
782  lhr(i,j) = undef
783  lhb(i,j) = undef
784  lhg(i,j) = undef
785  ghr(i,j) = undef
786  ghb(i,j) = undef
787  ghg(i,j) = undef
788  rnr(i,j) = undef
789  rnb(i,j) = undef
790  rng(i,j) = undef
791  rngrd(i,j) = undef
792  endif
793 
794  end do
795  end do
796  !$acc end parallel
797 
798  if ( .not. converged ) then
799  log_error("URBAN_DYN_kusaka01_SLC_main",*) "not converged"
800  call prc_abort
801  end if
802 
803  call put_history( uia, uja, &
804  shr(:,:), shb(:,:), shg(:,:), &
805  lhr(:,:), lhb(:,:), lhg(:,:), &
806  ghr(:,:), ghb(:,:), ghg(:,:), &
807  rnr(:,:), rnb(:,:), rng(:,:), &
808  rngrd(:,:) )
809 
810  !$acc end data
811 
812 #ifdef DEBUG_KUSAKA01
813  cnt_num1 = cnt_num1 + cn1
814  cnt_num2 = cnt_num2 + cn2
815  cnt_itr1 = cnt_itr1 + ci1
816  cnt_itr2 = cnt_itr2 + ci2
817  max_itr1 = max(max_itr1, mi1)
818  max_itr2 = max(max_itr2, mi2)
819 #endif
820 
821  return

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_const::const_undef, scale_file_cartesc::file_cartesc_close(), scale_file_cartesc::file_cartesc_flush(), scale_file_cartesc::file_cartesc_open(), 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_io::io_get_available_fid(), scale_io::io_get_fname(), scale_prc::prc_abort(), scale_prc::prc_myrank, and put_history().

Referenced by mod_urban_driver::urban_driver_calc_tendency().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ put_history()

subroutine scale_urban_dyn_kusaka01::put_history ( integer, intent(in)  UIA,
integer, intent(in)  UJA,
real(rp), dimension(uia,uja), intent(in)  SHR,
real(rp), dimension(uia,uja), intent(in)  SHB,
real(rp), dimension(uia,uja), intent(in)  SHG,
real(rp), dimension(uia,uja), intent(in)  LHR,
real(rp), dimension(uia,uja), intent(in)  LHB,
real(rp), dimension(uia,uja), intent(in)  LHG,
real(rp), dimension(uia,uja), intent(in)  GHR,
real(rp), dimension(uia,uja), intent(in)  GHB,
real(rp), dimension(uia,uja), intent(in)  GHG,
real(rp), dimension(uia,uja), intent(in)  RNR,
real(rp), dimension(uia,uja), intent(in)  RNB,
real(rp), dimension(uia,uja), intent(in)  RNG,
real(rp), dimension(uia,uja), intent(in)  RNgrd 
)

Definition at line 2551 of file scale_urban_dyn_kusaka01.F90.

2551  use scale_file_history, only: &
2552  file_history_put
2553  integer, intent(in) :: UIA, UJA
2554  real(RP), intent(in) :: SHR(UIA,UJA), SHB(UIA,UJA), SHG(UIA,UJA)
2555  real(RP), intent(in) :: LHR(UIA,UJA), LHB(UIA,UJA), LHG(UIA,UJA)
2556  real(RP), intent(in) :: GHR(UIA,UJA), GHB(UIA,UJA), GHG(UIA,UJA)
2557  real(RP), intent(in) :: RNR(UIA,UJA), RNB(UIA,UJA), RNG(UIA,UJA)
2558  real(RP), intent(in) :: RNgrd(UIA,UJA)
2559 
2560  call file_history_put( i_shr, shr(:,:) )
2561  call file_history_put( i_shb, shb(:,:) )
2562  call file_history_put( i_shg, shg(:,:) )
2563  call file_history_put( i_lhr, lhr(:,:) )
2564  call file_history_put( i_lhb, lhb(:,:) )
2565  call file_history_put( i_lhg, lhg(:,:) )
2566  call file_history_put( i_ghr, ghr(:,:) )
2567  call file_history_put( i_ghb, ghb(:,:) )
2568  call file_history_put( i_ghg, ghg(:,:) )
2569  call file_history_put( i_rnr, rnr(:,:) )
2570  call file_history_put( i_rnb, rnb(:,:) )
2571  call file_history_put( i_rng, rng(:,:) )
2572  call file_history_put( i_rngrd, rngrd(:,:) )
2573 
2574  return

Referenced by urban_dyn_kusaka01().

Here is the caller graph for this function:
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_cpl_sfc_index::i_r_direct
integer, parameter, public i_r_direct
Definition: scale_cpl_sfc_index.F90:37
scale_cpl_sfc_index::i_r_diffuse
integer, parameter, public i_r_diffuse
Definition: scale_cpl_sfc_index.F90:38
scale_prc::prc_local_comm_world
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:89
scale_const::const_rvap
real(rp), parameter, public const_rvap
specific gas constant (water vapor) [J/kg/K]
Definition: scale_const.F90:68
scale_cpl_sfc_index::i_r_ir
integer, parameter, public i_r_ir
Definition: scale_cpl_sfc_index.F90:29
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:35
scale_bulkflux
module Surface bulk flux
Definition: scale_bulkflux.F90:12
scale_prc::prc_myrank
integer, public prc_myrank
process num in local communicator
Definition: scale_prc.F90:91
scale_file_history
module file_history
Definition: scale_file_history.F90:15
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_cpl_sfc_index::i_r_nir
integer, parameter, public i_r_nir
Definition: scale_cpl_sfc_index.F90:30
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_precision::dp
integer, parameter, public dp
Definition: scale_precision.F90:32
scale_cpl_sfc_index::i_r_vis
integer, parameter, public i_r_vis
Definition: scale_cpl_sfc_index.F90:31
scale_const::const_rdry
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:59
scale_file_history::file_history_reg
subroutine, public file_history_reg(name, desc, unit, itemid, standard_name, ndims, dim_type, cell_measures, fill_halo)
Register/Append variable to history file.
Definition: scale_file_history.F90:685
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_prc::prc_ismaster
logical, public prc_ismaster
master process in local communicator?
Definition: scale_prc.F90:92
scale_bulkflux::bulkflux
procedure(bc), pointer, public bulkflux
Definition: scale_bulkflux.F90:84