SCALE-RM
mod_urban_phy_driver.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
9 !-------------------------------------------------------------------------------
11  !-----------------------------------------------------------------------------
12  !
13  !++ used modules
14  !
15  use scale_precision
16  use scale_stdio
17  use scale_prof
20 
21  use scale_const, only: &
22  i_sw => const_i_sw, &
23  i_lw => const_i_lw
24  !-----------------------------------------------------------------------------
25  implicit none
26  private
27  !-----------------------------------------------------------------------------
28  !
29  !++ Public procedure
30  !
31  public :: urban_phy_driver_setup
32  public :: urban_phy_driver_resume
33  public :: urban_phy_driver
34 
35  !-----------------------------------------------------------------------------
36  !
37  !++ Public parameters & variables
38  !
39  !-----------------------------------------------------------------------------
40  !
41  !++ Private procedure
42  !
43  !-----------------------------------------------------------------------------
44  !
45  !++ Private parameters & variables
46  !
47  !-----------------------------------------------------------------------------
48 contains
49  !-----------------------------------------------------------------------------
51  subroutine urban_phy_driver_setup
52  use scale_urban_phy, only: &
54  use mod_urban_admin, only: &
55  urban_type, &
56  urban_sw
57  use mod_urban_vars, only: &
58  urban_z0m, &
59  urban_z0h, &
60  urban_z0e
61  implicit none
62  !---------------------------------------------------------------------------
63 
64  if( io_l ) write(io_fid_log,*)
65  if( io_l ) write(io_fid_log,*) '++++++ Module[DRIVER] / Categ[URBAN PHY] / Origin[SCALE-RM]'
66 
67  if ( urban_sw ) then
68 
69  ! setup library component
70  call urban_phy_setup( urban_type, & ! [IN]
71  urban_z0m(:,:), & ! [OUT]
72  urban_z0h(:,:), & ! [OUT]
73  urban_z0e(:,:) ) ! [OUT]
74 
75  else
76  if( io_l ) write(io_fid_log,*) '*** this component is never called.'
77  endif
78 
79  return
80  end subroutine urban_phy_driver_setup
81 
82  !-----------------------------------------------------------------------------
84  subroutine urban_phy_driver_resume
85  use mod_admin_restart, only: &
87  use mod_urban_admin, only: &
88  urban_sw
89  implicit none
90 
91  if ( urban_sw ) then
92 
93  if ( .NOT. restart_run ) then ! tentative
94  ! run once (only for the diagnostic value)
95  call prof_rapstart('URB_Physics', 1)
96  call urban_phy_driver( update_flag = .true. )
97  call prof_rapend ('URB_Physics', 1)
98  end if
99 
100  end if
101 
102  return
103  end subroutine urban_phy_driver_resume
104  !-----------------------------------------------------------------------------
106  subroutine urban_phy_driver( update_flag )
108  hydrometeor_lhv => atmos_hydrometeor_lhv
109  use scale_time, only: &
110  nowdate => time_nowdate, &
111  dt => time_dtsec_urban
112  use scale_rm_statistics, only: &
114  stat_total
115  use scale_history, only: &
116  hist_in
117  use scale_mapproj, only: &
118  base_lon => mprj_basepoint_lon, &
119  base_lat => mprj_basepoint_lat
120  use scale_grid_real, only: &
121  real_z1
122  use scale_urban_phy, only: &
123  urban_phy
124  use mod_urban_vars, only: &
125  urban_tr, &
126  urban_tb, &
127  urban_tg, &
128  urban_tc, &
129  urban_qc, &
130  urban_uc, &
131  urban_trl, &
132  urban_tbl, &
133  urban_tgl, &
134  urban_rainr, &
135  urban_rainb, &
136  urban_raing, &
137  urban_roff, &
138  urban_tr_t, &
139  urban_tb_t, &
140  urban_tg_t, &
141  urban_tc_t, &
142  urban_qc_t, &
143  urban_uc_t, &
144  urban_trl_t, &
145  urban_tbl_t, &
146  urban_tgl_t, &
147  urban_rainr_t, &
148  urban_rainb_t, &
149  urban_raing_t, &
150  urban_roff_t, &
151  urban_sfc_temp, &
153  urban_sflx_mw, &
154  urban_sflx_mu, &
155  urban_sflx_mv, &
156  urban_sflx_sh, &
157  urban_sflx_lh, &
158  urban_sflx_gh, &
159  urban_sflx_evap, &
160  urban_z0m, &
161  urban_z0h, &
162  urban_z0e, &
163  urban_u10, &
164  urban_v10, &
165  urban_t2, &
166  urban_q2, &
167  atmos_temp, &
168  atmos_pres, &
169  atmos_w, &
170  atmos_u, &
171  atmos_v, &
172  atmos_dens, &
173  atmos_qv, &
174  atmos_pbl, &
175  atmos_sfc_pres, &
176  atmos_sflx_lw, &
177  atmos_sflx_sw, &
179  implicit none
180 
181  ! arguments
182  logical, intent(in) :: update_flag
183 
184  ! works
185  real(RP) :: lhv(ia,ja) ! latent heat of vaporization [J/kg]
186  real(RP) :: total ! dummy
187 
188  character(len=2) :: sk
189 
190  integer :: k, i, j
191  !---------------------------------------------------------------------------
192 
193  if ( update_flag ) then
194 
195  call urban_phy( urban_tr_t(:,:), & ! [OUT]
196  urban_tb_t(:,:), & ! [OUT]
197  urban_tg_t(:,:), & ! [OUT]
198  urban_tc_t(:,:), & ! [OUT]
199  urban_qc_t(:,:), & ! [OUT]
200  urban_uc_t(:,:), & ! [OUT]
201  urban_trl_t(:,:,:), & ! [OUT]
202  urban_tbl_t(:,:,:), & ! [OUT]
203  urban_tgl_t(:,:,:), & ! [OUT]
204  urban_rainr_t(:,:), & ! [OUT]
205  urban_rainb_t(:,:), & ! [OUT]
206  urban_raing_t(:,:), & ! [OUT]
207  urban_roff_t(:,:), & ! [OUT]
208  urban_sfc_temp(:,:), & ! [OUT]
209  urban_sfc_albedo(:,:,i_lw), & ! [OUT]
210  urban_sfc_albedo(:,:,i_sw), & ! [OUT]
211  urban_sflx_mw(:,:), & ! [OUT]
212  urban_sflx_mu(:,:), & ! [OUT]
213  urban_sflx_mv(:,:), & ! [OUT]
214  urban_sflx_sh(:,:), & ! [OUT]
215  urban_sflx_lh(:,:), & ! [OUT]
216  urban_sflx_gh(:,:), & ! [OUT]
217  urban_z0m(:,:), & ! [OUT]
218  urban_z0h(:,:), & ! [OUT]
219  urban_z0e(:,:), & ! [OUT]
220  urban_u10(:,:), & ! [OUT]
221  urban_v10(:,:), & ! [OUT]
222  urban_t2(:,:), & ! [OUT]
223  urban_q2(:,:), & ! [OUT]
224  atmos_temp(:,:), & ! [IN]
225  atmos_pres(:,:), & ! [IN]
226  atmos_w(:,:), & ! [IN]
227  atmos_u(:,:), & ! [IN]
228  atmos_v(:,:), & ! [IN]
229  atmos_dens(:,:), & ! [IN]
230  atmos_qv(:,:), & ! [IN]
231  real_z1(:,:), & ! [IN]
232  atmos_pbl(:,:), & ! [IN]
233  atmos_sfc_pres(:,:), & ! [IN]
234  atmos_sflx_lw(:,:,:), & ! [IN]
235  atmos_sflx_sw(:,:,:), & ! [IN]
236  atmos_sflx_prec(:,:), & ! [IN]
237  urban_tr(:,:), & ! [IN]
238  urban_tb(:,:), & ! [IN]
239  urban_tg(:,:), & ! [IN]
240  urban_tc(:,:), & ! [IN]
241  urban_qc(:,:), & ! [IN]
242  urban_uc(:,:), & ! [IN]
243  urban_trl(:,:,:), & ! [IN]
244  urban_tbl(:,:,:), & ! [IN]
245  urban_tgl(:,:,:), & ! [IN]
246  urban_rainr(:,:), & ! [IN]
247  urban_rainb(:,:), & ! [IN]
248  urban_raing(:,:), & ! [IN]
249  urban_roff(:,:), & ! [IN]
250  base_lon, & ! [IN]
251  base_lat, & ! [IN]
252  nowdate(:), & ! [IN]
253  dt ) ! [IN]
254 
255  call hydrometeor_lhv( lhv(:,:), atmos_temp(:,:) )
256 
257 !OCL XFILL
258  do j = js, je
259  do i = is, ie
260  urban_sflx_evap(i,j) = urban_sflx_lh(i,j) / lhv(i,j)
261  end do
262  end do
263 
264  call hist_in( urban_tr_t(:,:), 'URBAN_TR_t', 'tendency of URBAN_TR', 'K' )
265  call hist_in( urban_tb_t(:,:), 'URBAN_TB_t', 'tendency of URBAN_TB', 'K' )
266  call hist_in( urban_tg_t(:,:), 'URBAN_TG_t', 'tendency of URBAN_TG', 'K' )
267  call hist_in( urban_tc_t(:,:), 'URBAN_TC_t', 'tendency of URBAN_TC', 'K' )
268  call hist_in( urban_qc_t(:,:), 'URBAN_QC_t', 'tendency of URBAN_QC', 'kg/kg' )
269  call hist_in( urban_uc_t(:,:), 'URBAN_UC_t', 'tendency of URBAN_UC', 'm/s' )
270 
271  call hist_in( urban_trl_t(:,:,:), 'URBAN_TRL_t', 'tendency of URBAN_TRL', 'K', zdim='urban' )
272  call hist_in( urban_tbl_t(:,:,:), 'URBAN_TBL_t', 'tendency of URBAN_TBL', 'K', zdim='urban' )
273  call hist_in( urban_tgl_t(:,:,:), 'URBAN_TGL_t', 'tendency of URBAN_TGL', 'K', zdim='urban' )
274 
275  call hist_in( urban_rainr_t(:,:), 'URBAN_RAINR_t', 'tendency of URBAN_RAINR', 'K' )
276  call hist_in( urban_rainb_t(:,:), 'URBAN_RAINB_t', 'tendency of URBAN_RAINB', 'K' )
277  call hist_in( urban_raing_t(:,:), 'URBAN_RAING_t', 'tendency of URBAN_RAING', 'K' )
278  call hist_in( urban_roff_t(:,:), 'URBAN_ROFF_t', 'tendency of URBAN_ROFF', 'K' )
279 
280  endif
281 
282  if ( statistics_checktotal ) then
283  call stat_total( total, urban_tr_t(:,:), 'URBAN_TR_t' )
284  call stat_total( total, urban_tb_t(:,:), 'URBAN_TB_t' )
285  call stat_total( total, urban_tg_t(:,:), 'URBAN_TG_t' )
286  call stat_total( total, urban_tc_t(:,:), 'URBAN_TC_t' )
287  call stat_total( total, urban_qc_t(:,:), 'URBAN_QC_t' )
288  call stat_total( total, urban_uc_t(:,:), 'URBAN_UC_t' )
289 
290  do k = uks, uke
291  write(sk,'(I2.2)') k
292 
293  call stat_total( total, urban_trl_t(k,:,:), 'URBAN_TRL_t'//sk )
294  call stat_total( total, urban_tbl_t(k,:,:), 'URBAN_TBL_t'//sk )
295  call stat_total( total, urban_tgl_t(k,:,:), 'URBAN_TGL_t'//sk )
296  enddo
297 
298  call stat_total( total, urban_rainr_t(:,:), 'URBAN_RAINR_t' )
299  call stat_total( total, urban_rainb_t(:,:), 'URBAN_RAINB_t' )
300  call stat_total( total, urban_raing_t(:,:), 'URBAN_RAING_t' )
301  call stat_total( total, urban_roff_t(:,:), 'URBAN_ROFF_t' )
302  endif
303 
304  return
305  end subroutine urban_phy_driver
306 
307 end module mod_urban_phy_driver
real(rp), dimension(:,:), allocatable, public urban_qc_t
integer, public is
start point of inner domain: x, local
logical, public statistics_checktotal
calc&report variable totals to logfile?
logical, public urban_sw
integer, public je
end point of inner domain: y, local
integer, public const_i_lw
long-wave radiation index
Definition: scale_const.F90:95
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), public mprj_basepoint_lat
real(rp), dimension(:,:), allocatable, public urban_tg_t
real(rp), dimension(:,:), allocatable, public urban_u10
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:61
module URBAN / Physics
real(rp), dimension(:,:), allocatable, public urban_v10
real(rp), dimension(:,:), allocatable, public urban_raing_t
real(rp), dimension(:,:), allocatable, public urban_z0e
real(rp), dimension(:,:), allocatable, public urban_sflx_mu
real(rp), dimension(:,:), allocatable, public urban_sflx_evap
real(rp), dimension(:,:), allocatable, public urban_z0m
module STDIO
Definition: scale_stdio.F90:12
real(rp), dimension(:,:), allocatable, public urban_tb_t
real(rp), dimension(:,:), allocatable, public urban_t2
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), public mprj_basepoint_lon
procedure(urb), pointer, public urban_phy
real(rp), dimension(:,:), allocatable, public urban_uc
real(rp), dimension(:,:), allocatable, public atmos_sflx_prec
module Statistics
subroutine, public urban_phy_driver_setup
Setup.
module grid index
real(rp), dimension(:,:), allocatable, public atmos_pbl
real(rp), dimension(:,:), allocatable, public urban_sflx_sh
real(rp), dimension(:,:), allocatable, public urban_tr
integer, public ia
of whole cells: x, local, with HALO
real(rp), dimension(:,:,:), allocatable, public urban_tgl
real(rp), dimension(:,:), allocatable, public real_z1
Height of the lowermost grid from surface (cell center) [m].
module GRID (real space)
real(rp), dimension(:,:), allocatable, public atmos_qv
real(rp), dimension(:,:), allocatable, public atmos_sfc_pres
real(rp), dimension(:,:), allocatable, public urban_uc_t
subroutine, public urban_phy_driver(update_flag)
Driver.
real(rp), dimension(:,:), allocatable, public urban_z0h
real(rp), dimension(:,:,:), allocatable, public atmos_sflx_sw
real(rp), dimension(:,:), allocatable, public atmos_v
integer, public js
start point of inner domain: y, local
module TIME
Definition: scale_time.F90:15
real(rp), dimension(:,:), allocatable, public urban_roff
real(rp), dimension(:,:), allocatable, public atmos_dens
module administrator for restart
module CONSTANT
Definition: scale_const.F90:14
real(rp), dimension(:,:,:), allocatable, public urban_sfc_albedo
real(rp), dimension(:,:), allocatable, public urban_sflx_lh
character(len=h_short), public urban_type
subroutine, public urban_phy_driver_resume
Resume.
real(rp), dimension(:,:), allocatable, public urban_tc
real(rp), dimension(:,:), allocatable, public atmos_pres
real(rp), dimension(:,:), allocatable, public urban_rainr
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:156
module profiler
Definition: scale_prof.F90:10
integer, public ie
end point of inner domain: x, local
real(rp), dimension(:,:), allocatable, public atmos_temp
integer, public const_i_sw
short-wave radiation index
Definition: scale_const.F90:96
module PRECISION
subroutine, public urban_phy_setup(URBAN_TYPE, Z0M, Z0H, Z0E)
module HISTORY
real(rp), dimension(:,:), allocatable, public urban_q2
real(rp), dimension(:,:), allocatable, public urban_tg
real(rp), dimension(:,:), allocatable, public atmos_u
real(rp), dimension(:,:), allocatable, public urban_sfc_temp
integer, dimension(6), public time_nowdate
current time [YYYY MM DD HH MM SS]
Definition: scale_time.F90:65
real(rp), dimension(:,:,:), allocatable, public urban_trl
real(rp), dimension(:,:), allocatable, public urban_rainb_t
logical, public restart_run
Is this run restart?
real(rp), dimension(:,:,:), allocatable, public atmos_sflx_lw
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:204
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_sflx_gh
real(rp), dimension(:,:), allocatable, public urban_sflx_mv
real(rp), dimension(:,:), allocatable, public atmos_w
real(rp), dimension(:,:), allocatable, public urban_rainb
real(rp), dimension(:,:), allocatable, public urban_tc_t
real(rp), dimension(:,:), allocatable, public urban_roff_t
real(rp), dimension(:,:,:), allocatable, public urban_tbl_t
real(rp), dimension(:,:), allocatable, public urban_sflx_mw
module Map projection
module urban grid index
integer, public ja
of whole cells: y, local, with HALO