SCALE-RM
mod_urban_driver.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
9 !-------------------------------------------------------------------------------
10 #include "scalelib.h"
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_io
18  use scale_prof
20  use scale_tracer
22  !-----------------------------------------------------------------------------
23  implicit none
24  private
25  !-----------------------------------------------------------------------------
26  !
27  !++ Public procedure
28  !
29  public :: urban_driver_setup
31  public :: urban_driver_update
32  public :: urban_surface_get
33  public :: urban_surface_set
34 
35  !-----------------------------------------------------------------------------
36  !
37  !++ Public parameters & variables
38  !
39  !-----------------------------------------------------------------------------
40  !
41  !++ Private procedure
42  !
43  !-----------------------------------------------------------------------------
44  !
45  !++ Private parameters & variables
46  !
47  real(RP), private, allocatable :: AH_URB (:,:,:) ! urban grid average of anthropogenic sensible heat [W/m2]
48  real(RP), private, allocatable :: AHL_URB(:,:,:) ! urban grid average of anthropogenic latent heat [W/m2]
49  !-----------------------------------------------------------------------------
50 contains
51  !-----------------------------------------------------------------------------
53  subroutine urban_driver_setup
54  use scale_prc, only: &
55  prc_abort
56  use scale_const, only: &
57  undef => const_undef
58  use mod_urban_admin, only: &
59  urban_do, &
62  use scale_urban_dyn_kusaka01, only: &
64  use scale_landuse, only: &
66  use mod_urban_vars, only: &
67  urban_z0m, &
68  urban_z0h, &
69  urban_z0e, &
70  urban_zd
71  implicit none
72  !---------------------------------------------------------------------------
73 
74  log_newline
75  log_info("URBAN_driver_setup",*) 'Setup'
76 
77  if ( urban_do ) then
78 
79  allocate( ah_urb(uia,uja,1:24) )
80  allocate( ahl_urb(uia,uja,1:24) )
81  ah_urb(:,:,:) = undef
82  ahl_urb(:,:,:) = undef
83 
84  select case ( urban_dyn_type )
85  case ( 'KUSAKA01' )
86  call urban_dyn_kusaka01_setup( uia, uis, uie, uja, ujs, uje, & ! [IN]
87  landuse_fact_urban(:,:), & ! [IN]
88  urban_z0m(:,:), urban_z0h(:,:), urban_z0e(:,:), & ! [OUT]
89  urban_zd(:,:), & ! [OUT]
90  ah_urb(:,:,:), ahl_urb(:,:,:) ) ! [OUT]
91 
92  urban_sfc_type = 'KUSAKA01'
93  case default
94  log_error("URBAN_driver_setup",*) 'LAND_DYN_TYPE is invalid: ', trim(urban_dyn_type)
95  call prc_abort
96  end select
97 
98  select case ( urban_sfc_type )
99  case ( 'KUSAKA01' )
100  ! do nothing
101  case default
102  log_error("URBAN_driver_setup",*) 'LAND_SFC_TYPE is invalid: ', trim(urban_sfc_type)
103  call prc_abort
104  end select
105 
106  end if
107 
108  return
109  end subroutine urban_driver_setup
110 
111  !-----------------------------------------------------------------------------
113  subroutine urban_driver_calc_tendency( force )
114  use scale_const, only: &
115  d2r => const_d2r
116  use scale_statistics, only: &
118  statistics_total
119  use scale_urban_grid_cartesc_real, only: &
124  use scale_topography, only: &
125  tansl_x => topography_tansl_x, &
126  tansl_y => topography_tansl_y
127  use scale_file_history, only: &
128  file_history_in
129  use mod_atmos_admin, only: &
131  use mod_atmos_phy_ch_driver, only: &
133  use mod_urban_vars, only: &
134  atmos_temp, &
135  atmos_pres, &
136  atmos_u, &
137  atmos_v, &
138  atmos_dens, &
139  atmos_qv, &
140  atmos_pbl, &
141  atmos_sfc_dens, &
142  atmos_sfc_pres, &
143  atmos_sflx_lw, &
144  atmos_sflx_sw, &
145  atmos_cossza, &
147  atmos_sflx_engi, &
148  urban_trl_t, &
149  urban_tbl_t, &
150  urban_tgl_t, &
151  urban_tr_t, &
152  urban_tb_t, &
153  urban_tg_t, &
154  urban_tc_t, &
155  urban_qc_t, &
156  urban_uc_t, &
157  urban_rainr_t, &
158  urban_rainb_t, &
159  urban_raing_t, &
160  urban_sfc_temp, &
162  urban_sflx_mw, &
163  urban_sflx_mu, &
164  urban_sflx_mv, &
165  urban_sflx_sh, &
166  urban_sflx_lh, &
167  urban_sflx_shex, &
168  urban_sflx_qvex, &
169  urban_sflx_qtrc, &
170  urban_sflx_gh, &
171  urban_z0m, &
172  urban_z0h, &
173  urban_z0e, &
174  urban_zd, &
175  urban_ah, &
176  urban_ahl, &
177  urban_ustar, &
178  urban_tstar, &
179  urban_qstar, &
180  urban_wstar, &
181  urban_rlmo, &
182  urban_u10, &
183  urban_v10, &
184  urban_t2, &
185  urban_q2, &
186  urban_tr, &
187  urban_tb, &
188  urban_tg, &
189  urban_tc, &
190  urban_qc, &
191  urban_uc, &
192  urban_trl, &
193  urban_tbl, &
194  urban_tgl, &
195  urban_rainr, &
196  urban_rainb, &
197  urban_raing, &
198  urban_roff
199  use scale_atmos_hydrometeor, only: &
200  hydrometeor_lhv => atmos_hydrometeor_lhv, &
202  i_qv
203  use scale_time, only: &
204  dt => time_dtsec_urban, &
205  nowdate => time_nowdate
206  use scale_mapprojection, only: &
207  base_lon => mapprojection_basepoint_lon, &
208  base_lat => mapprojection_basepoint_lat
209  use scale_atmos_grid_cartesc_real, only: &
210  real_z1 => atmos_grid_cartesc_real_z1
211  use scale_urban_grid_cartesc, only: &
213  use scale_landuse, only: &
215  exists_urban => landuse_exists_urban
216  use mod_urban_admin, only: &
218  use scale_urban_dyn_kusaka01, only: &
220  implicit none
221  logical, intent(in) :: force
222 
223  real(rp) :: trl(uka,uia,uja), tbl(uka,uia,uja), tgl(uka,uia,uja)
224  real(rp) :: tr(uia,uja), tb(uia,uja), tg(uia,uja)
225  real(rp) :: tc(uia,uja), qc(uia,uja), uc(uia,uja)
226  real(rp) :: rainr(uia,uja), rainb(uia,uja), raing(uia,uja), roff(uia,uja)
227 
228  real(rp) :: lhv(uia,uja) ! latent heat of vaporization [J/kg]
229 
230  real(rp) :: urban_sflx_lhex(uia,uja)
231 
232  real(rp) :: lat, lon ! [deg]
233  integer :: tloc ! local time (1-24h)
234  real(rp) :: dsec ! second [s]
235 
236  integer :: k, i, j, iq
237  !---------------------------------------------------------------------------
238 
239  call prof_rapstart('URB_CalcTend', 1)
240 
241  !########## Get Surface Boundary from coupler ##########
242  call urban_surface_get
243 
244  !########## initialize tendency ##########
245 !OCL XFILL
246  do j = ujs, uje
247  do i = uis, uie
248  do k = uks, uke
249  urban_trl_t(k,i,j) = 0.0_rp
250  urban_tbl_t(k,i,j) = 0.0_rp
251  urban_tgl_t(k,i,j) = 0.0_rp
252  end do
253  end do
254  end do
255 
256 !OCL XFILL
257  do j = ujs, uje
258  do i = uis, uie
259  urban_tr_t(i,j) = 0.0_rp
260  urban_tb_t(i,j) = 0.0_rp
261  urban_tg_t(i,j) = 0.0_rp
262  urban_tc_t(i,j) = 0.0_rp
263  urban_qc_t(i,j) = 0.0_rp
264  urban_uc_t(i,j) = 0.0_rp
265 
266  urban_rainr_t(i,j) = 0.0_rp
267  urban_rainb_t(i,j) = 0.0_rp
268  urban_raing_t(i,j) = 0.0_rp
269  enddo
270  enddo
271 
272 !OCL XFILL
273  !$omp parallel do
274  do iq = 1, qa
275  do j = ujs, uje
276  do i = uis, uie
277  urban_sflx_qtrc(i,j,iq) = 0.0_rp
278  enddo
279  enddo
280  enddo
281 
282  select case ( urban_sfc_type )
283  case ( 'KUSAKA01' )
284 
285 !OCL XFILL
286  !$omp parallel do
287  do j = ujs, uje
288  do i = uis, uie
289  do k = uks, uke
290  trl(k,i,j) = urban_trl(k,i,j)
291  tbl(k,i,j) = urban_tbl(k,i,j)
292  tgl(k,i,j) = urban_tgl(k,i,j)
293  end do
294  end do
295  end do
296 
297 !OCL XFILL
298  !$omp parallel do
299  do j = ujs, uje
300  do i = uis, uie
301  tr(i,j) = urban_tr(i,j)
302  tb(i,j) = urban_tb(i,j)
303  tg(i,j) = urban_tg(i,j)
304  tc(i,j) = urban_tc(i,j)
305  qc(i,j) = urban_qc(i,j)
306  uc(i,j) = urban_uc(i,j)
307  rainr(i,j) = urban_rainr(i,j)
308  rainb(i,j) = urban_rainb(i,j)
309  raing(i,j) = urban_raing(i,j)
310  end do
311  end do
312 
313 
314  ! local time
315  lat = base_lat
316  lon = base_lon
317  if (lon < 0.0_rp ) lon = mod(lon, 360.0_rp) + 360.0_rp
318  if (lon > 360.0_rp ) lon = mod(lon, 360.0_rp)
319  tloc = mod( (nowdate(4) + int(lon/15.0_rp)),24 )
320  dsec = real( nowdate(5)*60.0_rp + nowdate(6), kind=rp ) / 3600.0_rp
321  if( tloc == 0 ) tloc = 24
322 
323  !--- Calculate AH at LST
324  if ( tloc == 24 ) then
325  do j = ujs, uje
326  do i = uis, uie
327  if ( exists_urban(i,j) ) then
328  urban_ah(i,j) = ( 1.0_rp-dsec ) * ah_urb(i,j,tloc ) &
329  + ( dsec ) * ah_urb(i,j,1 )
330  urban_ahl(i,j) = ( 1.0_rp-dsec ) * ahl_urb(i,j,tloc ) &
331  + ( dsec ) * ahl_urb(i,j,1 )
332  end if
333  enddo
334  enddo
335  else
336  do j = ujs, uje
337  do i = uis, uie
338  if ( exists_urban(i,j) ) then
339  urban_ah(i,j) = ( 1.0_rp-dsec ) * ah_urb(i,j,tloc ) &
340  + ( dsec ) * ah_urb(i,j,tloc+1)
341  urban_ahl(i,j) = ( 1.0_rp-dsec ) * ahl_urb(i,j,tloc ) &
342  + ( dsec ) * ahl_urb(i,j,tloc+1)
343  end if
344  enddo
345  enddo
346  endif
347 
348  call hydrometeor_lhv( uia, uis, uie, uja, ujs, uje, &
349  atmos_temp(:,:), lhv(:,:) )
350 
351  call urban_dyn_kusaka01( uka, uks, uke, uia, uis, uie, uja, ujs, uje, &
352  atmos_temp(:,:), atmos_pres(:,:), & ! [IN]
353  atmos_u(:,:), atmos_v(:,:), & ! [IN]
354  atmos_dens(:,:), atmos_qv(:,:), lhv(:,:), & ! [IN]
355  real_z1(:,:), atmos_pbl(:,:), & ! [IN]
356  atmos_sfc_dens(:,:), atmos_sfc_pres(:,:), & ! [IN]
357  atmos_sflx_lw(:,:,:), atmos_sflx_sw(:,:,:), & ! [IN]
358  atmos_sflx_water(:,:), atmos_sflx_engi(:,:), & ! [IN]
359  urban_z0m(:,:), urban_z0h(:,:), urban_z0e(:,:), & ! [IN]
360  urban_zd(:,:), & ! [IN]
361  cdz(:), & ! [IN]
362  tansl_x(:,:), tansl_y(:,:), & ! [IN]
363  landuse_fact_urban(:,:), & ! [IN]
364  dt, & ! [IN]
365  trl(:,:,:), tbl(:,:,:), tgl(:,:,:), & ! [INOUT]
366  tr(:,:), tb(:,:), tg(:,:), tc(:,:), qc(:,:), uc(:,:), & ! [INOUT]
367  rainr(:,:), rainb(:,:), raing(:,:), urban_roff(:,:), & ! [INOUT]
368  urban_sfc_temp(:,:), & ! [OUT]
369  urban_sfc_albedo(:,:,:,:), & ! [OUT]
370  urban_sflx_mw(:,:), urban_sflx_mu(:,:), urban_sflx_mv(:,:), & ! [OUT]
371  urban_sflx_sh(:,:), urban_sflx_lh(:,:), urban_sflx_gh(:,:), & ! [OUT]
372  urban_ustar(:,:), urban_tstar(:,:), urban_qstar(:,:), & ! [OUT]
373  urban_wstar(:,:), & ! [OUT]
374  urban_rlmo(:,:), & ! [OUT]
375  urban_u10(:,:), urban_v10(:,:), urban_t2(:,:), urban_q2(:,:) ) ! [OUT]
376 
377  !-----------------------------------------------------------
378  ! anthropogenic heat fluxes
379  !-----------------------------------------------------------
380  !$omp parallel do
381  do j = ujs, uje
382  do i = uis, uie
383  if ( exists_urban(i,j) ) then
384  urban_sflx_shex(i,j) = urban_ah(i,j) / landuse_fact_urban(i,j) ! Sensible heat flux [W/m2]
385  urban_sflx_lhex(i,j) = urban_ahl(i,j) / landuse_fact_urban(i,j) ! Latent heat flux [W/m2]
386  urban_sflx_sh(i,j) = urban_sflx_sh(i,j) + urban_sflx_shex(i,j)
387  urban_sflx_lh(i,j) = urban_sflx_lh(i,j) + urban_sflx_lhex(i,j)
388  end if
389  end do
390  end do
391 
392 !OCL XFILL
393  !$omp parallel do
394  do j = ujs, uje
395  do i = uis, uie
396  if ( exists_urban(i,j) ) then
397  do k = uks, uke
398  urban_trl_t(k,i,j) = ( trl(k,i,j) - urban_trl(k,i,j) ) / dt
399  urban_tbl_t(k,i,j) = ( tbl(k,i,j) - urban_tbl(k,i,j) ) / dt
400  urban_tgl_t(k,i,j) = ( tgl(k,i,j) - urban_tgl(k,i,j) ) / dt
401  end do
402  end if
403  end do
404  end do
405 
406 !OCL XFILL
407  !$omp parallel do
408  do j = ujs, uje
409  do i = uis, uie
410  if ( exists_urban(i,j) ) then
411  urban_tr_t(i,j) = ( tr(i,j) - urban_tr(i,j) ) / dt
412  urban_tb_t(i,j) = ( tb(i,j) - urban_tb(i,j) ) / dt
413  urban_tg_t(i,j) = ( tg(i,j) - urban_tg(i,j) ) / dt
414  urban_tc_t(i,j) = ( tc(i,j) - urban_tc(i,j) ) / dt
415  urban_qc_t(i,j) = ( qc(i,j) - urban_qc(i,j) ) / dt
416  urban_uc_t(i,j) = ( uc(i,j) - urban_uc(i,j) ) / dt
417  urban_rainr_t(i,j) = ( rainr(i,j) - urban_rainr(i,j) ) / dt
418  urban_rainb_t(i,j) = ( rainb(i,j) - urban_rainb(i,j) ) / dt
419  urban_raing_t(i,j) = ( raing(i,j) - urban_raing(i,j) ) / dt
420  end if
421  end do
422  end do
423 
424  if ( .NOT. atmos_hydrometeor_dry ) then
425  !$omp parallel do
426  do j = ujs, uje
427  do i = uis, uie
428  if ( exists_urban(i,j) ) then
429  urban_sflx_qtrc(i,j,i_qv) = urban_sflx_lh(i,j) / lhv(i,j)
430  urban_sflx_qvex(i,j) = urban_sflx_lhex(i,j) / lhv(i,j)
431  end if
432  enddo
433  enddo
434  endif
435 
436  end select
437 
438  ! Surface flux for chemical tracers
439  if ( atmos_sw_phy_ch ) then
440  call atmos_phy_ch_driver_urban_flux( urban_sflx_qtrc(:,:,:) ) ! [INOUT]
441  endif
442 
443  call file_history_in( urban_tr_t(:,:), 'URBAN_TR_t', 'tendency of URBAN_TR', 'K/s', dim_type='XY' )
444  call file_history_in( urban_tb_t(:,:), 'URBAN_TB_t', 'tendency of URBAN_TB', 'K/s', dim_type='XY' )
445  call file_history_in( urban_tg_t(:,:), 'URBAN_TG_t', 'tendency of URBAN_TG', 'K/s', dim_type='XY' )
446  call file_history_in( urban_tc_t(:,:), 'URBAN_TC_t', 'tendency of URBAN_TC', 'K/s', dim_type='XY' )
447  call file_history_in( urban_qc_t(:,:), 'URBAN_QC_t', 'tendency of URBAN_QC', 'kg/kg/s', dim_type='XY' )
448  call file_history_in( urban_uc_t(:,:), 'URBAN_UC_t', 'tendency of URBAN_UC', 'm/s2', dim_type='XY' )
449 
450  call file_history_in( urban_trl_t(:,:,:), 'URBAN_TRL_t', 'tendency of URBAN_TRL', 'K/s', dim_type='UXY' )
451  call file_history_in( urban_tbl_t(:,:,:), 'URBAN_TBL_t', 'tendency of URBAN_TBL', 'K/s', dim_type='UXY' )
452  call file_history_in( urban_tgl_t(:,:,:), 'URBAN_TGL_t', 'tendency of URBAN_TGL', 'K/s', dim_type='UXY' )
453 
454  call file_history_in( urban_rainr_t(:,:), 'URBAN_RAINR_t', 'tendency of URBAN_RAINR', 'kg/m2/s', dim_type='XY' )
455  call file_history_in( urban_rainb_t(:,:), 'URBAN_RAINB_t', 'tendency of URBAN_RAINB', 'kg/m2/s', dim_type='XY' )
456  call file_history_in( urban_raing_t(:,:), 'URBAN_RAING_t', 'tendency of URBAN_RAING', 'kg/m2/s', dim_type='XY' )
457  call file_history_in( urban_roff(:,:), 'URBAN_ROFF', 'urban runoff water', 'kg/m2/s', dim_type='XY' )
458 
459  if ( statistics_checktotal ) then
460 
461  call statistics_total( uka, uks, uke, uia, uis, uie, uja, ujs, uje, &
462  urban_trl_t(:,:,:), 'URBAN_TRL_t', &
465  call statistics_total( uka, uks, uke, uia, uis, uie, uja, ujs, uje, &
466  urban_tbl_t(:,:,:), 'URBAN_TBL_t', &
469  call statistics_total( uka, uks, uke, uia, uis, uie, uja, ujs, uje, &
470  urban_tgl_t(:,:,:), 'URBAN_TGL_t', &
473 
474  call statistics_total( uia, uis, uie, uja, ujs, uje, &
475  urban_tr_t(:,:), 'URBAN_TR_t', &
478  call statistics_total( uia, uis, uie, uja, ujs, uje, &
479  urban_tb_t(:,:), 'URBAN_TB_t', &
482  call statistics_total( uia, uis, uie, uja, ujs, uje, &
483  urban_tg_t(:,:), 'URBAN_TG_t', &
486  call statistics_total( uia, uis, uie, uja, ujs, uje, &
487  urban_tc_t(:,:), 'URBAN_TC_t', &
490  call statistics_total( uia, uis, uie, uja, ujs, uje, &
491  urban_qc_t(:,:), 'URBAN_QC_t', &
494  call statistics_total( uia, uis, uie, uja, ujs, uje, &
495  urban_uc_t(:,:), 'URBAN_UC_t', &
498 
499  call statistics_total( uia, uis, uie, uja, ujs, uje, &
500  urban_rainr_t(:,:), 'URBAN_RAINR_t', &
503  call statistics_total( uia, uis, uie, uja, ujs, uje, &
504  urban_rainb_t(:,:), 'URBAN_RAINB_t', &
507  call statistics_total( uia, uis, uie, uja, ujs, uje, &
508  urban_raing_t(:,:), 'URBAN_RAING_t', &
511  call statistics_total( uia, uis, uie, uja, ujs, uje, &
512  urban_roff(:,:), 'URBAN_ROFF', &
515  endif
516 
517 
518  call prof_rapend ('URB_CalcTend', 1)
519 
520  !########## Set Surface Boundary to coupler ##########
521  call urban_surface_set( countup=.true. )
522 
523  return
524  end subroutine urban_driver_calc_tendency
525 
526  !-----------------------------------------------------------------------------
528  subroutine urban_driver_update
529  use scale_time, only: &
530  dt => time_dtsec_urban
531  use mod_urban_vars, only: &
532  urban_trl_t, &
533  urban_tbl_t, &
534  urban_tgl_t, &
535  urban_tr_t, &
536  urban_tb_t, &
537  urban_tg_t, &
538  urban_tc_t, &
539  urban_qc_t, &
540  urban_uc_t, &
541  urban_rainr_t, &
542  urban_rainb_t, &
543  urban_raing_t, &
544  urban_tr, &
545  urban_tb, &
546  urban_tg, &
547  urban_tc, &
548  urban_qc, &
549  urban_uc, &
550  urban_trl, &
551  urban_tbl, &
552  urban_tgl, &
553  urban_rainr, &
554  urban_rainb, &
555  urban_raing, &
557  use mod_urban_admin, only: &
559  use scale_landuse, only: &
560  exists_urban => landuse_exists_urban
561  implicit none
562 
563  integer :: k, i, j
564  !---------------------------------------------------------------------------
565 
566  call prof_rapstart('URB_Update', 1)
567 
568  !########## Get Surface Boundary from coupler ##########
569  call urban_surface_get
570 
571  !########## Dynamics / Update variables ##########
572  select case ( urban_dyn_type )
573  case ( 'KUSAKA01' )
574 
575 !OCL XFILL
576  !$omp parallel do
577  do j = ujs, uje
578  do i = uis, uie
579  if ( exists_urban(i,j) ) then
580  do k = uks, uke
581  urban_trl(k,i,j) = urban_trl(k,i,j) + urban_trl_t(k,i,j) * dt
582  urban_tbl(k,i,j) = urban_tbl(k,i,j) + urban_tbl_t(k,i,j) * dt
583  urban_tgl(k,i,j) = urban_tgl(k,i,j) + urban_tgl_t(k,i,j) * dt
584  end do
585  end if
586  end do
587  end do
588 
589 !OCL XFILL
590  !$omp parallel do
591  do j = ujs, uje
592  do i = uis, uie
593  if ( exists_urban(i,j) ) then
594  urban_tr(i,j) = urban_tr(i,j) + urban_tr_t(i,j) * dt
595  urban_tb(i,j) = urban_tb(i,j) + urban_tb_t(i,j) * dt
596  urban_tg(i,j) = urban_tg(i,j) + urban_tg_t(i,j) * dt
597  urban_tc(i,j) = urban_tc(i,j) + urban_tc_t(i,j) * dt
598  urban_qc(i,j) = max( real(urban_qc(i,j) + urban_qc_t(i,j) * dt, rp), 0.0_rp )
599  urban_uc(i,j) = max( real(urban_uc(i,j) + urban_uc_t(i,j) * dt, rp), 0.0_rp )
600  urban_rainr(i,j) = max( real(urban_rainr(i,j) + urban_rainr_t(i,j) * dt, rp), 0.0_rp )
601  urban_rainb(i,j) = max( real(urban_rainb(i,j) + urban_rainb_t(i,j) * dt, rp), 0.0_rp )
602  urban_raing(i,j) = max( real(urban_raing(i,j) + urban_raing_t(i,j) * dt, rp), 0.0_rp )
603  end if
604  end do
605  end do
606 
607  end select
608 
609  call urban_vars_check
610 
611  call prof_rapend ('URB_Update', 1)
612 
613  return
614  end subroutine urban_driver_update
615 
616  !-----------------------------------------------------------------------------
618  subroutine urban_surface_get
619  use mod_urban_admin, only: &
620  urban_do
621  use mod_urban_vars, only: &
622  atmos_temp, &
623  atmos_pres, &
624  atmos_w, &
625  atmos_u, &
626  atmos_v, &
627  atmos_dens, &
628  atmos_qv, &
629  atmos_pbl, &
630  atmos_sfc_dens, &
631  atmos_sfc_pres, &
632  atmos_sflx_lw, &
633  atmos_sflx_sw, &
634  atmos_cossza, &
637  use mod_cpl_vars, only: &
639  implicit none
640 
641  real(rp) :: atmos_sflx_rad_dn(uia,uja,n_rad_dir,n_rad_rgn)
642 
643  integer :: i, j
644  !---------------------------------------------------------------------------
645 
646  call prof_rapstart('URB_SfcExch', 2)
647 
648  if ( urban_do ) then
649  call cpl_getatm_urb( atmos_temp(:,:), & ! [OUT]
650  atmos_pres(:,:), & ! [OUT]
651  atmos_w(:,:), & ! [OUT]
652  atmos_u(:,:), & ! [OUT]
653  atmos_v(:,:), & ! [OUT]
654  atmos_dens(:,:), & ! [OUT]
655  atmos_qv(:,:), & ! [OUT]
656  atmos_pbl(:,:), & ! [OUT]
657  atmos_sfc_dens(:,:), & ! [OUT]
658  atmos_sfc_pres(:,:), & ! [OUT]
659  atmos_sflx_rad_dn(:,:,:,:), & ! [OUT]
660  atmos_cossza(:,:), & ! [OUT]
661  atmos_sflx_water(:,:), & ! [OUT]
662  atmos_sflx_engi(:,:) ) ! [OUT]
663  endif
664 
665 !OCL XFILL
666  do j = ujs, uje
667  do i = uis, uie
668  atmos_sflx_lw(i,j,i_r_direct ) = atmos_sflx_rad_dn(i,j,i_r_direct ,i_r_ir) ! IR, direct
669  atmos_sflx_lw(i,j,i_r_diffuse) = atmos_sflx_rad_dn(i,j,i_r_diffuse,i_r_ir) ! IR, diffuse
670 
671  atmos_sflx_sw(i,j,i_r_direct ) = atmos_sflx_rad_dn(i,j,i_r_direct ,i_r_nir) & ! NIR, direct
672  + atmos_sflx_rad_dn(i,j,i_r_direct ,i_r_vis) ! VIS, direct
673  atmos_sflx_sw(i,j,i_r_diffuse) = atmos_sflx_rad_dn(i,j,i_r_diffuse,i_r_nir) & ! NIR, diffuse
674  + atmos_sflx_rad_dn(i,j,i_r_diffuse,i_r_vis) ! VIS, diffuse
675  enddo
676  enddo
677 
678  call prof_rapend ('URB_SfcExch', 2)
679 
680  return
681  end subroutine urban_surface_get
682 
683  !-----------------------------------------------------------------------------
685  subroutine urban_surface_set( countup )
686  use mod_urban_admin, only: &
687  urban_do
688  use mod_urban_vars, only: &
689  urban_sfc_temp, &
691  urban_sflx_mw, &
692  urban_sflx_mu, &
693  urban_sflx_mv, &
694  urban_sflx_sh, &
695  urban_sflx_lh, &
696  urban_sflx_shex, &
697  urban_sflx_qvex, &
698  urban_sflx_gh, &
699  urban_sflx_qtrc, &
700  urban_z0m, &
701  urban_z0h, &
702  urban_z0e, &
703  urban_u10, &
704  urban_v10, &
705  urban_t2, &
706  urban_q2
707  use mod_cpl_vars, only: &
708  cpl_puturb
709  use scale_landuse, only: &
710  exists_urban => landuse_exists_urban
711  implicit none
712 
713  ! arguments
714  logical, intent(in) :: countup
715  !---------------------------------------------------------------------------
716 
717  call prof_rapstart('URB_SfcExch', 2)
718 
719  if ( urban_do ) then
720  call cpl_puturb( urban_sfc_temp(:,:), & ! [IN]
721  urban_sfc_albedo(:,:,:,:), & ! [IN]
722  urban_z0m(:,:), & ! [IN]
723  urban_z0h(:,:), & ! [IN]
724  urban_z0e(:,:), & ! [IN]
725  urban_sflx_mw(:,:), & ! [IN]
726  urban_sflx_mu(:,:), & ! [IN]
727  urban_sflx_mv(:,:), & ! [IN]
728  urban_sflx_sh(:,:), & ! [IN]
729  urban_sflx_lh(:,:), & ! [IN]
730  urban_sflx_shex(:,:), & ! [IN]
731  urban_sflx_qvex(:,:), & ! [IN]
732  urban_sflx_gh(:,:), & ! [IN]
733  urban_sflx_qtrc(:,:,:), & ! [IN]
734  urban_u10(:,:), & ! [IN]
735  urban_v10(:,:), & ! [IN]
736  urban_t2(:,:), & ! [IN]
737  urban_q2(:,:), & ! [IN]
738  exists_urban(:,:), & ! [IN]
739  countup ) ! [IN]
740  endif
741 
742  call prof_rapend ('URB_SfcExch', 2)
743 
744  return
745  end subroutine urban_surface_set
746 
747 end module mod_urban_driver
mod_urban_vars::urban_q2
real(rp), dimension(:,:), allocatable, public urban_q2
Definition: mod_urban_vars.F90:120
mod_urban_vars::urban_trl_t
real(rp), dimension(:,:,:), allocatable, public urban_trl_t
Definition: mod_urban_vars.F90:78
mod_urban_vars::atmos_sflx_sw
real(rp), dimension(:,:,:), allocatable, public atmos_sflx_sw
Definition: mod_urban_vars.F90:134
scale_cpl_sfc_index::n_rad_dir
integer, parameter, public n_rad_dir
Definition: scale_cpl_sfc_index.F90:36
scale_statistics
module Statistics
Definition: scale_statistics.F90:11
mod_urban_vars::urban_z0m
real(rp), dimension(:,:), allocatable, public urban_z0m
Definition: mod_urban_vars.F90:104
mod_urban_vars::atmos_u
real(rp), dimension(:,:), allocatable, public atmos_u
Definition: mod_urban_vars.F90:126
scale_urban_grid_cartesc_index::uja
integer, public uja
Definition: scale_urban_grid_cartesC_index.F90:44
mod_urban_admin::urban_do
logical, public urban_do
Definition: mod_urban_admin.F90:32
mod_urban_vars::urban_tb_t
real(rp), dimension(:,:), allocatable, public urban_tb_t
Definition: mod_urban_vars.F90:85
mod_urban_vars::urban_sfc_temp
real(rp), dimension(:,:), allocatable, public urban_sfc_temp
Definition: mod_urban_vars.F90:74
scale_urban_grid_cartesc_real::urban_grid_cartesc_real_vol
real(rp), dimension(:,:,:), allocatable, public urban_grid_cartesc_real_vol
volume of grid cell
Definition: scale_urban_grid_cartesC_real.F90:38
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_cpl_sfc_index::i_r_direct
integer, parameter, public i_r_direct
Definition: scale_cpl_sfc_index.F90:37
scale_urban_grid_cartesc_real::urban_grid_cartesc_real_totvol
real(rp), public urban_grid_cartesc_real_totvol
total volume
Definition: scale_urban_grid_cartesC_real.F90:39
scale_tracer::qa
integer, public qa
Definition: scale_tracer.F90:34
scale_urban_grid_cartesc::urban_grid_cartesc_cdz
real(rp), dimension(:), allocatable, public urban_grid_cartesc_cdz
z-length of control volume [m]
Definition: scale_urban_grid_cartesC.F90:37
mod_urban_vars::urban_tb
real(rp), dimension(:,:), allocatable, public urban_tb
Definition: mod_urban_vars.F90:64
scale_urban_grid_cartesc_index::uia
integer, public uia
Definition: scale_urban_grid_cartesC_index.F90:40
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_z1
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_z1
Height of the lowermost grid from surface (cell center) [m].
Definition: scale_atmos_grid_cartesC_real.F90:59
mod_urban_vars::urban_raing_t
real(rp), dimension(:,:), allocatable, public urban_raing_t
Definition: mod_urban_vars.F90:89
scale_cpl_sfc_index::i_r_diffuse
integer, parameter, public i_r_diffuse
Definition: scale_cpl_sfc_index.F90:38
scale_topography::topography_tansl_y
real(rp), dimension(:,:), allocatable, public topography_tansl_y
tan(slope_y)
Definition: scale_topography.F90:40
mod_cpl_vars
module COUPLER Variables
Definition: mod_cpl_vars.F90:12
scale_precision
module PRECISION
Definition: scale_precision.F90:14
mod_atmos_admin
module ATMOS admin
Definition: mod_atmos_admin.F90:11
mod_urban_vars::urban_sfc_albedo
real(rp), dimension(:,:,:,:), allocatable, public urban_sfc_albedo
Definition: mod_urban_vars.F90:75
scale_urban_dyn_kusaka01
module urban / dynamics / Kusaka01
Definition: scale_urban_dyn_kusaka01.F90:12
mod_urban_vars::atmos_sflx_engi
real(rp), dimension(:,:), allocatable, public atmos_sflx_engi
Definition: mod_urban_vars.F90:137
scale_cpl_sfc_index::i_r_ir
integer, parameter, public i_r_ir
Definition: scale_cpl_sfc_index.F90:29
scale_prof::prof_rapstart
subroutine, public prof_rapstart(rapname_base, level, disable_barrier)
Start raptime.
Definition: scale_prof.F90:159
mod_urban_vars::urban_trl
real(rp), dimension(:,:,:), allocatable, public urban_trl
Definition: mod_urban_vars.F90:60
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
mod_urban_vars::urban_sflx_shex
real(rp), dimension(:,:), allocatable, public urban_sflx_shex
Definition: mod_urban_vars.F90:98
mod_urban_vars::urban_sflx_mu
real(rp), dimension(:,:), allocatable, public urban_sflx_mu
Definition: mod_urban_vars.F90:94
mod_cpl_vars::cpl_puturb
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_SHEX, SFLX_QVEX, SFLX_GH, SFLX_QTRC, U10, V10, T2, Q2, mask, countup)
Definition: mod_cpl_vars.F90:884
mod_urban_vars::urban_raing
real(rp), dimension(:,:), allocatable, public urban_raing
Definition: mod_urban_vars.F90:71
mod_urban_vars::urban_z0h
real(rp), dimension(:,:), allocatable, public urban_z0h
Definition: mod_urban_vars.F90:105
scale_topography
module TOPOGRAPHY
Definition: scale_topography.F90:11
scale_urban_dyn_kusaka01::urban_dyn_kusaka01_setup
subroutine, public urban_dyn_kusaka01_setup(UIA, UIS, UIE, UJA, UJS, UJE, fact_urban, Z0M, Z0H, Z0E, ZD, AH_URB, AHL_URB)
Setup.
Definition: scale_urban_dyn_kusaka01.F90:125
mod_urban_vars::urban_qc
real(rp), dimension(:,:), allocatable, public urban_qc
Definition: mod_urban_vars.F90:67
scale_atmos_hydrometeor::atmos_hydrometeor_dry
logical, public atmos_hydrometeor_dry
Definition: scale_atmos_hydrometeor.F90:97
scale_urban_grid_cartesc_index
module urban / grid / icosahedralA / index
Definition: scale_urban_grid_cartesC_index.F90:11
mod_urban_vars::atmos_pbl
real(rp), dimension(:,:), allocatable, public atmos_pbl
Definition: mod_urban_vars.F90:130
mod_urban_vars::urban_tbl_t
real(rp), dimension(:,:,:), allocatable, public urban_tbl_t
Definition: mod_urban_vars.F90:79
mod_urban_vars::urban_tgl_t
real(rp), dimension(:,:,:), allocatable, public urban_tgl_t
Definition: mod_urban_vars.F90:80
scale_atmos_grid_cartesc_real
module Atmosphere GRID CartesC Real(real space)
Definition: scale_atmos_grid_cartesC_real.F90:11
mod_urban_driver::urban_driver_update
subroutine, public urban_driver_update
Urban step.
Definition: mod_urban_driver.F90:529
mod_urban_vars::atmos_w
real(rp), dimension(:,:), allocatable, public atmos_w
Definition: mod_urban_vars.F90:125
scale_file_history
module file_history
Definition: scale_file_history.F90:15
mod_urban_vars::atmos_v
real(rp), dimension(:,:), allocatable, public atmos_v
Definition: mod_urban_vars.F90:127
mod_urban_driver::urban_driver_calc_tendency
subroutine, public urban_driver_calc_tendency(force)
Calclate tendency.
Definition: mod_urban_driver.F90:114
mod_urban_vars::atmos_sflx_water
real(rp), dimension(:,:), allocatable, public atmos_sflx_water
Definition: mod_urban_vars.F90:136
mod_urban_vars::atmos_dens
real(rp), dimension(:,:), allocatable, public atmos_dens
Definition: mod_urban_vars.F90:128
mod_urban_vars::urban_qstar
real(rp), dimension(:,:), allocatable, public urban_qstar
Definition: mod_urban_vars.F90:114
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_mapprojection::mapprojection_basepoint_lat
real(rp), public mapprojection_basepoint_lat
Definition: scale_mapprojection.F90:91
mod_urban_vars::urban_sflx_qtrc
real(rp), dimension(:,:,:), allocatable, public urban_sflx_qtrc
Definition: mod_urban_vars.F90:100
scale_mapprojection
module Map projection
Definition: scale_mapprojection.F90:12
scale_io
module STDIO
Definition: scale_io.F90:10
mod_urban_vars::urban_wstar
real(rp), dimension(:,:), allocatable, public urban_wstar
Definition: mod_urban_vars.F90:115
mod_urban_vars::urban_t2
real(rp), dimension(:,:), allocatable, public urban_t2
Definition: mod_urban_vars.F90:119
mod_urban_admin::urban_dyn_type
character(len=h_short), public urban_dyn_type
Definition: mod_urban_admin.F90:35
scale_urban_dyn_kusaka01::urban_dyn_kusaka01
subroutine, public urban_dyn_kusaka01(UKA, UKS, UKE, UIA, UIS, UIE, UJA, UJS, UJE, TMPA, PRSA, U1, V1, DENS, QA, LHV, Z1, PBL, 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.
Definition: scale_urban_dyn_kusaka01.F90:379
scale_cpl_sfc_index::i_r_nir
integer, parameter, public i_r_nir
Definition: scale_cpl_sfc_index.F90:30
scale_tracer::k
real(rp), public k
Definition: scale_tracer.F90:44
scale_const
module CONSTANT
Definition: scale_const.F90:11
mod_urban_vars::urban_sflx_mv
real(rp), dimension(:,:), allocatable, public urban_sflx_mv
Definition: mod_urban_vars.F90:95
mod_urban_vars::urban_tgl
real(rp), dimension(:,:,:), allocatable, public urban_tgl
Definition: mod_urban_vars.F90:62
scale_urban_grid_cartesc
module urban / grid / cartesianC
Definition: scale_urban_grid_cartesC.F90:12
mod_urban_vars::urban_v10
real(rp), dimension(:,:), allocatable, public urban_v10
Definition: mod_urban_vars.F90:118
mod_urban_vars::atmos_qv
real(rp), dimension(:,:), allocatable, public atmos_qv
Definition: mod_urban_vars.F90:129
mod_urban_vars::atmos_sfc_pres
real(rp), dimension(:,:), allocatable, public atmos_sfc_pres
Definition: mod_urban_vars.F90:132
scale_topography::topography_tansl_x
real(rp), dimension(:,:), allocatable, public topography_tansl_x
tan(slope_x)
Definition: scale_topography.F90:39
mod_urban_vars::atmos_sflx_lw
real(rp), dimension(:,:,:), allocatable, public atmos_sflx_lw
Definition: mod_urban_vars.F90:133
mod_urban_vars::urban_z0e
real(rp), dimension(:,:), allocatable, public urban_z0e
Definition: mod_urban_vars.F90:106
mod_cpl_vars::cpl_getatm_urb
subroutine, public cpl_getatm_urb(TEMP, PRES, W, U, V, DENS, QV, PBL, SFC_DENS, SFC_PRES, SFLX_rad_dn, cosSZA, SFLX_water, SFLX_ENGI)
Definition: mod_cpl_vars.F90:1286
mod_urban_vars::urban_u10
real(rp), dimension(:,:), allocatable, public urban_u10
Definition: mod_urban_vars.F90:117
mod_urban_driver
module URBAN driver
Definition: mod_urban_driver.F90:11
mod_urban_vars::atmos_cossza
real(rp), dimension(:,:), allocatable, public atmos_cossza
Definition: mod_urban_vars.F90:135
mod_urban_driver::urban_surface_get
subroutine, public urban_surface_get
Get surface boundary.
Definition: mod_urban_driver.F90:619
scale_prof
module profiler
Definition: scale_prof.F90:11
mod_urban_driver::urban_surface_set
subroutine, public urban_surface_set(countup)
Set surface boundary to other model.
Definition: mod_urban_driver.F90:686
mod_atmos_phy_ch_driver::atmos_phy_ch_driver_urban_flux
subroutine, public atmos_phy_ch_driver_urban_flux(SFLX_QTRC)
Driver.
Definition: mod_atmos_phy_ch_driver.F90:288
mod_atmos_admin::atmos_sw_phy_ch
logical, public atmos_sw_phy_ch
Definition: mod_atmos_admin.F90:54
scale_mapprojection::mapprojection_basepoint_lon
real(rp), public mapprojection_basepoint_lon
Definition: scale_mapprojection.F90:90
scale_time
module TIME
Definition: scale_time.F90:11
mod_urban_vars::urban_vars_check
subroutine, public urban_vars_check(force)
Budget monitor for urban.
Definition: mod_urban_vars.F90:940
mod_urban_vars::atmos_pres
real(rp), dimension(:,:), allocatable, public atmos_pres
Definition: mod_urban_vars.F90:124
mod_urban_vars::urban_tstar
real(rp), dimension(:,:), allocatable, public urban_tstar
Definition: mod_urban_vars.F90:113
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
scale_time::time_dtsec_urban
real(dp), public time_dtsec_urban
time interval of urban step [sec]
Definition: scale_time.F90:48
mod_urban_vars::urban_tc
real(rp), dimension(:,:), allocatable, public urban_tc
Definition: mod_urban_vars.F90:66
mod_urban_vars::urban_tg
real(rp), dimension(:,:), allocatable, public urban_tg
Definition: mod_urban_vars.F90:65
scale_atmos_hydrometeor::i_qv
integer, public i_qv
Definition: scale_atmos_hydrometeor.F90:77
scale_urban_grid_cartesc_index::uks
integer, public uks
Definition: scale_urban_grid_cartesC_index.F90:37
mod_urban_vars::urban_tg_t
real(rp), dimension(:,:), allocatable, public urban_tg_t
Definition: mod_urban_vars.F90:86
mod_urban_vars::urban_uc
real(rp), dimension(:,:), allocatable, public urban_uc
Definition: mod_urban_vars.F90:68
mod_urban_vars::urban_tr_t
real(rp), dimension(:,:), allocatable, public urban_tr_t
Definition: mod_urban_vars.F90:84
mod_atmos_phy_ch_driver
module ATMOSPHERE / Physics Chemistry
Definition: mod_atmos_phy_ch_driver.F90:12
mod_urban_driver::urban_driver_setup
subroutine, public urban_driver_setup
Setup.
Definition: mod_urban_driver.F90:54
scale_statistics::statistics_checktotal
logical, public statistics_checktotal
calc&report variable totals to logfile?
Definition: scale_statistics.F90:64
mod_urban_vars::urban_sflx_lh
real(rp), dimension(:,:), allocatable, public urban_sflx_lh
Definition: mod_urban_vars.F90:97
mod_urban_admin::urban_sfc_type
character(len=h_short), public urban_sfc_type
Definition: mod_urban_admin.F90:39
scale_landuse::landuse_exists_urban
logical, dimension(:,:), allocatable, public landuse_exists_urban
urban calculation flag
Definition: scale_landuse.F90:51
scale_cpl_sfc_index
module coupler / surface-atmospehre
Definition: scale_cpl_sfc_index.F90:11
scale_time::time_nowdate
integer, dimension(6), public time_nowdate
current time [YYYY MM DD HH MM SS]
Definition: scale_time.F90:66
mod_urban_vars::urban_rainb
real(rp), dimension(:,:), allocatable, public urban_rainb
Definition: mod_urban_vars.F90:70
scale_cpl_sfc_index::i_r_vis
integer, parameter, public i_r_vis
Definition: scale_cpl_sfc_index.F90:31
mod_urban_vars::urban_ah
real(rp), dimension(:,:), allocatable, public urban_ah
Definition: mod_urban_vars.F90:108
mod_urban_vars::urban_ahl
real(rp), dimension(:,:), allocatable, public urban_ahl
Definition: mod_urban_vars.F90:109
mod_urban_vars::urban_sflx_sh
real(rp), dimension(:,:), allocatable, public urban_sflx_sh
Definition: mod_urban_vars.F90:96
mod_urban_vars::urban_qc_t
real(rp), dimension(:,:), allocatable, public urban_qc_t
Definition: mod_urban_vars.F90:83
mod_urban_vars::urban_tr
real(rp), dimension(:,:), allocatable, public urban_tr
Definition: mod_urban_vars.F90:63
scale_urban_grid_cartesc_index::uka
integer, public uka
Definition: scale_urban_grid_cartesC_index.F90:36
mod_urban_vars::urban_roff
real(rp), dimension(:,:), allocatable, public urban_roff
Definition: mod_urban_vars.F90:91
mod_urban_vars::urban_rainb_t
real(rp), dimension(:,:), allocatable, public urban_rainb_t
Definition: mod_urban_vars.F90:88
scale_urban_grid_cartesc_index::uis
integer, public uis
Definition: scale_urban_grid_cartesC_index.F90:41
mod_urban_vars::urban_ustar
real(rp), dimension(:,:), allocatable, public urban_ustar
Definition: mod_urban_vars.F90:112
mod_urban_vars::urban_zd
real(rp), dimension(:,:), allocatable, public urban_zd
Definition: mod_urban_vars.F90:107
scale_prof::prof_rapend
subroutine, public prof_rapend(rapname_base, level, disable_barrier)
Save raptime.
Definition: scale_prof.F90:217
scale_const::const_d2r
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:32
mod_urban_vars::urban_sflx_qvex
real(rp), dimension(:,:), allocatable, public urban_sflx_qvex
Definition: mod_urban_vars.F90:99
scale_landuse::landuse_fact_urban
real(rp), dimension(:,:), allocatable, public landuse_fact_urban
urban factor
Definition: scale_landuse.F90:46
scale_urban_grid_cartesc_real::urban_grid_cartesc_real_totarea
real(rp), public urban_grid_cartesc_real_totarea
total area
Definition: scale_urban_grid_cartesC_real.F90:37
mod_urban_vars::urban_uc_t
real(rp), dimension(:,:), allocatable, public urban_uc_t
Definition: mod_urban_vars.F90:82
mod_urban_vars::urban_rainr_t
real(rp), dimension(:,:), allocatable, public urban_rainr_t
Definition: mod_urban_vars.F90:87
mod_urban_vars::urban_tbl
real(rp), dimension(:,:,:), allocatable, public urban_tbl
Definition: mod_urban_vars.F90:61
scale_landuse
module LANDUSE
Definition: scale_landuse.F90:19
mod_urban_admin
module Urban admin
Definition: mod_urban_admin.F90:11
mod_urban_vars::urban_sflx_mw
real(rp), dimension(:,:), allocatable, public urban_sflx_mw
Definition: mod_urban_vars.F90:93
scale_cpl_sfc_index::n_rad_rgn
integer, parameter, public n_rad_rgn
Definition: scale_cpl_sfc_index.F90:28
mod_urban_vars::atmos_sfc_dens
real(rp), dimension(:,:), allocatable, public atmos_sfc_dens
Definition: mod_urban_vars.F90:131
scale_urban_grid_cartesc_real
module urban / grid / cartesianC / real
Definition: scale_urban_grid_cartesC_real.F90:12
mod_urban_vars::urban_sflx_gh
real(rp), dimension(:,:), allocatable, public urban_sflx_gh
Definition: mod_urban_vars.F90:101
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:41
scale_urban_grid_cartesc_index::uje
integer, public uje
Definition: scale_urban_grid_cartesC_index.F90:46
mod_urban_vars::urban_rlmo
real(rp), dimension(:,:), allocatable, public urban_rlmo
Definition: mod_urban_vars.F90:116
scale_urban_grid_cartesc_index::ujs
integer, public ujs
Definition: scale_urban_grid_cartesC_index.F90:45
mod_urban_vars::urban_rainr
real(rp), dimension(:,:), allocatable, public urban_rainr
Definition: mod_urban_vars.F90:69
mod_urban_vars
module URBAN Variables
Definition: mod_urban_vars.F90:12
scale_urban_grid_cartesc_index::uie
integer, public uie
Definition: scale_urban_grid_cartesC_index.F90:42
scale_urban_grid_cartesc_index::uke
integer, public uke
Definition: scale_urban_grid_cartesC_index.F90:38
mod_urban_vars::atmos_temp
real(rp), dimension(:,:), allocatable, public atmos_temp
Definition: mod_urban_vars.F90:123
mod_urban_vars::urban_tc_t
real(rp), dimension(:,:), allocatable, public urban_tc_t
Definition: mod_urban_vars.F90:81
scale_urban_grid_cartesc_real::urban_grid_cartesc_real_area
real(rp), dimension(:,:), allocatable, public urban_grid_cartesc_real_area
area of grid cell
Definition: scale_urban_grid_cartesC_real.F90:36