SCALE-RM
Functions/Subroutines
mod_urban_driver Module Reference

module URBAN driver More...

Functions/Subroutines

subroutine, public urban_driver_setup
 Setup. More...
 
subroutine, public urban_driver_calc_tendency (force)
 Calclate tendency. More...
 
subroutine, public urban_driver_update
 Urban step. More...
 
subroutine, public urban_surface_get
 Get surface boundary. More...
 
subroutine, public urban_surface_set (countup)
 Set surface boundary to other model. More...
 

Detailed Description

module URBAN driver

Description
Urban module driver
Author
Team SCALE
NAMELIST
  • No namelist group
History Output
namedescriptionunitvariable
URBAN_QC_t tendency of URBAN_QC kg/kg URBAN_QC_t
URBAN_RAINB_t tendency of URBAN_RAINB K URBAN_RAINB_t
URBAN_RAING_t tendency of URBAN_RAING K URBAN_RAING_t
URBAN_RAINR_t tendency of URBAN_RAINR K URBAN_RAINR_t
URBAN_ROFF_t tendency of URBAN_ROFF K URBAN_ROFF_t
URBAN_TBL_t tendency of URBAN_TBL K URBAN_TBL_t
URBAN_TB_t tendency of URBAN_TB K URBAN_TB_t
URBAN_TC_t tendency of URBAN_TC K URBAN_TC_t
URBAN_TGL_t tendency of URBAN_TGL K URBAN_TGL_t
URBAN_TG_t tendency of URBAN_TG K URBAN_TG_t
URBAN_TRL_t tendency of URBAN_TRL K URBAN_TRL_t
URBAN_TR_t tendency of URBAN_TR K URBAN_TR_t
URBAN_UC_t tendency of URBAN_UC m/s URBAN_UC_t

Function/Subroutine Documentation

◆ urban_driver_setup()

subroutine, public mod_urban_driver::urban_driver_setup ( )

Setup.

Definition at line 52 of file mod_urban_driver.F90.

References scale_prc::prc_abort(), scale_urban_grid_cartesc_index::uia, scale_urban_grid_cartesc_index::uie, scale_urban_grid_cartesc_index::uis, scale_urban_grid_cartesc_index::uja, scale_urban_grid_cartesc_index::uje, scale_urban_grid_cartesc_index::ujs, mod_urban_admin::urban_do, scale_urban_dyn_kusaka01::urban_dyn_kusaka01_setup(), mod_urban_admin::urban_dyn_type, mod_urban_admin::urban_sfc_type, mod_urban_vars::urban_z0e, mod_urban_vars::urban_z0h, and mod_urban_vars::urban_z0m.

Referenced by mod_rm_driver::rm_driver().

52  use scale_prc, only: &
53  prc_abort
54  use mod_urban_admin, only: &
55  urban_do, &
58  use scale_urban_dyn_kusaka01, only: &
60  use mod_urban_vars, only: &
61  urban_z0m, &
62  urban_z0h, &
63  urban_z0e
64  implicit none
65  !---------------------------------------------------------------------------
66 
67  log_newline
68  log_info("URBAN_driver_setup",*) 'Setup'
69 
70  if ( urban_do ) then
71 
72  select case ( urban_dyn_type )
73  case ( 'KUSAKA01' )
74  call urban_dyn_kusaka01_setup( uia, uis, uie, uja, ujs, uje, &
75  urban_z0m(:,:), urban_z0h(:,:), urban_z0e(:,:) ) ! [OUT]
76 
77  urban_sfc_type = 'KUSAKA01'
78  case default
79  log_error("URBAN_driver_setup",*) 'LAND_DYN_TYPE is invalid: ', trim(urban_dyn_type)
80  call prc_abort
81  end select
82 
83  select case ( urban_sfc_type )
84  case ( 'KUSAKA01' )
85  ! do nothing
86  case default
87  log_error("URBAN_driver_setup",*) 'LAND_SFC_TYPE is invalid: ', trim(urban_sfc_type)
88  call prc_abort
89  end select
90 
91  end if
92 
93  return
real(rp), dimension(:,:), allocatable, public urban_z0e
real(rp), dimension(:,:), allocatable, public urban_z0m
module URBAN Variables
module PROCESS
Definition: scale_prc.F90:11
real(rp), dimension(:,:), allocatable, public urban_z0h
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
module urban / dynamics / Kusaka01
character(len=h_short), public urban_dyn_type
subroutine, public urban_dyn_kusaka01_setup(UIA, UIS, UIE, UJA, UJS, UJE, Z0M, Z0H, Z0E)
Setup.
character(len=h_short), public urban_sfc_type
module Urban admin
logical, public urban_do
Here is the call graph for this function:
Here is the caller graph for this function:

◆ urban_driver_calc_tendency()

subroutine, public mod_urban_driver::urban_driver_calc_tendency ( logical, intent(in)  force)

Calclate tendency.

Definition at line 99 of file mod_urban_driver.F90.

References mod_urban_vars::atmos_cossza, mod_urban_vars::atmos_dens, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_z1, scale_atmos_hydrometeor::atmos_hydrometeor_dry, mod_urban_vars::atmos_pbl, mod_atmos_phy_ch_driver::atmos_phy_ch_driver_urban_flux(), mod_urban_vars::atmos_pres, mod_urban_vars::atmos_qv, mod_urban_vars::atmos_sfc_dens, mod_urban_vars::atmos_sfc_pres, mod_urban_vars::atmos_sflx_lw, mod_urban_vars::atmos_sflx_rain, mod_urban_vars::atmos_sflx_snow, mod_urban_vars::atmos_sflx_sw, mod_atmos_admin::atmos_sw_phy_ch, mod_urban_vars::atmos_temp, mod_urban_vars::atmos_u, mod_urban_vars::atmos_v, mod_urban_vars::atmos_w, scale_const::const_d2r, scale_atmos_hydrometeor::i_qv, scale_landuse::landuse_fact_urban, scale_mapprojection::mapprojection_basepoint_lat, scale_mapprojection::mapprojection_basepoint_lon, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), scale_tracer::qa, scale_statistics::statistics_checktotal, scale_time::time_dtsec_urban, scale_time::time_nowdate, scale_urban_grid_cartesc_index::uia, scale_urban_grid_cartesc_index::uie, scale_urban_grid_cartesc_index::uis, scale_urban_grid_cartesc_index::uja, scale_urban_grid_cartesc_index::uje, scale_urban_grid_cartesc_index::ujs, scale_urban_grid_cartesc_index::uka, scale_urban_grid_cartesc_index::uke, scale_urban_grid_cartesc_index::uks, scale_urban_dyn_kusaka01::urban_dyn_kusaka01(), scale_urban_grid_cartesc::urban_grid_cartesc_cdz, scale_urban_grid_cartesc_real::urban_grid_cartesc_real_area, scale_urban_grid_cartesc_real::urban_grid_cartesc_real_totarea, scale_urban_grid_cartesc_real::urban_grid_cartesc_real_totvol, scale_urban_grid_cartesc_real::urban_grid_cartesc_real_vol, mod_urban_vars::urban_q2, mod_urban_vars::urban_qc, mod_urban_vars::urban_qc_t, mod_urban_vars::urban_rainb, mod_urban_vars::urban_rainb_t, mod_urban_vars::urban_raing, mod_urban_vars::urban_raing_t, mod_urban_vars::urban_rainr, mod_urban_vars::urban_rainr_t, mod_urban_vars::urban_roff, mod_urban_vars::urban_roff_t, mod_urban_vars::urban_sfc_albedo, mod_urban_vars::urban_sfc_temp, mod_urban_admin::urban_sfc_type, mod_urban_vars::urban_sflx_gh, mod_urban_vars::urban_sflx_lh, mod_urban_vars::urban_sflx_mu, mod_urban_vars::urban_sflx_mv, mod_urban_vars::urban_sflx_mw, mod_urban_vars::urban_sflx_qtrc, mod_urban_vars::urban_sflx_sh, urban_surface_get(), urban_surface_set(), mod_urban_vars::urban_t2, mod_urban_vars::urban_tb, mod_urban_vars::urban_tb_t, mod_urban_vars::urban_tbl, mod_urban_vars::urban_tbl_t, mod_urban_vars::urban_tc, mod_urban_vars::urban_tc_t, mod_urban_vars::urban_tg, mod_urban_vars::urban_tg_t, mod_urban_vars::urban_tgl, mod_urban_vars::urban_tgl_t, mod_urban_vars::urban_tr, mod_urban_vars::urban_tr_t, mod_urban_vars::urban_trl, mod_urban_vars::urban_trl_t, mod_urban_vars::urban_u10, mod_urban_vars::urban_uc, mod_urban_vars::urban_uc_t, mod_urban_vars::urban_v10, mod_urban_vars::urban_z0e, mod_urban_vars::urban_z0h, and mod_urban_vars::urban_z0m.

Referenced by mod_rm_driver::restart_read(), and mod_rm_driver::rm_driver().

99  use scale_const, only: &
100  d2r => const_d2r
101  use scale_statistics, only: &
103  statistics_total
104  use scale_urban_grid_cartesc_real, only: &
109  use scale_file_history, only: &
110  file_history_in
111  use mod_atmos_admin, only: &
113  use mod_atmos_phy_ch_driver, only: &
115  use mod_urban_vars, only: &
116  atmos_temp, &
117  atmos_pres, &
118  atmos_w, &
119  atmos_u, &
120  atmos_v, &
121  atmos_dens, &
122  atmos_qv, &
123  atmos_pbl, &
124  atmos_sfc_dens, &
125  atmos_sfc_pres, &
126  atmos_sflx_lw, &
127  atmos_sflx_sw, &
128  atmos_cossza, &
129  atmos_sflx_rain, &
130  atmos_sflx_snow, &
131  urban_trl_t, &
132  urban_tbl_t, &
133  urban_tgl_t, &
134  urban_tr_t, &
135  urban_tb_t, &
136  urban_tg_t, &
137  urban_tc_t, &
138  urban_qc_t, &
139  urban_uc_t, &
140  urban_rainr_t, &
141  urban_rainb_t, &
142  urban_raing_t, &
143  urban_roff_t, &
144  urban_sfc_temp, &
146  urban_sflx_mw, &
147  urban_sflx_mu, &
148  urban_sflx_mv, &
149  urban_sflx_sh, &
150  urban_sflx_lh, &
151  urban_sflx_qtrc, &
152  urban_sflx_gh, &
153  urban_z0m, &
154  urban_z0h, &
155  urban_z0e, &
156  urban_u10, &
157  urban_v10, &
158  urban_t2, &
159  urban_q2, &
160  urban_tr, &
161  urban_tb, &
162  urban_tg, &
163  urban_tc, &
164  urban_qc, &
165  urban_uc, &
166  urban_trl, &
167  urban_tbl, &
168  urban_tgl, &
169  urban_rainr, &
170  urban_rainb, &
171  urban_raing, &
172  urban_roff
173  use scale_atmos_hydrometeor, only: &
174  hydrometeor_lhv => atmos_hydrometeor_lhv, &
176  i_qv
177  use scale_time, only: &
178  dt => time_dtsec_urban, &
179  nowdate => time_nowdate
180  use scale_mapprojection, only: &
181  base_lon => mapprojection_basepoint_lon, &
182  base_lat => mapprojection_basepoint_lat
183  use scale_atmos_grid_cartesc_real, only: &
184  real_z1 => atmos_grid_cartesc_real_z1
185  use scale_urban_grid_cartesc, only: &
187  use scale_landuse, only: &
189  use mod_urban_admin, only: &
191  use scale_urban_dyn_kusaka01, only: &
193  implicit none
194  logical, intent(in) :: force
195 
196  real(RP) :: trl(uka,uia,uja), tbl(uka,uia,uja), tgl(uka,uia,uja)
197  real(RP) :: tr(uia,uja), tb(uia,uja), tg(uia,uja)
198  real(RP) :: tc(uia,uja), qc(uia,uja), uc(uia,uja)
199  real(RP) :: rainr(uia,uja), rainb(uia,uja), raing(uia,uja), roff(uia,uja)
200 
201  real(RP) :: lhv(uia,uja) ! latent heat of vaporization [J/kg]
202 
203  real(RP) :: lat, lon ! [deg]
204  integer :: tloc ! local time (1-24h)
205  real(RP) :: dsec ! second [s]
206 
207  integer :: k, i, j, iq
208  !---------------------------------------------------------------------------
209 
210  call prof_rapstart('URB_CalcTend', 1)
211 
212  !########## Get Surface Boundary from coupler ##########
213  call urban_surface_get
214 
215  !########## initialize tendency ##########
216 !OCL XFILL
217  do j = ujs, uje
218  do i = uis, uie
219  do k = uks, uke
220  urban_trl_t(k,i,j) = 0.0_rp
221  urban_tbl_t(k,i,j) = 0.0_rp
222  urban_tgl_t(k,i,j) = 0.0_rp
223  end do
224  end do
225  end do
226 
227 !OCL XFILL
228  do j = ujs, uje
229  do i = uis, uie
230  urban_tr_t(i,j) = 0.0_rp
231  urban_tb_t(i,j) = 0.0_rp
232  urban_tg_t(i,j) = 0.0_rp
233  urban_tc_t(i,j) = 0.0_rp
234  urban_qc_t(i,j) = 0.0_rp
235  urban_uc_t(i,j) = 0.0_rp
236 
237  urban_rainr_t(i,j) = 0.0_rp
238  urban_rainb_t(i,j) = 0.0_rp
239  urban_raing_t(i,j) = 0.0_rp
240  urban_roff_t(i,j) = 0.0_rp
241  enddo
242  enddo
243 
244 !OCL XFILL
245  !$omp parallel do
246  do iq = 1, qa
247  do j = ujs, uje
248  do i = uis, uie
249  urban_sflx_qtrc(i,j,iq) = 0.0_rp
250  enddo
251  enddo
252  enddo
253 
254  select case ( urban_sfc_type )
255  case ( 'KUSAKA01' )
256 
257 !OCL XFILL
258  !$omp parallel do
259  do j = ujs, uje
260  do i = uis, uie
261  do k = uks, uke
262  trl(k,i,j) = urban_trl(k,i,j)
263  tbl(k,i,j) = urban_tbl(k,i,j)
264  tgl(k,i,j) = urban_tgl(k,i,j)
265  end do
266  end do
267  end do
268 
269 !OCL XFILL
270  !$omp parallel do
271  do j = ujs, uje
272  do i = uis, uie
273  tr(i,j) = urban_tr(i,j)
274  tb(i,j) = urban_tb(i,j)
275  tg(i,j) = urban_tg(i,j)
276  tc(i,j) = urban_tc(i,j)
277  qc(i,j) = urban_qc(i,j)
278  uc(i,j) = urban_uc(i,j)
279  rainr(i,j) = urban_rainr(i,j)
280  rainb(i,j) = urban_rainb(i,j)
281  raing(i,j) = urban_raing(i,j)
282  roff(i,j) = urban_roff(i,j)
283  end do
284  end do
285 
286  ! local time
287  lat = base_lat
288  lon = base_lon
289  if (lon < 0.0_rp ) lon = mod(lon, 360.0_rp) + 360.0_rp
290  if (lon > 360.0_rp ) lon = mod(lon, 360.0_rp)
291  tloc = mod( (nowdate(4) + int(lon/15.0_rp)),24 )
292  dsec = real( NOWDATE(5)*60.0_RP + NOWDATE(6), kind=RP ) / 3600.0_rp
293  if( tloc == 0 ) tloc = 24
294 
295  call hydrometeor_lhv( uia, uis, uie, uja, ujs, uje, &
296  atmos_temp(:,:), lhv(:,:) )
297 
298  call urban_dyn_kusaka01( uka, uks, uke, uia, uis, uie, uja, ujs, uje, &
299  atmos_temp(:,:), atmos_pres(:,:), & ! [IN]
300  atmos_w(:,:), atmos_u(:,:), atmos_v(:,:), & ! [IN]
301  atmos_dens(:,:), atmos_qv(:,:), lhv(:,:), & ! [IN]
302  real_z1(:,:), atmos_pbl(:,:), & ! [IN]
303  atmos_sfc_dens(:,:), atmos_sfc_pres(:,:), & ! [IN]
304  atmos_sflx_lw(:,:,:), atmos_sflx_sw(:,:,:), & ! [IN]
305  atmos_sflx_rain(:,:), atmos_sflx_snow(:,:), & ! [IN]
306  cdz(:), & ! [IN]
307  landuse_fact_urban(:,:), & ! [IN]
308  tloc, dsec, dt, & ! [IN]
309  trl(:,:,:), tbl(:,:,:), tgl(:,:,:), & ! [INOUT]
310  tr(:,:), tb(:,:), tg(:,:), tc(:,:), qc(:,:), uc(:,:), & ! [INOUT]
311  rainr(:,:), rainb(:,:), raing(:,:), roff(:,:), & ! [INOUT]
312  urban_sfc_temp(:,:), & ! [OUT]
313  urban_sfc_albedo(:,:,:,:), & ! [OUT]
314  urban_sflx_mw(:,:), urban_sflx_mu(:,:), urban_sflx_mv(:,:), & ! [OUT]
315  urban_sflx_sh(:,:), urban_sflx_lh(:,:), urban_sflx_gh(:,:), & ! [OUT]
316  urban_z0m(:,:), urban_z0h(:,:), urban_z0e(:,:), & ! [OUT]
317  urban_u10(:,:), urban_v10(:,:), urban_t2(:,:), urban_q2(:,:) ) ! [OUT]
318 
319 !OCL XFILL
320  !$omp parallel do
321  do j = ujs, uje
322  do i = uis, uie
323  do k = uks, uke
324  urban_trl_t(k,i,j) = ( trl(k,i,j) - urban_trl(k,i,j) ) / dt
325  urban_tbl_t(k,i,j) = ( tbl(k,i,j) - urban_tbl(k,i,j) ) / dt
326  urban_tgl_t(k,i,j) = ( tgl(k,i,j) - urban_tgl(k,i,j) ) / dt
327  end do
328  end do
329  end do
330 
331 !OCL XFILL
332  !$omp parallel do
333  do j = ujs, uje
334  do i = uis, uie
335  urban_tr_t(i,j) = ( tr(i,j) - urban_tr(i,j) ) / dt
336  urban_tb_t(i,j) = ( tb(i,j) - urban_tb(i,j) ) / dt
337  urban_tg_t(i,j) = ( tg(i,j) - urban_tg(i,j) ) / dt
338  urban_tc_t(i,j) = ( tc(i,j) - urban_tc(i,j) ) / dt
339  urban_qc_t(i,j) = ( qc(i,j) - urban_qc(i,j) ) / dt
340  urban_uc_t(i,j) = ( uc(i,j) - urban_uc(i,j) ) / dt
341  urban_rainr_t(i,j) = ( rainr(i,j) - urban_rainr(i,j) ) / dt
342  urban_rainb_t(i,j) = ( rainb(i,j) - urban_rainb(i,j) ) / dt
343  urban_raing_t(i,j) = ( raing(i,j) - urban_raing(i,j) ) / dt
344  urban_roff_t(i,j) = ( roff(i,j) - urban_roff(i,j) ) / dt
345  end do
346  end do
347 
348  if ( .NOT. atmos_hydrometeor_dry ) then
349  !$omp parallel do
350  do j = ujs, uje
351  do i = uis, uie
352  urban_sflx_qtrc(i,j,i_qv) = urban_sflx_lh(i,j) / lhv(i,j)
353  enddo
354  enddo
355  endif
356 
357  end select
358 
359  ! Surface flux for chemical tracers
360  if ( atmos_sw_phy_ch ) then
361  call atmos_phy_ch_driver_urban_flux( urban_sflx_qtrc(:,:,:) ) ! [INOUT]
362  endif
363 
364  !########## Set Surface Boundary to coupler ##########
365  call urban_surface_set( countup=.true. )
366 
367 
368  call file_history_in( urban_tr_t(:,:), 'URBAN_TR_t', 'tendency of URBAN_TR', 'K', dim_type='XY' )
369  call file_history_in( urban_tb_t(:,:), 'URBAN_TB_t', 'tendency of URBAN_TB', 'K', dim_type='XY' )
370  call file_history_in( urban_tg_t(:,:), 'URBAN_TG_t', 'tendency of URBAN_TG', 'K', dim_type='XY' )
371  call file_history_in( urban_tc_t(:,:), 'URBAN_TC_t', 'tendency of URBAN_TC', 'K', dim_type='XY' )
372  call file_history_in( urban_qc_t(:,:), 'URBAN_QC_t', 'tendency of URBAN_QC', 'kg/kg', dim_type='XY' )
373  call file_history_in( urban_uc_t(:,:), 'URBAN_UC_t', 'tendency of URBAN_UC', 'm/s', dim_type='XY' )
374 
375  call file_history_in( urban_trl_t(:,:,:), 'URBAN_TRL_t', 'tendency of URBAN_TRL', 'K', dim_type='UXY' )
376  call file_history_in( urban_tbl_t(:,:,:), 'URBAN_TBL_t', 'tendency of URBAN_TBL', 'K', dim_type='UXY' )
377  call file_history_in( urban_tgl_t(:,:,:), 'URBAN_TGL_t', 'tendency of URBAN_TGL', 'K', dim_type='UXY' )
378 
379  call file_history_in( urban_rainr_t(:,:), 'URBAN_RAINR_t', 'tendency of URBAN_RAINR', 'K', dim_type='XY' )
380  call file_history_in( urban_rainb_t(:,:), 'URBAN_RAINB_t', 'tendency of URBAN_RAINB', 'K', dim_type='XY' )
381  call file_history_in( urban_raing_t(:,:), 'URBAN_RAING_t', 'tendency of URBAN_RAING', 'K', dim_type='XY' )
382  call file_history_in( urban_roff_t(:,:), 'URBAN_ROFF_t', 'tendency of URBAN_ROFF', 'K', dim_type='XY' )
383 
384 
385  if ( statistics_checktotal ) then
386 
387  call statistics_total( uka, uks, uke, uia, uis, uie, uja, ujs, uje, &
388  urban_trl_t(:,:,:), 'URBAN_TRL_t', &
391  call statistics_total( uka, uks, uke, uia, uis, uie, uja, ujs, uje, &
392  urban_tbl_t(:,:,:), 'URBAN_TBL_t', &
395  call statistics_total( uka, uks, uke, uia, uis, uie, uja, ujs, uje, &
396  urban_tgl_t(:,:,:), 'URBAN_TGL_t', &
399 
400  call statistics_total( uia, uis, uie, uja, ujs, uje, &
401  urban_tr_t(:,:), 'URBAN_TR_t', &
404  call statistics_total( uia, uis, uie, uja, ujs, uje, &
405  urban_tb_t(:,:), 'URBAN_TB_t', &
408  call statistics_total( uia, uis, uie, uja, ujs, uje, &
409  urban_tg_t(:,:), 'URBAN_TG_t', &
412  call statistics_total( uia, uis, uie, uja, ujs, uje, &
413  urban_tc_t(:,:), 'URBAN_TC_t', &
416  call statistics_total( uia, uis, uie, uja, ujs, uje, &
417  urban_qc_t(:,:), 'URBAN_QC_t', &
420  call statistics_total( uia, uis, uie, uja, ujs, uje, &
421  urban_uc_t(:,:), 'URBAN_UC_t', &
424 
425  call statistics_total( uia, uis, uie, uja, ujs, uje, &
426  urban_rainr_t(:,:), 'URBAN_RAINR_t', &
429  call statistics_total( uia, uis, uie, uja, ujs, uje, &
430  urban_rainb_t(:,:), 'URBAN_RAINB_t', &
433  call statistics_total( uia, uis, uie, uja, ujs, uje, &
434  urban_raing_t(:,:), 'URBAN_RAING_t', &
437  call statistics_total( uia, uis, uie, uja, ujs, uje, &
438  urban_roff_t(:,:), 'URBAN_ROFF_t', &
441  endif
442 
443 
444  call prof_rapend ('URB_CalcTend', 1)
445 
446  return
module ATMOS admin
real(rp), dimension(:,:), allocatable, public urban_qc_t
real(rp), dimension(:,:,:), allocatable, public urban_sflx_qtrc
subroutine, public atmos_phy_ch_driver_urban_flux(SFLX_QTRC)
Driver.
real(rp), dimension(:,:), allocatable, public urban_qc
real(rp), dimension(:,:), allocatable, public urban_rainr_t
real(rp), dimension(:,:), allocatable, public landuse_fact_urban
urban factor
real(rp), dimension(:,:,:), allocatable, public urban_trl_t
real(rp), dimension(:,:), allocatable, public urban_tg_t
real(rp), dimension(:,:), allocatable, public urban_u10
real(rp), dimension(:,:), allocatable, public urban_v10
real(rp), public urban_grid_cartesc_real_totarea
total area
real(rp), dimension(:,:), allocatable, public urban_raing_t
real(rp), dimension(:,:), allocatable, public urban_z0e
real(rp), dimension(:,:), allocatable, public urban_sflx_mu
real(rp), dimension(:,:), allocatable, public urban_z0m
real(rp), dimension(:,:), allocatable, public urban_tb_t
module ATMOSPHERE / Physics Chemistry
real(rp), dimension(:,:), allocatable, public urban_t2
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:32
real(rp), dimension(:,:), allocatable, public urban_tb
module URBAN Variables
real(rp), dimension(:,:), allocatable, public urban_raing
logical, public statistics_checktotal
calc&report variable totals to logfile?
real(rp), dimension(:,:), allocatable, public urban_uc
real(rp), dimension(:,:,:,:), allocatable, public urban_sfc_albedo
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_z1
Height of the lowermost grid from surface (cell center) [m].
subroutine, public urban_surface_set(countup)
Set surface boundary to other model.
real(rp), dimension(:,:,:), allocatable, public urban_grid_cartesc_real_vol
volume of grid cell
real(rp), dimension(:,:), allocatable, public atmos_pbl
real(rp), dimension(:,:), allocatable, public urban_sflx_sh
real(rp), dimension(:,:), allocatable, public urban_tr
real(rp), dimension(:,:,:), allocatable, public urban_tgl
real(rp), public urban_grid_cartesc_real_totvol
total volume
module atmosphere / hydrometeor
module LANDUSE
real(rp), dimension(:,:), allocatable, public atmos_qv
real(rp), dimension(:,:), allocatable, public atmos_sfc_pres
real(rp), dimension(:), allocatable, public urban_grid_cartesc_cdz
z-length of control volume [m]
real(rp), dimension(:,:), allocatable, public atmos_cossza
real(rp), dimension(:,:), allocatable, public urban_uc_t
real(rp), dimension(:,:), allocatable, public urban_z0h
real(rp), dimension(:,:,:), allocatable, public atmos_sflx_sw
real(rp), dimension(:,:), allocatable, public atmos_v
module TIME
Definition: scale_time.F90:16
real(rp), dimension(:,:), allocatable, public urban_roff
module urban / grid / cartesianC
real(rp), dimension(:,:), allocatable, public atmos_dens
logical, public atmos_sw_phy_ch
real(rp), dimension(:,:), allocatable, public urban_grid_cartesc_real_area
area of grid cell
real(rp), dimension(:,:), allocatable, public atmos_sflx_rain
module CONSTANT
Definition: scale_const.F90:11
real(rp), dimension(:,:), allocatable, public urban_sflx_lh
real(rp), dimension(:,:), allocatable, public urban_tc
real(rp), dimension(:,:), allocatable, public atmos_pres
module urban / dynamics / Kusaka01
real(rp), public mapprojection_basepoint_lon
real(rp), dimension(:,:), allocatable, public urban_rainr
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.
module Atmosphere GRID CartesC Real(real space)
module Map projection
real(rp), dimension(:,:), allocatable, public atmos_temp
real(rp), dimension(:,:), allocatable, public urban_q2
module urban / grid / cartesianC / real
real(rp), dimension(:,:), allocatable, public urban_tg
real(rp), dimension(:,:), allocatable, public atmos_u
character(len=h_short), public urban_sfc_type
real(rp), dimension(:,:), allocatable, public urban_sfc_temp
module Statistics
integer, dimension(6), public time_nowdate
current time [YYYY MM DD HH MM SS]
Definition: scale_time.F90:69
real(rp), dimension(:,:,:), allocatable, public urban_trl
real(rp), dimension(:,:), allocatable, public urban_rainb_t
real(rp), dimension(:,:,:), allocatable, public atmos_sflx_lw
real(rp), dimension(:,:), allocatable, public atmos_sflx_snow
module Urban admin
real(rp), dimension(:,:,:), allocatable, public urban_tbl
real(dp), public time_dtsec_urban
time interval of urban step [sec]
Definition: scale_time.F90:52
real(rp), public mapprojection_basepoint_lat
real(rp), dimension(:,:), allocatable, public urban_tr_t
real(rp), dimension(:,:,:), allocatable, public urban_tgl_t
real(rp), dimension(:,:), allocatable, public urban_sflx_gh
real(rp), dimension(:,:), allocatable, public urban_sflx_mv
real(rp), dimension(:,:), allocatable, public atmos_w
real(rp), dimension(:,:), allocatable, public urban_rainb
real(rp), dimension(:,:), allocatable, public urban_tc_t
real(rp), dimension(:,:), allocatable, public urban_roff_t
real(rp), dimension(:,:,:), allocatable, public urban_tbl_t
module file_history
real(rp), dimension(:,:), allocatable, public urban_sflx_mw
real(rp), dimension(:,:), allocatable, public atmos_sfc_dens
Here is the call graph for this function:
Here is the caller graph for this function:

◆ urban_driver_update()

subroutine, public mod_urban_driver::urban_driver_update ( )

Urban step.

Definition at line 452 of file mod_urban_driver.F90.

References mod_urban_vars::atmos_cossza, mod_urban_vars::atmos_dens, mod_urban_vars::atmos_pbl, mod_urban_vars::atmos_pres, mod_urban_vars::atmos_qv, mod_urban_vars::atmos_sfc_dens, mod_urban_vars::atmos_sfc_pres, mod_urban_vars::atmos_sflx_lw, mod_urban_vars::atmos_sflx_rain, mod_urban_vars::atmos_sflx_snow, mod_urban_vars::atmos_sflx_sw, mod_urban_vars::atmos_temp, mod_urban_vars::atmos_u, mod_urban_vars::atmos_v, mod_urban_vars::atmos_w, scale_landuse::landuse_fact_urban, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), scale_time::time_dtsec_urban, scale_urban_grid_cartesc_index::uie, scale_urban_grid_cartesc_index::uis, scale_urban_grid_cartesc_index::uje, scale_urban_grid_cartesc_index::ujs, scale_urban_grid_cartesc_index::uke, scale_urban_grid_cartesc_index::uks, mod_urban_admin::urban_dyn_type, mod_urban_vars::urban_q2, mod_urban_vars::urban_qc, mod_urban_vars::urban_qc_t, mod_urban_vars::urban_rainb, mod_urban_vars::urban_rainb_t, mod_urban_vars::urban_raing, mod_urban_vars::urban_raing_t, mod_urban_vars::urban_rainr, mod_urban_vars::urban_rainr_t, mod_urban_vars::urban_roff, mod_urban_vars::urban_roff_t, mod_urban_vars::urban_sfc_temp, mod_urban_vars::urban_sflx_gh, mod_urban_vars::urban_sflx_lh, mod_urban_vars::urban_sflx_mu, mod_urban_vars::urban_sflx_mv, mod_urban_vars::urban_sflx_mw, mod_urban_vars::urban_sflx_qtrc, mod_urban_vars::urban_sflx_sh, urban_surface_get(), mod_urban_vars::urban_t2, mod_urban_vars::urban_tb, mod_urban_vars::urban_tb_t, mod_urban_vars::urban_tbl, mod_urban_vars::urban_tbl_t, mod_urban_vars::urban_tc, mod_urban_vars::urban_tc_t, mod_urban_vars::urban_tg, mod_urban_vars::urban_tg_t, mod_urban_vars::urban_tgl, mod_urban_vars::urban_tgl_t, mod_urban_vars::urban_tr, mod_urban_vars::urban_tr_t, mod_urban_vars::urban_trl, mod_urban_vars::urban_trl_t, mod_urban_vars::urban_u10, mod_urban_vars::urban_uc, mod_urban_vars::urban_uc_t, mod_urban_vars::urban_v10, mod_urban_vars::urban_vars_history(), mod_urban_vars::urban_vars_total(), mod_urban_vars::urban_z0e, mod_urban_vars::urban_z0h, and mod_urban_vars::urban_z0m.

Referenced by mod_rm_driver::rm_driver().

452  use scale_time, only: &
453  dt => time_dtsec_urban
454  use mod_urban_vars, only: &
455  atmos_temp, &
456  atmos_pres, &
457  atmos_w, &
458  atmos_u, &
459  atmos_v, &
460  atmos_dens, &
461  atmos_qv, &
462  atmos_pbl, &
463  atmos_sfc_dens, &
464  atmos_sfc_pres, &
465  atmos_sflx_lw, &
466  atmos_sflx_sw, &
467  atmos_cossza, &
468  atmos_sflx_rain, &
469  atmos_sflx_snow, &
470  urban_trl_t, &
471  urban_tbl_t, &
472  urban_tgl_t, &
473  urban_tr_t, &
474  urban_tb_t, &
475  urban_tg_t, &
476  urban_tc_t, &
477  urban_qc_t, &
478  urban_uc_t, &
479  urban_rainr_t, &
480  urban_rainb_t, &
481  urban_raing_t, &
482  urban_roff_t, &
483  urban_sflx_mw, &
484  urban_sflx_mu, &
485  urban_sflx_mv, &
486  urban_sflx_sh, &
487  urban_sflx_lh, &
488  urban_sflx_qtrc, &
489  urban_sflx_gh, &
490  urban_z0m, &
491  urban_z0h, &
492  urban_z0e, &
493  urban_u10, &
494  urban_v10, &
495  urban_t2, &
496  urban_q2, &
497  urban_tr, &
498  urban_tb, &
499  urban_tg, &
500  urban_tc, &
501  urban_qc, &
502  urban_uc, &
503  urban_trl, &
504  urban_tbl, &
505  urban_tgl, &
506  urban_rainr, &
507  urban_rainb, &
508  urban_raing, &
509  urban_roff, &
510  urban_sfc_temp, &
513  use scale_landuse, only: &
515  use mod_urban_admin, only: &
517  implicit none
518 
519  integer :: k, i, j
520  !---------------------------------------------------------------------------
521 
522  call prof_rapstart('URB_Update', 1)
523 
524  !########## Get Surface Boundary from coupler ##########
525  call urban_surface_get
526 
527  !########## Dynamics / Update variables ##########
528  select case ( urban_dyn_type )
529  case ( 'KUSAKA01' )
530 
531 !OCL XFILL
532  !$omp parallel do
533  do j = ujs, uje
534  do i = uis, uie
535  do k = uks, uke
536  urban_trl(k,i,j) = urban_trl(k,i,j) + urban_trl_t(k,i,j) * dt
537  urban_tbl(k,i,j) = urban_tbl(k,i,j) + urban_tbl_t(k,i,j) * dt
538  urban_tgl(k,i,j) = urban_tgl(k,i,j) + urban_tgl_t(k,i,j) * dt
539  end do
540  end do
541  end do
542 
543 !OCL XFILL
544  !$omp parallel do
545  do j = ujs, uje
546  do i = uis, uie
547  urban_tr(i,j) = urban_tr(i,j) + urban_tr_t(i,j) * dt
548  urban_tb(i,j) = urban_tb(i,j) + urban_tb_t(i,j) * dt
549  urban_tg(i,j) = urban_tg(i,j) + urban_tg_t(i,j) * dt
550  urban_tc(i,j) = urban_tc(i,j) + urban_tc_t(i,j) * dt
551  urban_qc(i,j) = urban_qc(i,j) + urban_qc_t(i,j) * dt
552  urban_uc(i,j) = urban_uc(i,j) + urban_uc_t(i,j) * dt
553  urban_rainr(i,j) = urban_rainr(i,j) + urban_rainr_t(i,j) * dt
554  urban_rainb(i,j) = urban_rainb(i,j) + urban_rainb_t(i,j) * dt
555  urban_raing(i,j) = urban_raing(i,j) + urban_raing_t(i,j) * dt
556  urban_roff(i,j) = urban_roff(i,j) + urban_roff_t(i,j) * dt
557  end do
558  end do
559 
560  end select
561 
562  call urban_vars_total
563 
564  !########## History & Monitor ##########
565  call urban_vars_history
566 
567 
568  call prof_rapend ('URB_Update', 1)
569 
570  return
real(rp), dimension(:,:), allocatable, public urban_qc_t
real(rp), dimension(:,:,:), allocatable, public urban_sflx_qtrc
real(rp), dimension(:,:), allocatable, public urban_qc
real(rp), dimension(:,:), allocatable, public urban_rainr_t
real(rp), dimension(:,:), allocatable, public landuse_fact_urban
urban factor
real(rp), dimension(:,:,:), allocatable, public urban_trl_t
real(rp), dimension(:,:), allocatable, public urban_tg_t
real(rp), dimension(:,:), allocatable, public urban_u10
real(rp), dimension(:,:), allocatable, public urban_v10
real(rp), dimension(:,:), allocatable, public urban_raing_t
real(rp), dimension(:,:), allocatable, public urban_z0e
real(rp), dimension(:,:), allocatable, public urban_sflx_mu
real(rp), dimension(:,:), allocatable, public urban_z0m
subroutine, public urban_vars_total
Budget monitor for urban.
real(rp), dimension(:,:), allocatable, public urban_tb_t
real(rp), dimension(:,:), allocatable, public urban_t2
subroutine, public urban_vars_history
History output set for urban variables.
real(rp), dimension(:,:), allocatable, public urban_tb
module URBAN Variables
real(rp), dimension(:,:), allocatable, public urban_raing
real(rp), dimension(:,:), allocatable, public urban_uc
real(rp), dimension(:,:), allocatable, public atmos_pbl
real(rp), dimension(:,:), allocatable, public urban_sflx_sh
real(rp), dimension(:,:), allocatable, public urban_tr
real(rp), dimension(:,:,:), allocatable, public urban_tgl
module LANDUSE
real(rp), dimension(:,:), allocatable, public atmos_qv
real(rp), dimension(:,:), allocatable, public atmos_sfc_pres
real(rp), dimension(:,:), allocatable, public atmos_cossza
real(rp), dimension(:,:), allocatable, public urban_uc_t
real(rp), dimension(:,:), allocatable, public urban_z0h
real(rp), dimension(:,:,:), allocatable, public atmos_sflx_sw
real(rp), dimension(:,:), allocatable, public atmos_v
module TIME
Definition: scale_time.F90:16
real(rp), dimension(:,:), allocatable, public urban_roff
real(rp), dimension(:,:), allocatable, public atmos_dens
real(rp), dimension(:,:), allocatable, public atmos_sflx_rain
real(rp), dimension(:,:), allocatable, public urban_sflx_lh
real(rp), dimension(:,:), allocatable, public urban_tc
real(rp), dimension(:,:), allocatable, public atmos_pres
real(rp), dimension(:,:), allocatable, public urban_rainr
character(len=h_short), public urban_dyn_type
real(rp), dimension(:,:), allocatable, public atmos_temp
real(rp), dimension(:,:), allocatable, public urban_q2
real(rp), dimension(:,:), allocatable, public urban_tg
real(rp), dimension(:,:), allocatable, public atmos_u
real(rp), dimension(:,:), allocatable, public urban_sfc_temp
real(rp), dimension(:,:,:), allocatable, public urban_trl
real(rp), dimension(:,:), allocatable, public urban_rainb_t
real(rp), dimension(:,:,:), allocatable, public atmos_sflx_lw
real(rp), dimension(:,:), allocatable, public atmos_sflx_snow
module Urban admin
real(rp), dimension(:,:,:), allocatable, public urban_tbl
real(dp), public time_dtsec_urban
time interval of urban step [sec]
Definition: scale_time.F90:52
real(rp), dimension(:,:), allocatable, public urban_tr_t
real(rp), dimension(:,:,:), allocatable, public urban_tgl_t
real(rp), dimension(:,:), allocatable, public urban_sflx_gh
real(rp), dimension(:,:), allocatable, public urban_sflx_mv
real(rp), dimension(:,:), allocatable, public atmos_w
real(rp), dimension(:,:), allocatable, public urban_rainb
real(rp), dimension(:,:), allocatable, public urban_tc_t
real(rp), dimension(:,:), allocatable, public urban_roff_t
real(rp), dimension(:,:,:), allocatable, public urban_tbl_t
real(rp), dimension(:,:), allocatable, public urban_sflx_mw
real(rp), dimension(:,:), allocatable, public atmos_sfc_dens
Here is the call graph for this function:
Here is the caller graph for this function:

◆ urban_surface_get()

subroutine, public mod_urban_driver::urban_surface_get ( )

Get surface boundary.

Definition at line 576 of file mod_urban_driver.F90.

References mod_urban_vars::atmos_cossza, mod_urban_vars::atmos_dens, mod_urban_vars::atmos_pbl, mod_urban_vars::atmos_pres, mod_urban_vars::atmos_qv, mod_urban_vars::atmos_sfc_dens, mod_urban_vars::atmos_sfc_pres, mod_urban_vars::atmos_sflx_lw, mod_urban_vars::atmos_sflx_rain, mod_urban_vars::atmos_sflx_snow, mod_urban_vars::atmos_sflx_sw, mod_urban_vars::atmos_temp, mod_urban_vars::atmos_u, mod_urban_vars::atmos_v, mod_urban_vars::atmos_w, mod_cpl_vars::cpl_getatm_urb(), 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_prof::prof_rapend(), scale_prof::prof_rapstart(), scale_urban_grid_cartesc_index::uie, scale_urban_grid_cartesc_index::uis, scale_urban_grid_cartesc_index::uje, scale_urban_grid_cartesc_index::ujs, and mod_urban_admin::urban_do.

Referenced by urban_driver_calc_tendency(), and urban_driver_update().

576  use mod_urban_admin, only: &
577  urban_do
578  use mod_urban_vars, only: &
579  atmos_temp, &
580  atmos_pres, &
581  atmos_w, &
582  atmos_u, &
583  atmos_v, &
584  atmos_dens, &
585  atmos_qv, &
586  atmos_pbl, &
587  atmos_sfc_dens, &
588  atmos_sfc_pres, &
589  atmos_sflx_lw, &
590  atmos_sflx_sw, &
591  atmos_cossza, &
592  atmos_sflx_rain, &
594  use mod_cpl_vars, only: &
596  implicit none
597 
598  real(RP) :: atmos_sflx_rad_dn(uia,uja,n_rad_dir,n_rad_rgn)
599 
600  integer :: i, j
601  !---------------------------------------------------------------------------
602 
603  call prof_rapstart('URB_SfcExch', 2)
604 
605  if ( urban_do ) then
606  call cpl_getatm_urb( atmos_temp(:,:), & ! [OUT]
607  atmos_pres(:,:), & ! [OUT]
608  atmos_w(:,:), & ! [OUT]
609  atmos_u(:,:), & ! [OUT]
610  atmos_v(:,:), & ! [OUT]
611  atmos_dens(:,:), & ! [OUT]
612  atmos_qv(:,:), & ! [OUT]
613  atmos_pbl(:,:), & ! [OUT]
614  atmos_sfc_dens(:,:), & ! [OUT]
615  atmos_sfc_pres(:,:), & ! [OUT]
616  atmos_sflx_rad_dn(:,:,:,:), & ! [OUT]
617  atmos_cossza(:,:), & ! [OUT]
618  atmos_sflx_rain(:,:), & ! [OUT]
619  atmos_sflx_snow(:,:) ) ! [OUT]
620  endif
621 
622 !OCL XFILL
623  do j = ujs, uje
624  do i = uis, uie
625  atmos_sflx_lw(i,j,i_r_direct ) = atmos_sflx_rad_dn(i,j,i_r_direct ,i_r_ir) ! IR, direct
626  atmos_sflx_lw(i,j,i_r_diffuse) = atmos_sflx_rad_dn(i,j,i_r_diffuse,i_r_ir) ! IR, diffuse
627 
628  atmos_sflx_sw(i,j,i_r_direct ) = atmos_sflx_rad_dn(i,j,i_r_direct ,i_r_nir) & ! NIR, direct
629  + atmos_sflx_rad_dn(i,j,i_r_direct ,i_r_vis) ! VIS, direct
630  atmos_sflx_sw(i,j,i_r_diffuse) = atmos_sflx_rad_dn(i,j,i_r_diffuse,i_r_nir) & ! NIR, diffuse
631  + atmos_sflx_rad_dn(i,j,i_r_diffuse,i_r_vis) ! VIS, diffuse
632  enddo
633  enddo
634 
635  call prof_rapend ('URB_SfcExch', 2)
636 
637  return
module URBAN Variables
real(rp), dimension(:,:), allocatable, public atmos_pbl
real(rp), dimension(:,:), allocatable, public atmos_qv
real(rp), dimension(:,:), allocatable, public atmos_sfc_pres
real(rp), dimension(:,:), allocatable, public atmos_cossza
module COUPLER Variables
real(rp), dimension(:,:,:), allocatable, public atmos_sflx_sw
real(rp), dimension(:,:), allocatable, public atmos_v
real(rp), dimension(:,:), allocatable, public atmos_dens
real(rp), dimension(:,:), allocatable, public atmos_sflx_rain
real(rp), dimension(:,:), allocatable, public atmos_pres
real(rp), dimension(:,:), allocatable, public atmos_temp
real(rp), dimension(:,:), allocatable, public atmos_u
real(rp), dimension(:,:,:), allocatable, public atmos_sflx_lw
real(rp), dimension(:,:), allocatable, public atmos_sflx_snow
module Urban admin
logical, public urban_do
subroutine, public cpl_getatm_urb(TEMP, PRES, W, U, V, DENS, QV, PBL, SFC_DENS, SFC_PRES, SFLX_rad_dn, cosSZA, SFLX_rain, SFLX_snow)
real(rp), dimension(:,:), allocatable, public atmos_w
real(rp), dimension(:,:), allocatable, public atmos_sfc_dens
Here is the call graph for this function:
Here is the caller graph for this function:

◆ urban_surface_set()

subroutine, public mod_urban_driver::urban_surface_set ( logical, intent(in)  countup)

Set surface boundary to other model.

Definition at line 643 of file mod_urban_driver.F90.

References mod_cpl_vars::cpl_puturb(), scale_prof::prof_rapend(), scale_prof::prof_rapstart(), mod_urban_admin::urban_do, mod_urban_vars::urban_q2, mod_urban_vars::urban_sfc_albedo, mod_urban_vars::urban_sfc_temp, mod_urban_vars::urban_sflx_gh, mod_urban_vars::urban_sflx_lh, mod_urban_vars::urban_sflx_mu, mod_urban_vars::urban_sflx_mv, mod_urban_vars::urban_sflx_mw, mod_urban_vars::urban_sflx_qtrc, mod_urban_vars::urban_sflx_sh, mod_urban_vars::urban_t2, mod_urban_vars::urban_u10, mod_urban_vars::urban_v10, mod_urban_vars::urban_z0e, mod_urban_vars::urban_z0h, and mod_urban_vars::urban_z0m.

Referenced by mod_rm_driver::restart_read(), mod_rm_prep::rm_prep(), and urban_driver_calc_tendency().

643  use mod_urban_admin, only: &
644  urban_do
645  use mod_urban_vars, only: &
646  urban_sfc_temp, &
648  urban_sflx_mw, &
649  urban_sflx_mu, &
650  urban_sflx_mv, &
651  urban_sflx_sh, &
652  urban_sflx_lh, &
653  urban_sflx_gh, &
654  urban_sflx_qtrc, &
655  urban_z0m, &
656  urban_z0h, &
657  urban_z0e, &
658  urban_u10, &
659  urban_v10, &
660  urban_t2, &
661  urban_q2
662  use mod_cpl_vars, only: &
663  cpl_puturb
664  implicit none
665 
666  ! arguments
667  logical, intent(in) :: countup
668  !---------------------------------------------------------------------------
669 
670  call prof_rapstart('URB_SfcExch', 2)
671 
672  if ( urban_do ) then
673  call cpl_puturb( urban_sfc_temp(:,:), & ! [IN]
674  urban_sfc_albedo(:,:,:,:), & ! [IN]
675  urban_z0m(:,:), & ! [IN]
676  urban_z0h(:,:), & ! [IN]
677  urban_z0e(:,:), & ! [IN]
678  urban_sflx_mw(:,:), & ! [IN]
679  urban_sflx_mu(:,:), & ! [IN]
680  urban_sflx_mv(:,:), & ! [IN]
681  urban_sflx_sh(:,:), & ! [IN]
682  urban_sflx_lh(:,:), & ! [IN]
683  urban_sflx_gh(:,:), & ! [IN]
684  urban_sflx_qtrc(:,:,:), & ! [IN]
685  urban_u10(:,:), & ! [IN]
686  urban_v10(:,:), & ! [IN]
687  urban_t2(:,:), & ! [IN]
688  urban_q2(:,:), & ! [IN]
689  countup ) ! [IN]
690  endif
691 
692  call prof_rapend ('URB_SfcExch', 2)
693 
694  return
real(rp), dimension(:,:,:), allocatable, public urban_sflx_qtrc
real(rp), dimension(:,:), allocatable, public urban_u10
real(rp), dimension(:,:), allocatable, public urban_v10
real(rp), dimension(:,:), allocatable, public urban_z0e
real(rp), dimension(:,:), allocatable, public urban_sflx_mu
real(rp), dimension(:,:), allocatable, public urban_z0m
real(rp), dimension(:,:), allocatable, public urban_t2
module URBAN Variables
real(rp), dimension(:,:,:,:), allocatable, public urban_sfc_albedo
real(rp), dimension(:,:), allocatable, public urban_sflx_sh
module COUPLER Variables
real(rp), dimension(:,:), allocatable, public urban_z0h
subroutine, public cpl_puturb(SFC_TEMP, SFC_albedo, SFC_Z0M, SFC_Z0H, SFC_Z0E, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_LH, SFLX_G, SFLX_QTRC, U10, V10, T2, Q2, countup)
real(rp), dimension(:,:), allocatable, public urban_sflx_lh
real(rp), dimension(:,:), allocatable, public urban_q2
real(rp), dimension(:,:), allocatable, public urban_sfc_temp
module Urban admin
logical, public urban_do
real(rp), dimension(:,:), allocatable, public urban_sflx_gh
real(rp), dimension(:,:), allocatable, public urban_sflx_mv
real(rp), dimension(:,:), allocatable, public urban_sflx_mw
Here is the call graph for this function:
Here is the caller graph for this function: