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_resume
 Resume. More...
 
subroutine, public urban_driver
 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

Function/Subroutine Documentation

◆ urban_driver_setup()

subroutine, public mod_urban_driver::urban_driver_setup ( )

Setup.

Definition at line 55 of file mod_urban_driver.f90.

References scale_stdio::io_fid_log, scale_stdio::io_l, and mod_urban_phy_driver::urban_phy_driver_setup().

Referenced by mod_rm_driver::scalerm().

55  use mod_urban_phy_driver, only: &
57  implicit none
58  !---------------------------------------------------------------------------
59 
60  if( io_l ) write(io_fid_log,*)
61  if( io_l ) write(io_fid_log,*) '++++++ Module[DRIVER] / Categ[URBAN] / Origin[SCALE-RM]'
62 
64 
65  return
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
module URBAN / Physics Urban Canopy Model (UCM)
subroutine, public urban_phy_driver_setup
Setup.
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
Here is the call graph for this function:
Here is the caller graph for this function:

◆ urban_driver_resume()

subroutine, public mod_urban_driver::urban_driver_resume ( )

Resume.

Definition at line 71 of file mod_urban_driver.f90.

References scale_stdio::io_fid_log, scale_stdio::io_l, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), mod_urban_phy_driver::urban_phy_driver_resume(), urban_surface_get(), urban_surface_set(), mod_urban_admin::urban_sw, and mod_urban_vars::urban_vars_history().

Referenced by mod_rm_driver::resume_state().

71  use mod_urban_phy_driver, only: &
73  use mod_urban_vars, only: &
75  use mod_urban_admin, only: &
76  urban_sw
77  implicit none
78  !---------------------------------------------------------------------------
79 
80  if( io_l ) write(io_fid_log,*)
81  if( io_l ) write(io_fid_log,*) '++++++ Module[DRIVER] / Categ[URBAN] / Origin[SCALE-RM]'
82 
83  !########## Get Surface Boundary from coupler ##########
84  call urban_surface_get
85 
87 
88  !########## Set Surface Boundary to coupler ##########
89  call urban_surface_set( countup=.true. )
90 
91  !########## History & Monitor ##########
92  if ( urban_sw ) then
93  call prof_rapstart('URB_History', 1)
95  call prof_rapend ('URB_History', 1)
96  endif
97 
98  return
logical, public urban_sw
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
subroutine, public urban_vars_history
History output set for urban variables.
module URBAN / Physics Urban Canopy Model (UCM)
module URBAN Variables
subroutine, public urban_surface_set(countup)
Set surface boundary to other model.
subroutine, public urban_phy_driver_resume
Resume.
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:132
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:178
module Urban admin
Here is the call graph for this function:
Here is the caller graph for this function:

◆ urban_driver()

subroutine, public mod_urban_driver::urban_driver ( )

Urban step.

Definition at line 104 of file mod_urban_driver.f90.

References scale_grid_index::ie, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, scale_prof::prof_rapend(), scale_prof::prof_rapstart(), scale_time::time_dtsec_urban, scale_urban_grid_index::uke, scale_urban_grid_index::uks, mod_urban_phy_driver::urban_phy_driver(), 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, urban_surface_get(), urban_surface_set(), mod_urban_admin::urban_sw, 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_uc, mod_urban_vars::urban_uc_t, mod_urban_vars::urban_vars_history(), and mod_urban_vars::urban_vars_total().

Referenced by mod_rm_driver::scalerm().

104  use scale_time, only: &
105  dt => time_dtsec_urban
106  use mod_urban_admin, only: &
107  urban_sw
108  use mod_urban_vars, only: &
109  urban_tr_t, &
110  urban_tb_t, &
111  urban_tg_t, &
112  urban_tc_t, &
113  urban_qc_t, &
114  urban_uc_t, &
115  urban_trl_t, &
116  urban_tbl_t, &
117  urban_tgl_t, &
118  urban_rainr_t, &
119  urban_rainb_t, &
120  urban_raing_t, &
121  urban_roff_t, &
122  urban_tr, &
123  urban_tb, &
124  urban_tg, &
125  urban_tc, &
126  urban_qc, &
127  urban_uc, &
128  urban_trl, &
129  urban_tbl, &
130  urban_tgl, &
131  urban_rainr, &
132  urban_rainb, &
133  urban_raing, &
134  urban_roff, &
137  use mod_urban_phy_driver, only: &
139  implicit none
140 
141  integer :: k, i, j
142  !---------------------------------------------------------------------------
143 
144  !########## Get Surface Boundary from coupler ##########
145  call prof_rapstart('URB_SfcExch', 2)
146  call urban_surface_get
147  call prof_rapend ('URB_SfcExch', 2)
148 
149  !########## Physics ##########
150  if ( urban_sw ) then
151  call prof_rapstart('URB_Physics', 1)
152  call urban_phy_driver( update_flag = .true. )
153  call prof_rapend ('URB_Physics', 1)
154  endif
155 
156  !########## Update ##########
157  do j = js, je
158  do i = is, ie
159  do k = uks, uke
160  urban_trl(k,i,j) = urban_trl(k,i,j) + urban_trl_t(k,i,j) * dt
161  urban_tbl(k,i,j) = urban_tbl(k,i,j) + urban_tbl_t(k,i,j) * dt
162  urban_tgl(k,i,j) = urban_tgl(k,i,j) + urban_tgl_t(k,i,j) * dt
163  end do
164  end do
165  end do
166  do j = js, je
167  do i = is, ie
168  urban_tr(i,j) = urban_tr(i,j) + urban_tr_t(i,j) * dt
169  urban_tb(i,j) = urban_tb(i,j) + urban_tb_t(i,j) * dt
170  urban_tg(i,j) = urban_tg(i,j) + urban_tg_t(i,j) * dt
171  urban_tc(i,j) = urban_tc(i,j) + urban_tc_t(i,j) * dt
172  urban_qc(i,j) = urban_qc(i,j) + urban_qc_t(i,j) * dt
173  urban_uc(i,j) = urban_uc(i,j) + urban_uc_t(i,j) * dt
174 
175  urban_rainr(i,j) = urban_rainr(i,j) + urban_rainr_t(i,j) * dt
176  urban_rainb(i,j) = urban_rainb(i,j) + urban_rainb_t(i,j) * dt
177  urban_raing(i,j) = urban_raing(i,j) + urban_raing_t(i,j) * dt
178  urban_roff(i,j) = urban_roff(i,j) + urban_roff_t(i,j) * dt
179  end do
180  end do
181 
182  call urban_vars_total
183 
184  !########## Set Surface Boundary to coupler ##########
185  call prof_rapstart('URB_SfcExch', 2)
186  call urban_surface_set( countup=.true. )
187  call prof_rapend ('URB_SfcExch', 2)
188 
189  !########## reset tendencies ##########
190 !OCL XFILL
191  do j = js, je
192  do i = is, ie
193  do k = uks, uke
194  urban_trl_t(k,i,j) = 0.0_rp
195  urban_tbl_t(k,i,j) = 0.0_rp
196  urban_tgl_t(k,i,j) = 0.0_rp
197  end do
198  end do
199  end do
200 
201 !OCL XFILL
202  do j = js, je
203  do i = is, ie
204  urban_tr_t(i,j) = 0.0_rp
205  urban_tb_t(i,j) = 0.0_rp
206  urban_tg_t(i,j) = 0.0_rp
207  urban_tc_t(i,j) = 0.0_rp
208  urban_qc_t(i,j) = 0.0_rp
209  urban_uc_t(i,j) = 0.0_rp
210 
211  urban_rainr_t(i,j) = 0.0_rp
212  urban_rainb_t(i,j) = 0.0_rp
213  urban_raing_t(i,j) = 0.0_rp
214  urban_roff_t(i,j) = 0.0_rp
215  enddo
216  enddo
217 
218  !########## History & Monitor ##########
219  call prof_rapstart('URB_History', 1)
220  call urban_vars_history
221  call prof_rapend ('URB_History', 1)
222 
223  return
real(rp), dimension(:,:), allocatable, public urban_qc_t
logical, public urban_sw
real(rp), dimension(:,:), allocatable, public urban_qc
real(rp), dimension(:,:), allocatable, public urban_rainr_t
real(rp), dimension(:,:,:), allocatable, public urban_trl_t
real(rp), dimension(:,:), allocatable, public urban_tg_t
real(rp), dimension(:,:), allocatable, public urban_raing_t
subroutine, public urban_vars_total
Budget monitor for urban.
real(rp), dimension(:,:), allocatable, public urban_tb_t
subroutine, public urban_vars_history
History output set for urban variables.
real(rp), dimension(:,:), allocatable, public urban_tb
module URBAN / Physics Urban Canopy Model (UCM)
module URBAN Variables
real(rp), dimension(:,:), allocatable, public urban_raing
real(rp), dimension(:,:), allocatable, public urban_uc
subroutine, public urban_surface_set(countup)
Set surface boundary to other model.
real(rp), dimension(:,:), allocatable, public urban_tr
real(rp), dimension(:,:,:), allocatable, public urban_tgl
real(rp), dimension(:,:), allocatable, public urban_uc_t
subroutine, public urban_phy_driver(update_flag)
Driver.
module TIME
Definition: scale_time.F90:15
real(rp), dimension(:,:), allocatable, public urban_roff
real(rp), dimension(:,:), allocatable, public urban_tc
real(rp), dimension(:,:), allocatable, public urban_rainr
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:132
real(rp), dimension(:,:), allocatable, public urban_tg
real(rp), dimension(:,:,:), allocatable, public urban_trl
real(rp), dimension(:,:), allocatable, public urban_rainb_t
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:178
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:49
real(rp), dimension(:,:), allocatable, public urban_tr_t
real(rp), dimension(:,:,:), allocatable, public urban_tgl_t
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
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 229 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_pres, mod_urban_vars::atmos_sflx_lw, mod_urban_vars::atmos_sflx_prec, 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(), and mod_urban_admin::urban_sw.

Referenced by urban_driver(), and urban_driver_resume().

229  use mod_urban_admin, only: &
230  urban_sw
231  use mod_urban_vars, only: &
232  atmos_temp, &
233  atmos_pres, &
234  atmos_w, &
235  atmos_u, &
236  atmos_v, &
237  atmos_dens, &
238  atmos_qv, &
239  atmos_pbl, &
240  atmos_sfc_pres, &
241  atmos_sflx_lw, &
242  atmos_sflx_sw, &
243  atmos_cossza, &
245  use mod_cpl_vars, only: &
247  implicit none
248 
249  real(RP) :: atmos_sflx_rad_dn(ia,ja,2,2)
250  real(RP) :: atmos_sflx_rain (ia,ja)
251  real(RP) :: atmos_sflx_snow (ia,ja)
252  !---------------------------------------------------------------------------
253 
254  if ( urban_sw ) then
255  call cpl_getatm_urb( atmos_temp(:,:), & ! [OUT]
256  atmos_pres(:,:), & ! [OUT]
257  atmos_w(:,:), & ! [OUT]
258  atmos_u(:,:), & ! [OUT]
259  atmos_v(:,:), & ! [OUT]
260  atmos_dens(:,:), & ! [OUT]
261  atmos_qv(:,:), & ! [OUT]
262  atmos_pbl(:,:), & ! [OUT]
263  atmos_sfc_pres(:,:), & ! [OUT]
264  atmos_sflx_rad_dn(:,:,:,:), & ! [OUT]
265  atmos_cossza(:,:), & ! [OUT]
266  atmos_sflx_rain(:,:), & ! [OUT]
267  atmos_sflx_snow(:,:) ) ! [OUT]
268  endif
269 
270  atmos_sflx_sw(:,:,:) = atmos_sflx_rad_dn(:,:,i_sw,:) ! direct/diffuse
271  atmos_sflx_lw(:,:,:) = atmos_sflx_rad_dn(:,:,i_lw,:) ! direct/diffuse
272 
273  atmos_sflx_prec(:,:) = atmos_sflx_rain(:,:) + atmos_sflx_snow(:,:) ! liquid+ice
274 
275  return
logical, public urban_sw
module URBAN Variables
subroutine, public cpl_getatm_urb(TEMP, PRES, W, U, V, DENS, QV, PBL, SFC_PRES, SFLX_rad_dn, cosSZA, SFLX_rain, SFLX_snow)
real(rp), dimension(:,:), allocatable, public atmos_sflx_prec
integer, parameter, public i_lw
integer, parameter, public i_sw
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_pres
real(rp), dimension(:,:), allocatable, public atmos_temp
real(rp), dimension(:,:), allocatable, public atmos_u
real(rp), dimension(:,:,:), allocatable, public atmos_sflx_lw
module Urban admin
real(rp), dimension(:,:), allocatable, public atmos_w
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 281 of file mod_urban_driver.f90.

References mod_cpl_vars::cpl_puturb(), mod_urban_vars::urban_q2, mod_urban_vars::urban_sfc_albedo, mod_urban_vars::urban_sfc_temp, mod_urban_vars::urban_sflx_evap, 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_sh, mod_urban_admin::urban_sw, 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_mkinit::mkinit(), mod_rm_driver::resume_state(), urban_driver(), and urban_driver_resume().

281  use mod_urban_admin, only: &
282  urban_sw
283  use mod_urban_vars, only: &
284  urban_sfc_temp, &
286  urban_sflx_mw, &
287  urban_sflx_mu, &
288  urban_sflx_mv, &
289  urban_sflx_sh, &
290  urban_sflx_lh, &
291  urban_sflx_gh, &
292  urban_sflx_evap, &
293  urban_z0m, &
294  urban_z0h, &
295  urban_z0e, &
296  urban_u10, &
297  urban_v10, &
298  urban_t2, &
299  urban_q2
300  use mod_cpl_vars, only: &
301  cpl_puturb
302  implicit none
303 
304  ! arguments
305  logical, intent(in) :: countup
306  !---------------------------------------------------------------------------
307 
308  if ( urban_sw ) then
309  call cpl_puturb( urban_sfc_temp(:,:), & ! [IN]
310  urban_sfc_albedo(:,:,:), & ! [IN]
311  urban_z0m(:,:), & ! [IN]
312  urban_z0h(:,:), & ! [IN]
313  urban_z0e(:,:), & ! [IN]
314  urban_sflx_mw(:,:), & ! [IN]
315  urban_sflx_mu(:,:), & ! [IN]
316  urban_sflx_mv(:,:), & ! [IN]
317  urban_sflx_sh(:,:), & ! [IN]
318  urban_sflx_lh(:,:), & ! [IN]
319  urban_sflx_gh(:,:), & ! [IN]
320  urban_sflx_evap(:,:), & ! [IN]
321  urban_u10(:,:), & ! [IN]
322  urban_v10(:,:), & ! [IN]
323  urban_t2(:,:), & ! [IN]
324  urban_q2(:,:), & ! [IN]
325  countup ) ! [IN]
326  endif
327 
328  return
logical, public urban_sw
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_sflx_evap
real(rp), dimension(:,:), allocatable, public urban_z0m
real(rp), dimension(:,:), allocatable, public urban_t2
module URBAN Variables
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_GH, SFLX_evap, U10, V10, T2, Q2, countup)
real(rp), dimension(:,:), allocatable, public urban_sflx_sh
module COUPLER Variables
real(rp), dimension(:,:), allocatable, public urban_z0h
real(rp), dimension(:,:,:), allocatable, public urban_sfc_albedo
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
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: