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
30  public :: urban_driver_finalize
32  public :: urban_driver_update
33  public :: urban_surface_get
34  public :: urban_surface_set
35 
36  !-----------------------------------------------------------------------------
37  !
38  !++ Public parameters & variables
39  !
40  !-----------------------------------------------------------------------------
41  !
42  !++ Private procedure
43  !
44  !-----------------------------------------------------------------------------
45  !
46  !++ Private parameters & variables
47  !
48  real(RP), private, allocatable :: AH_URB (:,:,:) ! urban grid average of anthropogenic sensible heat [W/m2]
49  real(RP), private, allocatable :: AHL_URB(:,:,:) ! urban grid average of anthropogenic latent heat [W/m2]
50  real(RP), private :: AH_TOFFSET ! time offset for AH [Hour]
51  !-----------------------------------------------------------------------------
52 contains
53  !-----------------------------------------------------------------------------
55  subroutine urban_driver_setup
56  use scale_prc, only: &
57  prc_abort
58  use scale_const, only: &
59  undef => const_undef
60  use mod_urban_admin, only: &
61  urban_do, &
64  use scale_urban_dyn_kusaka01, only: &
66  use scale_landuse, only: &
68  use mod_urban_vars, only: &
69  urban_z0m, &
70  urban_z0h, &
71  urban_z0e, &
72  urban_zd
73  implicit none
74  !---------------------------------------------------------------------------
75 
76  log_newline
77  log_info("URBAN_driver_setup",*) 'Setup'
78 
79  if ( urban_do ) then
80 
81  allocate( ah_urb(uia,uja,1:24) )
82  allocate( ahl_urb(uia,uja,1:24) )
83  ah_urb(:,:,:) = undef
84  ahl_urb(:,:,:) = undef
85  !$acc enter data create(AH_URB,AHL_URB)
86 
87  select case ( urban_dyn_type )
88  case ( 'KUSAKA01' )
89  call urban_dyn_kusaka01_setup( uia, uis, uie, uja, ujs, uje, & ! [IN]
90  landuse_fact_urban(:,:), & ! [IN]
91  urban_z0m(:,:), urban_z0h(:,:), urban_z0e(:,:), & ! [OUT]
92  urban_zd(:,:), & ! [OUT]
93  ah_urb(:,:,:), ahl_urb(:,:,:), ah_toffset ) ! [OUT]
94 
95  urban_sfc_type = 'KUSAKA01'
96  case default
97  log_error("URBAN_driver_setup",*) 'LAND_DYN_TYPE is invalid: ', trim(urban_dyn_type)
98  call prc_abort
99  end select
100 
101  select case ( urban_sfc_type )
102  case ( 'KUSAKA01' )
103  ! do nothing
104  case default
105  log_error("URBAN_driver_setup",*) 'LAND_SFC_TYPE is invalid: ', trim(urban_sfc_type)
106  call prc_abort
107  end select
108 
109  end if
110 
111  return
112  end subroutine urban_driver_setup
113 
114  !-----------------------------------------------------------------------------
116  subroutine urban_driver_finalize
117  use mod_urban_admin, only: &
118  urban_do, &
119  urban_dyn_type, &
121  use scale_urban_dyn_kusaka01, only: &
123  implicit none
124  !---------------------------------------------------------------------------
125 
126  log_newline
127  log_info("URBAN_driver_finalize",*) 'Finalize'
128 
129  if ( urban_do ) then
130 
131  select case ( urban_dyn_type )
132  case ( 'KUSAKA01' )
134  end select
135 
136  select case ( urban_sfc_type )
137  case ( 'KUSAKA01' )
138  end select
139 
140  !$acc exit data delete(AH_URB,AHL_URB)
141  deallocate( ah_urb )
142  deallocate( ahl_urb )
143 
144  end if
145 
146  return
147  end subroutine urban_driver_finalize
148 
149  !-----------------------------------------------------------------------------
151  subroutine urban_driver_calc_tendency( force )
152  use scale_const, only: &
153  d2r => const_d2r
154  use scale_statistics, only: &
156  statistics_total
157  use scale_urban_grid_cartesc_real, only: &
162  use scale_topography, only: &
163  tansl_x => topography_tansl_x, &
164  tansl_y => topography_tansl_y
165  use scale_file_history, only: &
166  file_history_in
167  use mod_atmos_admin, only: &
169  use mod_atmos_phy_ch_driver, only: &
171  use mod_urban_vars, only: &
172  atmos_temp, &
173  atmos_pres, &
174  atmos_u, &
175  atmos_v, &
176  atmos_dens, &
177  atmos_qv, &
178  atmos_sfc_dens, &
179  atmos_sfc_pres, &
180  atmos_sflx_lw, &
181  atmos_sflx_sw, &
182  atmos_cossza, &
184  atmos_sflx_engi, &
185  urban_trl_t, &
186  urban_tbl_t, &
187  urban_tgl_t, &
188  urban_tr_t, &
189  urban_tb_t, &
190  urban_tg_t, &
191  urban_tc_t, &
192  urban_qc_t, &
193  urban_uc_t, &
194  urban_rainr_t, &
195  urban_rainb_t, &
196  urban_raing_t, &
197  urban_sfc_temp, &
199  urban_sflx_mw, &
200  urban_sflx_mu, &
201  urban_sflx_mv, &
202  urban_sflx_sh, &
203  urban_sflx_lh, &
204  urban_sflx_shex, &
205  urban_sflx_lhex, &
206  urban_sflx_qvex, &
207  urban_sflx_qtrc, &
208  urban_sflx_gh, &
209  urban_z0m, &
210  urban_z0h, &
211  urban_z0e, &
212  urban_zd, &
213  urban_ah, &
214  urban_ahl, &
215  urban_ustar, &
216  urban_tstar, &
217  urban_qstar, &
218  urban_wstar, &
219  urban_rlmo, &
220  urban_u10, &
221  urban_v10, &
222  urban_t2, &
223  urban_q2, &
224  urban_tr, &
225  urban_tb, &
226  urban_tg, &
227  urban_tc, &
228  urban_qc, &
229  urban_uc, &
230  urban_trl, &
231  urban_tbl, &
232  urban_tgl, &
233  urban_rainr, &
234  urban_rainb, &
235  urban_raing, &
236  urban_roff
237  use scale_atmos_hydrometeor, only: &
238  hydrometeor_lhv => atmos_hydrometeor_lhv, &
240  i_qv
241  use scale_time, only: &
242  dt => time_dtsec_urban, &
243  nowdate => time_nowdate
244  use scale_atmos_grid_cartesc_real, only: &
245  real_z1 => atmos_grid_cartesc_real_z1
246  use scale_urban_grid_cartesc, only: &
248  use scale_landuse, only: &
250  exists_urban => landuse_exists_urban
251  use mod_urban_admin, only: &
253  use scale_urban_dyn_kusaka01, only: &
255  implicit none
256  logical, intent(in) :: force
257 
258  real(rp) :: trl(uka,uia,uja), tbl(uka,uia,uja), tgl(uka,uia,uja)
259  real(rp) :: tr(uia,uja), tb(uia,uja), tg(uia,uja)
260  real(rp) :: tc(uia,uja), qc(uia,uja), uc(uia,uja)
261  real(rp) :: rainr(uia,uja), rainb(uia,uja), raing(uia,uja)
262 
263  real(rp) :: lhv(uia,uja) ! latent heat of vaporization [J/kg]
264 
265  ! real(RP) :: URBAN_SFLX_LHEX(UIA,UJA)
266 
267  integer :: tloc, tloc_next ! universal time (1-24h)
268  real(rp) :: dsec ! second [s]
269 
270  integer :: k, i, j, iq
271  !---------------------------------------------------------------------------
272 
273  call prof_rapstart('URB_CalcTend', 1)
274 
275  !$acc data create(TRL,TBL,TGL,TR,TB,TG,TC,QC,UC,RAINR,RAINB,RAING,LHV)
276 
277  !########## Get Surface Boundary from coupler ##########
278  call urban_surface_get
279 
280  !########## initialize tendency ##########
281 !OCL XFILL
282  !$acc kernels
283  do j = ujs, uje
284  do i = uis, uie
285  do k = uks, uke
286  urban_trl_t(k,i,j) = 0.0_rp
287  urban_tbl_t(k,i,j) = 0.0_rp
288  urban_tgl_t(k,i,j) = 0.0_rp
289  end do
290  end do
291  end do
292  !$acc end kernels
293 
294 !OCL XFILL
295  !$acc kernels
296  do j = ujs, uje
297  do i = uis, uie
298  urban_tr_t(i,j) = 0.0_rp
299  urban_tb_t(i,j) = 0.0_rp
300  urban_tg_t(i,j) = 0.0_rp
301  urban_tc_t(i,j) = 0.0_rp
302  urban_qc_t(i,j) = 0.0_rp
303  urban_uc_t(i,j) = 0.0_rp
304 
305  urban_rainr_t(i,j) = 0.0_rp
306  urban_rainb_t(i,j) = 0.0_rp
307  urban_raing_t(i,j) = 0.0_rp
308  enddo
309  enddo
310  !$acc end kernels
311 
312 !OCL XFILL
313  !$omp parallel do
314  !$acc kernels
315  do iq = 1, qa
316  do j = ujs, uje
317  do i = uis, uie
318  urban_sflx_qtrc(i,j,iq) = 0.0_rp
319  enddo
320  enddo
321  enddo
322  !$acc end kernels
323 
324  select case ( urban_sfc_type )
325  case ( 'KUSAKA01' )
326 
327 !OCL XFILL
328  !$omp parallel do
329  !$acc kernels
330  do j = ujs, uje
331  do i = uis, uie
332  do k = uks, uke
333  trl(k,i,j) = urban_trl(k,i,j)
334  tbl(k,i,j) = urban_tbl(k,i,j)
335  tgl(k,i,j) = urban_tgl(k,i,j)
336  end do
337  end do
338  end do
339  !$acc end kernels
340 
341 !OCL XFILL
342  !$omp parallel do
343  !$acc kernels
344  do j = ujs, uje
345  do i = uis, uie
346  tr(i,j) = urban_tr(i,j)
347  tb(i,j) = urban_tb(i,j)
348  tg(i,j) = urban_tg(i,j)
349  tc(i,j) = urban_tc(i,j)
350  qc(i,j) = urban_qc(i,j)
351  uc(i,j) = urban_uc(i,j)
352  rainr(i,j) = urban_rainr(i,j)
353  rainb(i,j) = urban_rainb(i,j)
354  raing(i,j) = urban_raing(i,j)
355  end do
356  end do
357  !$acc end kernels
358 
359 
360  ! universal time
361  dsec = real( nowdate(5)*60.0_rp + nowdate(6), kind=rp ) / 3600.0_rp ! [hour]
362  tloc = nowdate(4) ! [hour]
363  tloc = modulo(tloc-1,24)+1
364  if ( tloc == 24 ) then
365  tloc_next = 1
366  else
367  tloc_next = tloc + 1
368  end if
369  !--- Calculate AH at UTC
370  !$acc kernels
371  do j = ujs, uje
372  do i = uis, uie
373  if ( exists_urban(i,j) ) then
374  urban_ah(i,j) = ( 1.0_rp-dsec ) * ah_urb(i,j, tloc) &
375  + ( dsec ) * ah_urb(i,j, tloc_next)
376  urban_ahl(i,j) = ( 1.0_rp-dsec ) * ahl_urb(i,j, tloc) &
377  + ( dsec ) * ahl_urb(i,j, tloc_next)
378  end if
379  enddo
380  enddo
381  !$acc end kernels
382 
383  call hydrometeor_lhv( uia, uis, uie, uja, ujs, uje, &
384  atmos_temp(:,:), lhv(:,:) )
385 
386  call urban_dyn_kusaka01( uka, uks, uke, uia, uis, uie, uja, ujs, uje, &
387  atmos_temp(:,:), atmos_pres(:,:), & ! [IN]
388  atmos_u(:,:), atmos_v(:,:), & ! [IN]
389  atmos_dens(:,:), atmos_qv(:,:), lhv(:,:), & ! [IN]
390  real_z1(:,:), & ! [IN]
391  atmos_sfc_dens(:,:), atmos_sfc_pres(:,:), & ! [IN]
392  atmos_sflx_lw(:,:,:), atmos_sflx_sw(:,:,:), & ! [IN]
393  atmos_sflx_water(:,:), atmos_sflx_engi(:,:), & ! [IN]
394  urban_z0m(:,:), urban_z0h(:,:), urban_z0e(:,:), & ! [IN]
395  urban_zd(:,:), & ! [IN]
396  cdz(:), & ! [IN]
397  tansl_x(:,:), tansl_y(:,:), & ! [IN]
398  landuse_fact_urban(:,:), & ! [IN]
399  dt, & ! [IN]
400  trl(:,:,:), tbl(:,:,:), tgl(:,:,:), & ! [INOUT]
401  tr(:,:), tb(:,:), tg(:,:), tc(:,:), qc(:,:), uc(:,:), & ! [INOUT]
402  rainr(:,:), rainb(:,:), raing(:,:), urban_roff(:,:), & ! [INOUT]
403  urban_sfc_temp(:,:), & ! [OUT]
404  urban_sfc_albedo(:,:,:,:), & ! [OUT]
405  urban_sflx_mw(:,:), urban_sflx_mu(:,:), urban_sflx_mv(:,:), & ! [OUT]
406  urban_sflx_sh(:,:), urban_sflx_lh(:,:), urban_sflx_gh(:,:), & ! [OUT]
407  urban_ustar(:,:), urban_tstar(:,:), urban_qstar(:,:), & ! [OUT]
408  urban_wstar(:,:), & ! [OUT]
409  urban_rlmo(:,:), & ! [OUT]
410  urban_u10(:,:), urban_v10(:,:), urban_t2(:,:), urban_q2(:,:) ) ! [OUT]
411 
412  !-----------------------------------------------------------
413  ! anthropogenic heat fluxes
414  !-----------------------------------------------------------
415  !$omp parallel do
416  !$acc kernels
417  do j = ujs, uje
418  do i = uis, uie
419  if ( exists_urban(i,j) ) then
420  !URBAN_SFLX_SHEX(i,j) = URBAN_AH (i,j) / LANDUSE_fact_urban(i,j) ! Sensible heat flux [W/m2]
421  !URBAN_SFLX_LHEX(i,j) = URBAN_AHL(i,j) / LANDUSE_fact_urban(i,j) ! Latent heat flux [W/m2]
422  !URBAN_SFLX_SH (i,j) = URBAN_SFLX_SH(i,j) + URBAN_SFLX_SHEX(i,j)
423  !URBAN_SFLX_LH (i,j) = URBAN_SFLX_LH(i,j) + URBAN_SFLX_LHEX(i,j)
424  urban_sflx_shex(i,j) = urban_ah(i,j) ! Sensible anthropogenic heat flux [W/m2]
425  urban_sflx_lhex(i,j) = urban_ahl(i,j) ! Latent anthropogenic heat flux [W/m2]
426  end if
427  end do
428  end do
429  !$acc end kernels
430 
431 !OCL XFILL
432  !$omp parallel do
433  !$acc kernels
434  do j = ujs, uje
435  do i = uis, uie
436  if ( exists_urban(i,j) ) then
437  do k = uks, uke
438  urban_trl_t(k,i,j) = ( trl(k,i,j) - urban_trl(k,i,j) ) / dt
439  urban_tbl_t(k,i,j) = ( tbl(k,i,j) - urban_tbl(k,i,j) ) / dt
440  urban_tgl_t(k,i,j) = ( tgl(k,i,j) - urban_tgl(k,i,j) ) / dt
441  end do
442  end if
443  end do
444  end do
445  !$acc end kernels
446 
447 !OCL XFILL
448  !$omp parallel do
449  !$acc kernels
450  do j = ujs, uje
451  do i = uis, uie
452  if ( exists_urban(i,j) ) then
453  urban_tr_t(i,j) = ( tr(i,j) - urban_tr(i,j) ) / dt
454  urban_tb_t(i,j) = ( tb(i,j) - urban_tb(i,j) ) / dt
455  urban_tg_t(i,j) = ( tg(i,j) - urban_tg(i,j) ) / dt
456  urban_tc_t(i,j) = ( tc(i,j) - urban_tc(i,j) ) / dt
457  urban_qc_t(i,j) = ( qc(i,j) - urban_qc(i,j) ) / dt
458  urban_uc_t(i,j) = ( uc(i,j) - urban_uc(i,j) ) / dt
459  urban_rainr_t(i,j) = ( rainr(i,j) - urban_rainr(i,j) ) / dt
460  urban_rainb_t(i,j) = ( rainb(i,j) - urban_rainb(i,j) ) / dt
461  urban_raing_t(i,j) = ( raing(i,j) - urban_raing(i,j) ) / dt
462  end if
463  end do
464  end do
465  !$acc end kernels
466 
467  if ( .NOT. atmos_hydrometeor_dry ) then
468  !$omp parallel do
469  !$acc kernels
470  do j = ujs, uje
471  do i = uis, uie
472  if ( exists_urban(i,j) ) then
473  urban_sflx_qtrc(i,j,i_qv) = urban_sflx_lh(i,j) / lhv(i,j)
474  urban_sflx_qvex(i,j) = urban_sflx_lhex(i,j) / lhv(i,j)
475  end if
476  enddo
477  enddo
478  !$acc end kernels
479  endif
480 
481  end select
482 
483  ! Surface flux for chemical tracers
484  if ( atmos_sw_phy_ch ) then
485  call atmos_phy_ch_driver_urban_flux( urban_sflx_qtrc(:,:,:) ) ! [INOUT]
486  endif
487 
488  call file_history_in( urban_tr_t(:,:), 'URBAN_TR_t', 'tendency of URBAN_TR', 'K/s', dim_type='XY' )
489  call file_history_in( urban_tb_t(:,:), 'URBAN_TB_t', 'tendency of URBAN_TB', 'K/s', dim_type='XY' )
490  call file_history_in( urban_tg_t(:,:), 'URBAN_TG_t', 'tendency of URBAN_TG', 'K/s', dim_type='XY' )
491  call file_history_in( urban_tc_t(:,:), 'URBAN_TC_t', 'tendency of URBAN_TC', 'K/s', dim_type='XY' )
492  call file_history_in( urban_qc_t(:,:), 'URBAN_QC_t', 'tendency of URBAN_QC', 'kg/kg/s', dim_type='XY' )
493  call file_history_in( urban_uc_t(:,:), 'URBAN_UC_t', 'tendency of URBAN_UC', 'm/s2', dim_type='XY' )
494 
495  call file_history_in( urban_trl_t(:,:,:), 'URBAN_TRL_t', 'tendency of URBAN_TRL', 'K/s', dim_type='UXY' )
496  call file_history_in( urban_tbl_t(:,:,:), 'URBAN_TBL_t', 'tendency of URBAN_TBL', 'K/s', dim_type='UXY' )
497  call file_history_in( urban_tgl_t(:,:,:), 'URBAN_TGL_t', 'tendency of URBAN_TGL', 'K/s', dim_type='UXY' )
498 
499  call file_history_in( urban_rainr_t(:,:), 'URBAN_RAINR_t', 'tendency of URBAN_RAINR', 'kg/m2/s', dim_type='XY' )
500  call file_history_in( urban_rainb_t(:,:), 'URBAN_RAINB_t', 'tendency of URBAN_RAINB', 'kg/m2/s', dim_type='XY' )
501  call file_history_in( urban_raing_t(:,:), 'URBAN_RAING_t', 'tendency of URBAN_RAING', 'kg/m2/s', dim_type='XY' )
502  call file_history_in( urban_roff(:,:), 'URBAN_ROFF', 'urban runoff water', 'kg/m2/s', dim_type='XY' )
503 
504  if ( statistics_checktotal ) then
505 
506  call statistics_total( uka, uks, uke, uia, uis, uie, uja, ujs, uje, &
507  urban_trl_t(:,:,:), 'URBAN_TRL_t', &
510  call statistics_total( uka, uks, uke, uia, uis, uie, uja, ujs, uje, &
511  urban_tbl_t(:,:,:), 'URBAN_TBL_t', &
514  call statistics_total( uka, uks, uke, uia, uis, uie, uja, ujs, uje, &
515  urban_tgl_t(:,:,:), 'URBAN_TGL_t', &
518 
519  call statistics_total( uia, uis, uie, uja, ujs, uje, &
520  urban_tr_t(:,:), 'URBAN_TR_t', &
523  call statistics_total( uia, uis, uie, uja, ujs, uje, &
524  urban_tb_t(:,:), 'URBAN_TB_t', &
527  call statistics_total( uia, uis, uie, uja, ujs, uje, &
528  urban_tg_t(:,:), 'URBAN_TG_t', &
531  call statistics_total( uia, uis, uie, uja, ujs, uje, &
532  urban_tc_t(:,:), 'URBAN_TC_t', &
535  call statistics_total( uia, uis, uie, uja, ujs, uje, &
536  urban_qc_t(:,:), 'URBAN_QC_t', &
539  call statistics_total( uia, uis, uie, uja, ujs, uje, &
540  urban_uc_t(:,:), 'URBAN_UC_t', &
543 
544  call statistics_total( uia, uis, uie, uja, ujs, uje, &
545  urban_rainr_t(:,:), 'URBAN_RAINR_t', &
548  call statistics_total( uia, uis, uie, uja, ujs, uje, &
549  urban_rainb_t(:,:), 'URBAN_RAINB_t', &
552  call statistics_total( uia, uis, uie, uja, ujs, uje, &
553  urban_raing_t(:,:), 'URBAN_RAING_t', &
556  call statistics_total( uia, uis, uie, uja, ujs, uje, &
557  urban_roff(:,:), 'URBAN_ROFF', &
560  endif
561 
562 
563  !########## Set Surface Boundary to coupler ##########
564  call urban_surface_set( countup=.true. )
565 
566  !$acc end data
567 
568  call prof_rapend ('URB_CalcTend', 1)
569 
570  return
571  end subroutine urban_driver_calc_tendency
572 
573  !-----------------------------------------------------------------------------
575  subroutine urban_driver_update
576  use scale_time, only: &
577  dt => time_dtsec_urban
578  use mod_urban_vars, only: &
579  urban_trl_t, &
580  urban_tbl_t, &
581  urban_tgl_t, &
582  urban_tr_t, &
583  urban_tb_t, &
584  urban_tg_t, &
585  urban_tc_t, &
586  urban_qc_t, &
587  urban_uc_t, &
588  urban_rainr_t, &
589  urban_rainb_t, &
590  urban_raing_t, &
591  urban_tr, &
592  urban_tb, &
593  urban_tg, &
594  urban_tc, &
595  urban_qc, &
596  urban_uc, &
597  urban_trl, &
598  urban_tbl, &
599  urban_tgl, &
600  urban_rainr, &
601  urban_rainb, &
602  urban_raing, &
604  use mod_urban_admin, only: &
606  use scale_landuse, only: &
607  exists_urban => landuse_exists_urban
608  implicit none
609 
610  integer :: k, i, j
611  !---------------------------------------------------------------------------
612 
613  call prof_rapstart('URB_Update', 1)
614 
615  !########## Get Surface Boundary from coupler ##########
616  call urban_surface_get
617 
618  !########## Dynamics / Update variables ##########
619  select case ( urban_dyn_type )
620  case ( 'KUSAKA01' )
621 
622 !OCL XFILL
623  !$omp parallel do
624  !$acc kernels
625  do j = ujs, uje
626  do i = uis, uie
627  if ( exists_urban(i,j) ) then
628  do k = uks, uke
629  urban_trl(k,i,j) = urban_trl(k,i,j) + urban_trl_t(k,i,j) * dt
630  urban_tbl(k,i,j) = urban_tbl(k,i,j) + urban_tbl_t(k,i,j) * dt
631  urban_tgl(k,i,j) = urban_tgl(k,i,j) + urban_tgl_t(k,i,j) * dt
632  end do
633  end if
634  end do
635  end do
636  !$acc end kernels
637 
638 !OCL XFILL
639  !$omp parallel do
640  !$acc kernels
641  do j = ujs, uje
642  do i = uis, uie
643  if ( exists_urban(i,j) ) then
644  urban_tr(i,j) = urban_tr(i,j) + urban_tr_t(i,j) * dt
645  urban_tb(i,j) = urban_tb(i,j) + urban_tb_t(i,j) * dt
646  urban_tg(i,j) = urban_tg(i,j) + urban_tg_t(i,j) * dt
647  urban_tc(i,j) = urban_tc(i,j) + urban_tc_t(i,j) * dt
648  urban_qc(i,j) = max( real(urban_qc(i,j) + urban_qc_t(i,j) * dt, rp), 0.0_rp )
649  urban_uc(i,j) = max( real(urban_uc(i,j) + urban_uc_t(i,j) * dt, rp), 0.0_rp )
650  urban_rainr(i,j) = max( real(urban_rainr(i,j) + urban_rainr_t(i,j) * dt, rp), 0.0_rp )
651  urban_rainb(i,j) = max( real(urban_rainb(i,j) + urban_rainb_t(i,j) * dt, rp), 0.0_rp )
652  urban_raing(i,j) = max( real(urban_raing(i,j) + urban_raing_t(i,j) * dt, rp), 0.0_rp )
653  end if
654  end do
655  end do
656  !$acc end kernels
657 
658  end select
659 
660  call urban_vars_check
661 
662  call prof_rapend ('URB_Update', 1)
663 
664  return
665  end subroutine urban_driver_update
666 
667  !-----------------------------------------------------------------------------
669  subroutine urban_surface_get
670  use mod_urban_admin, only: &
671  urban_do
672  use mod_urban_vars, only: &
673  atmos_temp, &
674  atmos_pres, &
675  atmos_w, &
676  atmos_u, &
677  atmos_v, &
678  atmos_dens, &
679  atmos_qv, &
680  atmos_pbl, &
681  atmos_sfc_dens, &
682  atmos_sfc_pres, &
683  atmos_sflx_lw, &
684  atmos_sflx_sw, &
685  atmos_cossza, &
688  use mod_cpl_vars, only: &
690  implicit none
691 
692  real(rp) :: atmos_sflx_rad_dn(uia,uja,n_rad_dir,n_rad_rgn)
693 
694  integer :: i, j
695  !---------------------------------------------------------------------------
696 
697  call prof_rapstart('URB_SfcExch', 3)
698 
699  !$acc data create(ATMOS_SFLX_rad_dn)
700 
701  if ( urban_do ) then
702  call cpl_getatm_urb( atmos_temp(:,:), & ! [OUT]
703  atmos_pres(:,:), & ! [OUT]
704  atmos_w(:,:), & ! [OUT]
705  atmos_u(:,:), & ! [OUT]
706  atmos_v(:,:), & ! [OUT]
707  atmos_dens(:,:), & ! [OUT]
708  atmos_qv(:,:), & ! [OUT]
709  atmos_pbl(:,:), & ! [OUT]
710  atmos_sfc_dens(:,:), & ! [OUT]
711  atmos_sfc_pres(:,:), & ! [OUT]
712  atmos_sflx_rad_dn(:,:,:,:), & ! [OUT]
713  atmos_cossza(:,:), & ! [OUT]
714  atmos_sflx_water(:,:), & ! [OUT]
715  atmos_sflx_engi(:,:) ) ! [OUT]
716  endif
717 
718 !OCL XFILL
719  !$omp parallel do
720  !$acc kernels
721  do j = ujs, uje
722  do i = uis, uie
723  atmos_sflx_lw(i,j,i_r_direct ) = atmos_sflx_rad_dn(i,j,i_r_direct ,i_r_ir) ! IR, direct
724  atmos_sflx_lw(i,j,i_r_diffuse) = atmos_sflx_rad_dn(i,j,i_r_diffuse,i_r_ir) ! IR, diffuse
725 
726  atmos_sflx_sw(i,j,i_r_direct ) = atmos_sflx_rad_dn(i,j,i_r_direct ,i_r_nir) & ! NIR, direct
727  + atmos_sflx_rad_dn(i,j,i_r_direct ,i_r_vis) ! VIS, direct
728  atmos_sflx_sw(i,j,i_r_diffuse) = atmos_sflx_rad_dn(i,j,i_r_diffuse,i_r_nir) & ! NIR, diffuse
729  + atmos_sflx_rad_dn(i,j,i_r_diffuse,i_r_vis) ! VIS, diffuse
730  enddo
731  enddo
732  !$acc end kernels
733 
734  !$acc end data
735 
736  call prof_rapend ('URB_SfcExch', 3)
737 
738  return
739  end subroutine urban_surface_get
740 
741  !-----------------------------------------------------------------------------
743  subroutine urban_surface_set( countup )
744  use mod_urban_admin, only: &
745  urban_do
746  use mod_urban_vars, only: &
747  urban_sfc_temp, &
749  urban_sflx_mw, &
750  urban_sflx_mu, &
751  urban_sflx_mv, &
752  urban_sflx_sh, &
753  urban_sflx_lh, &
754  urban_sflx_shex, &
755  urban_sflx_lhex, &
756  urban_sflx_qvex, &
757  urban_sflx_gh, &
758  urban_sflx_qtrc, &
759  urban_z0m, &
760  urban_z0h, &
761  urban_z0e, &
762  urban_u10, &
763  urban_v10, &
764  urban_t2, &
765  urban_q2
766  use mod_cpl_vars, only: &
767  cpl_puturb
768  use scale_landuse, only: &
769  exists_urban => landuse_exists_urban
770  implicit none
771 
772  ! arguments
773  logical, intent(in) :: countup
774  !---------------------------------------------------------------------------
775 
776  call prof_rapstart('URB_SfcExch', 3)
777 
778  if ( urban_do ) then
779  call cpl_puturb( urban_sfc_temp(:,:), & ! [IN]
780  urban_sfc_albedo(:,:,:,:), & ! [IN]
781  urban_z0m(:,:), & ! [IN]
782  urban_z0h(:,:), & ! [IN]
783  urban_z0e(:,:), & ! [IN]
784  urban_sflx_mw(:,:), & ! [IN]
785  urban_sflx_mu(:,:), & ! [IN]
786  urban_sflx_mv(:,:), & ! [IN]
787  urban_sflx_sh(:,:), & ! [IN]
788  urban_sflx_lh(:,:), & ! [IN]
789  urban_sflx_shex(:,:), & ! [IN]
790  urban_sflx_lhex(:,:), & ! [IN]
791  urban_sflx_qvex(:,:), & ! [IN]
792  urban_sflx_gh(:,:), & ! [IN]
793  urban_sflx_qtrc(:,:,:), & ! [IN]
794  urban_u10(:,:), & ! [IN]
795  urban_v10(:,:), & ! [IN]
796  urban_t2(:,:), & ! [IN]
797  urban_q2(:,:), & ! [IN]
798  exists_urban(:,:), & ! [IN]
799  countup ) ! [IN]
800  endif
801 
802  call prof_rapend ('URB_SfcExch', 3)
803 
804  return
805  end subroutine urban_surface_set
806 
807 end module mod_urban_driver
mod_urban_vars::urban_q2
real(rp), dimension(:,:), allocatable, public urban_q2
Definition: mod_urban_vars.F90:122
mod_urban_vars::urban_trl_t
real(rp), dimension(:,:,:), allocatable, public urban_trl_t
Definition: mod_urban_vars.F90:79
mod_urban_vars::atmos_sflx_sw
real(rp), dimension(:,:,:), allocatable, public atmos_sflx_sw
Definition: mod_urban_vars.F90:136
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:106
mod_urban_vars::atmos_u
real(rp), dimension(:,:), allocatable, public atmos_u
Definition: mod_urban_vars.F90:128
scale_urban_grid_cartesc_index::uja
integer, public uja
Definition: scale_urban_grid_cartesC_index.F90:45
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:86
mod_urban_vars::urban_sfc_temp
real(rp), dimension(:,:), allocatable, public urban_sfc_temp
Definition: mod_urban_vars.F90:75
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:39
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_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:40
scale_tracer::qa
integer, public qa
Definition: scale_tracer.F90:35
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:38
mod_urban_vars::urban_tb
real(rp), dimension(:,:), allocatable, public urban_tb
Definition: mod_urban_vars.F90:65
scale_urban_grid_cartesc_index::uia
integer, public uia
Definition: scale_urban_grid_cartesC_index.F90:41
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:60
mod_urban_vars::urban_raing_t
real(rp), dimension(:,:), allocatable, public urban_raing_t
Definition: mod_urban_vars.F90:90
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:41
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:76
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:139
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:174
mod_urban_vars::urban_trl
real(rp), dimension(:,:,:), allocatable, public urban_trl
Definition: mod_urban_vars.F90:61
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:99
mod_urban_vars::urban_sflx_mu
real(rp), dimension(:,:), allocatable, public urban_sflx_mu
Definition: mod_urban_vars.F90:95
mod_urban_vars::urban_raing
real(rp), dimension(:,:), allocatable, public urban_raing
Definition: mod_urban_vars.F90:72
mod_urban_vars::urban_z0h
real(rp), dimension(:,:), allocatable, public urban_z0h
Definition: mod_urban_vars.F90:107
scale_topography
module TOPOGRAPHY
Definition: scale_topography.F90:11
mod_urban_vars::urban_qc
real(rp), dimension(:,:), allocatable, public urban_qc
Definition: mod_urban_vars.F90:68
scale_atmos_hydrometeor::atmos_hydrometeor_dry
logical, public atmos_hydrometeor_dry
Definition: scale_atmos_hydrometeor.F90:114
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:132
mod_urban_vars::urban_tbl_t
real(rp), dimension(:,:,:), allocatable, public urban_tbl_t
Definition: mod_urban_vars.F90:80
mod_urban_vars::urban_tgl_t
real(rp), dimension(:,:,:), allocatable, public urban_tgl_t
Definition: mod_urban_vars.F90:81
scale_atmos_grid_cartesc_real
module Atmosphere GRID CartesC Real(real space)
Definition: scale_atmos_grid_cartesC_real.F90:11
mod_urban_vars::urban_sflx_lhex
real(rp), dimension(:,:), allocatable, public urban_sflx_lhex
Definition: mod_urban_vars.F90:100
mod_urban_driver::urban_driver_update
subroutine, public urban_driver_update
Urban step.
Definition: mod_urban_driver.F90:576
mod_urban_vars::atmos_w
real(rp), dimension(:,:), allocatable, public atmos_w
Definition: mod_urban_vars.F90:127
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:129
mod_urban_driver::urban_driver_calc_tendency
subroutine, public urban_driver_calc_tendency(force)
Calclate tendency.
Definition: mod_urban_driver.F90:152
mod_urban_vars::atmos_sflx_water
real(rp), dimension(:,:), allocatable, public atmos_sflx_water
Definition: mod_urban_vars.F90:138
mod_urban_vars::atmos_dens
real(rp), dimension(:,:), allocatable, public atmos_dens
Definition: mod_urban_vars.F90:130
mod_urban_vars::urban_qstar
real(rp), dimension(:,:), allocatable, public urban_qstar
Definition: mod_urban_vars.F90:116
scale_prc
module PROCESS
Definition: scale_prc.F90:11
mod_urban_vars::urban_sflx_qtrc
real(rp), dimension(:,:,:), allocatable, public urban_sflx_qtrc
Definition: mod_urban_vars.F90:102
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:117
mod_urban_vars::urban_t2
real(rp), dimension(:,:), allocatable, public urban_t2
Definition: mod_urban_vars.F90:121
mod_urban_admin::urban_dyn_type
character(len=h_short), public urban_dyn_type
Definition: mod_urban_admin.F90:35
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:45
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:96
mod_urban_vars::urban_tgl
real(rp), dimension(:,:,:), allocatable, public urban_tgl
Definition: mod_urban_vars.F90:63
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:120
mod_urban_vars::atmos_qv
real(rp), dimension(:,:), allocatable, public atmos_qv
Definition: mod_urban_vars.F90:131
mod_urban_vars::atmos_sfc_pres
real(rp), dimension(:,:), allocatable, public atmos_sfc_pres
Definition: mod_urban_vars.F90:134
scale_topography::topography_tansl_x
real(rp), dimension(:,:), allocatable, public topography_tansl_x
tan(slope_x)
Definition: scale_topography.F90:40
mod_urban_vars::atmos_sflx_lw
real(rp), dimension(:,:,:), allocatable, public atmos_sflx_lw
Definition: mod_urban_vars.F90:135
mod_urban_vars::urban_z0e
real(rp), dimension(:,:), allocatable, public urban_z0e
Definition: mod_urban_vars.F90:108
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:1450
mod_urban_vars::urban_u10
real(rp), dimension(:,:), allocatable, public urban_u10
Definition: mod_urban_vars.F90:119
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:137
mod_urban_driver::urban_surface_get
subroutine, public urban_surface_get
Get surface boundary.
Definition: mod_urban_driver.F90:670
scale_prof
module profiler
Definition: scale_prof.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, AH_TOFFSET)
Setup.
Definition: scale_urban_dyn_kusaka01.F90:139
mod_urban_driver::urban_surface_set
subroutine, public urban_surface_set(countup)
Set surface boundary to other model.
Definition: mod_urban_driver.F90:744
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:328
mod_atmos_admin::atmos_sw_phy_ch
logical, public atmos_sw_phy_ch
Definition: mod_atmos_admin.F90:54
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:1059
mod_urban_vars::atmos_pres
real(rp), dimension(:,:), allocatable, public atmos_pres
Definition: mod_urban_vars.F90:126
mod_urban_vars::urban_tstar
real(rp), dimension(:,:), allocatable, public urban_tstar
Definition: mod_urban_vars.F90:115
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:67
mod_urban_vars::urban_tg
real(rp), dimension(:,:), allocatable, public urban_tg
Definition: mod_urban_vars.F90:66
scale_atmos_hydrometeor::i_qv
integer, public i_qv
Definition: scale_atmos_hydrometeor.F90:93
scale_urban_grid_cartesc_index::uks
integer, public uks
Definition: scale_urban_grid_cartesC_index.F90:38
mod_urban_vars::urban_tg_t
real(rp), dimension(:,:), allocatable, public urban_tg_t
Definition: mod_urban_vars.F90:87
mod_urban_vars::urban_uc
real(rp), dimension(:,:), allocatable, public urban_uc
Definition: mod_urban_vars.F90:69
mod_urban_vars::urban_tr_t
real(rp), dimension(:,:), allocatable, public urban_tr_t
Definition: mod_urban_vars.F90:85
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:56
scale_statistics::statistics_checktotal
logical, public statistics_checktotal
calc&report variable totals to logfile?
Definition: scale_statistics.F90:109
mod_urban_vars::urban_sflx_lh
real(rp), dimension(:,:), allocatable, public urban_sflx_lh
Definition: mod_urban_vars.F90:98
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:52
scale_urban_dyn_kusaka01::urban_dyn_kusaka01_finalize
subroutine, public urban_dyn_kusaka01_finalize
Finalize.
Definition: scale_urban_dyn_kusaka01.F90:399
scale_cpl_sfc_index
module coupler / surface-atmospehre
Definition: scale_cpl_sfc_index.F90:11
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_LHEX, SFLX_QVEX, SFLX_GH, SFLX_QTRC, U10, V10, T2, Q2, mask, countup)
Definition: mod_cpl_vars.F90:1025
scale_time::time_nowdate
integer, dimension(6), public time_nowdate
current time [YYYY MM DD HH MM SS]
Definition: scale_time.F90:68
mod_urban_vars::urban_rainb
real(rp), dimension(:,:), allocatable, public urban_rainb
Definition: mod_urban_vars.F90:71
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:110
mod_urban_vars::urban_ahl
real(rp), dimension(:,:), allocatable, public urban_ahl
Definition: mod_urban_vars.F90:111
mod_urban_vars::urban_sflx_sh
real(rp), dimension(:,:), allocatable, public urban_sflx_sh
Definition: mod_urban_vars.F90:97
mod_urban_vars::urban_qc_t
real(rp), dimension(:,:), allocatable, public urban_qc_t
Definition: mod_urban_vars.F90:84
mod_urban_vars::urban_tr
real(rp), dimension(:,:), allocatable, public urban_tr
Definition: mod_urban_vars.F90:64
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, 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:456
scale_urban_grid_cartesc_index::uka
integer, public uka
Definition: scale_urban_grid_cartesC_index.F90:37
mod_urban_vars::urban_roff
real(rp), dimension(:,:), allocatable, public urban_roff
Definition: mod_urban_vars.F90:92
mod_urban_vars::urban_rainb_t
real(rp), dimension(:,:), allocatable, public urban_rainb_t
Definition: mod_urban_vars.F90:89
scale_urban_grid_cartesc_index::uis
integer, public uis
Definition: scale_urban_grid_cartesC_index.F90:42
mod_urban_vars::urban_ustar
real(rp), dimension(:,:), allocatable, public urban_ustar
Definition: mod_urban_vars.F90:114
mod_urban_vars::urban_zd
real(rp), dimension(:,:), allocatable, public urban_zd
Definition: mod_urban_vars.F90:109
scale_prof::prof_rapend
subroutine, public prof_rapend(rapname_base, level, disable_barrier)
Save raptime.
Definition: scale_prof.F90:246
scale_const::const_d2r
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:33
mod_urban_vars::urban_sflx_qvex
real(rp), dimension(:,:), allocatable, public urban_sflx_qvex
Definition: mod_urban_vars.F90:101
scale_landuse::landuse_fact_urban
real(rp), dimension(:,:), allocatable, public landuse_fact_urban
urban factor
Definition: scale_landuse.F90:47
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:38
mod_urban_vars::urban_uc_t
real(rp), dimension(:,:), allocatable, public urban_uc_t
Definition: mod_urban_vars.F90:83
mod_urban_vars::urban_rainr_t
real(rp), dimension(:,:), allocatable, public urban_rainr_t
Definition: mod_urban_vars.F90:88
mod_urban_vars::urban_tbl
real(rp), dimension(:,:,:), allocatable, public urban_tbl
Definition: mod_urban_vars.F90:62
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:94
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:133
scale_urban_grid_cartesc_real
module urban / grid / cartesianC / real
Definition: scale_urban_grid_cartesC_real.F90:12
mod_urban_driver::urban_driver_finalize
subroutine, public urban_driver_finalize
Finalize.
Definition: mod_urban_driver.F90:117
mod_urban_vars::urban_sflx_gh
real(rp), dimension(:,:), allocatable, public urban_sflx_gh
Definition: mod_urban_vars.F90:103
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_urban_grid_cartesc_index::uje
integer, public uje
Definition: scale_urban_grid_cartesC_index.F90:47
mod_urban_vars::urban_rlmo
real(rp), dimension(:,:), allocatable, public urban_rlmo
Definition: mod_urban_vars.F90:118
scale_urban_grid_cartesc_index::ujs
integer, public ujs
Definition: scale_urban_grid_cartesC_index.F90:46
mod_urban_vars::urban_rainr
real(rp), dimension(:,:), allocatable, public urban_rainr
Definition: mod_urban_vars.F90:70
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:43
scale_urban_grid_cartesc_index::uke
integer, public uke
Definition: scale_urban_grid_cartesC_index.F90:39
mod_urban_vars::atmos_temp
real(rp), dimension(:,:), allocatable, public atmos_temp
Definition: mod_urban_vars.F90:125
mod_urban_vars::urban_tc_t
real(rp), dimension(:,:), allocatable, public urban_tc_t
Definition: mod_urban_vars.F90:82
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:37