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  implicit none
58  !---------------------------------------------------------------------------
59 
60  if( io_l ) write(io_fid_log,*)
61  if( io_l ) write(io_fid_log,*) '++++++ Module[DRIVER] / Categ[URBAN PHY] / Origin[SCALE-RM]'
62 
63  if ( urban_sw ) then
64 
65  ! setup library component
67 
68  else
69  if( io_l ) write(io_fid_log,*) '*** this component is never called.'
70  endif
71 
72  return
73  end subroutine urban_phy_driver_setup
74 
75  !-----------------------------------------------------------------------------
77  subroutine urban_phy_driver_resume
78  use mod_admin_restart, only: &
80  use mod_urban_admin, only: &
81  urban_sw
82  implicit none
83 
84  if ( urban_sw ) then
85 
86  if ( .NOT. restart_run ) then ! tentative
87  ! run once (only for the diagnostic value)
88  call prof_rapstart('URB_Physics', 1)
89  call urban_phy_driver( update_flag = .true. )
90  call prof_rapend ('URB_Physics', 1)
91  end if
92 
93  end if
94 
95  return
96  end subroutine urban_phy_driver_resume
97  !-----------------------------------------------------------------------------
99  subroutine urban_phy_driver( update_flag )
101  atmos_thermodyn_templhv
102  use scale_time, only: &
103  nowdate => time_nowdate, &
104  dt => time_dtsec_urban
105  use scale_rm_statistics, only: &
107  stat_total
108  use scale_history, only: &
109  hist_in
110  use scale_mapproj, only: &
111  base_lon => mprj_basepoint_lon, &
112  base_lat => mprj_basepoint_lat
113  use scale_grid_real, only: &
114  real_z1
115  use scale_urban_phy, only: &
116  urban_phy
117  use mod_urban_vars, only: &
118  urban_tr, &
119  urban_tb, &
120  urban_tg, &
121  urban_tc, &
122  urban_qc, &
123  urban_uc, &
124  urban_trl, &
125  urban_tbl, &
126  urban_tgl, &
127  urban_rainr, &
128  urban_rainb, &
129  urban_raing, &
130  urban_roff, &
131  urban_tr_t, &
132  urban_tb_t, &
133  urban_tg_t, &
134  urban_tc_t, &
135  urban_qc_t, &
136  urban_uc_t, &
137  urban_trl_t, &
138  urban_tbl_t, &
139  urban_tgl_t, &
140  urban_rainr_t, &
141  urban_rainb_t, &
142  urban_raing_t, &
143  urban_roff_t, &
144  urban_sfc_temp, &
146  urban_sflx_mw, &
147  urban_sflx_mu, &
148  urban_sflx_mv, &
149  urban_sflx_sh, &
150  urban_sflx_lh, &
151  urban_sflx_gh, &
152  urban_sflx_evap, &
153  urban_z0m, &
154  urban_z0h, &
155  urban_z0e, &
156  urban_u10, &
157  urban_v10, &
158  urban_t2, &
159  urban_q2, &
160  atmos_temp, &
161  atmos_pres, &
162  atmos_w, &
163  atmos_u, &
164  atmos_v, &
165  atmos_dens, &
166  atmos_qv, &
167  atmos_pbl, &
168  atmos_sfc_pres, &
169  atmos_sflx_lw, &
170  atmos_sflx_sw, &
172  implicit none
173 
174  ! arguments
175  logical, intent(in) :: update_flag
176 
177  ! works
178  real(RP) :: total ! dummy
179  real(RP) :: lhv(ia,ja)
180 
181  character(len=2) :: sk
182 
183  integer :: k, i, j
184  !---------------------------------------------------------------------------
185 
186  if ( update_flag ) then
187 
188  call urban_phy( urban_tr_t(:,:), & ! [OUT]
189  urban_tb_t(:,:), & ! [OUT]
190  urban_tg_t(:,:), & ! [OUT]
191  urban_tc_t(:,:), & ! [OUT]
192  urban_qc_t(:,:), & ! [OUT]
193  urban_uc_t(:,:), & ! [OUT]
194  urban_trl_t(:,:,:), & ! [OUT]
195  urban_tbl_t(:,:,:), & ! [OUT]
196  urban_tgl_t(:,:,:), & ! [OUT]
197  urban_rainr_t(:,:), & ! [OUT]
198  urban_rainb_t(:,:), & ! [OUT]
199  urban_raing_t(:,:), & ! [OUT]
200  urban_roff_t(:,:), & ! [OUT]
201  urban_sfc_temp(:,:), & ! [OUT]
202  urban_sfc_albedo(:,:,i_lw), & ! [OUT]
203  urban_sfc_albedo(:,:,i_sw), & ! [OUT]
204  urban_sflx_mw(:,:), & ! [OUT]
205  urban_sflx_mu(:,:), & ! [OUT]
206  urban_sflx_mv(:,:), & ! [OUT]
207  urban_sflx_sh(:,:), & ! [OUT]
208  urban_sflx_lh(:,:), & ! [OUT]
209  urban_sflx_gh(:,:), & ! [OUT]
210  urban_z0m(:,:), & ! [OUT]
211  urban_z0h(:,:), & ! [OUT]
212  urban_z0e(:,:), & ! [OUT]
213  urban_u10(:,:), & ! [OUT]
214  urban_v10(:,:), & ! [OUT]
215  urban_t2(:,:), & ! [OUT]
216  urban_q2(:,:), & ! [OUT]
217  atmos_temp(:,:), & ! [IN]
218  atmos_pres(:,:), & ! [IN]
219  atmos_w(:,:), & ! [IN]
220  atmos_u(:,:), & ! [IN]
221  atmos_v(:,:), & ! [IN]
222  atmos_dens(:,:), & ! [IN]
223  atmos_qv(:,:), & ! [IN]
224  real_z1(:,:), & ! [IN]
225  atmos_pbl(:,:), & ! [IN]
226  atmos_sfc_pres(:,:), & ! [IN]
227  atmos_sflx_lw(:,:,:), & ! [IN]
228  atmos_sflx_sw(:,:,:), & ! [IN]
229  atmos_sflx_prec(:,:), & ! [IN]
230  urban_tr(:,:), & ! [IN]
231  urban_tb(:,:), & ! [IN]
232  urban_tg(:,:), & ! [IN]
233  urban_tc(:,:), & ! [IN]
234  urban_qc(:,:), & ! [IN]
235  urban_uc(:,:), & ! [IN]
236  urban_trl(:,:,:), & ! [IN]
237  urban_tbl(:,:,:), & ! [IN]
238  urban_tgl(:,:,:), & ! [IN]
239  urban_rainr(:,:), & ! [IN]
240  urban_rainb(:,:), & ! [IN]
241  urban_raing(:,:), & ! [IN]
242  urban_roff(:,:), & ! [IN]
243  base_lon, & ! [IN]
244  base_lat, & ! [IN]
245  nowdate(:), & ! [IN]
246  dt ) ! [IN]
247 
248  call atmos_thermodyn_templhv( lhv, atmos_temp )
249 
250 !OCL XFILL
251  do j = js, je
252  do i = is, ie
253  urban_sflx_evap(i,j) = urban_sflx_lh(i,j) / lhv(i,j)
254  end do
255  end do
256 
257  call hist_in( urban_tr_t(:,:), 'URBAN_TR_t', 'tendency of URBAN_TR', 'K' )
258  call hist_in( urban_tb_t(:,:), 'URBAN_TB_t', 'tendency of URBAN_TB', 'K' )
259  call hist_in( urban_tg_t(:,:), 'URBAN_TG_t', 'tendency of URBAN_TG', 'K' )
260  call hist_in( urban_tc_t(:,:), 'URBAN_TC_t', 'tendency of URBAN_TC', 'K' )
261  call hist_in( urban_qc_t(:,:), 'URBAN_QC_t', 'tendency of URBAN_QC', 'kg/kg' )
262  call hist_in( urban_uc_t(:,:), 'URBAN_UC_t', 'tendency of URBAN_UC', 'm/s' )
263 
264  call hist_in( urban_trl_t(:,:,:), 'URBAN_TRL_t', 'tendency of URBAN_TRL', 'K', zdim='urban' )
265  call hist_in( urban_tbl_t(:,:,:), 'URBAN_TBL_t', 'tendency of URBAN_TBL', 'K', zdim='urban' )
266  call hist_in( urban_tgl_t(:,:,:), 'URBAN_TGL_t', 'tendency of URBAN_TGL', 'K', zdim='urban' )
267 
268  call hist_in( urban_rainr_t(:,:), 'URBAN_RAINR_t', 'tendency of URBAN_RAINR', 'K' )
269  call hist_in( urban_rainb_t(:,:), 'URBAN_RAINB_t', 'tendency of URBAN_RAINB', 'K' )
270  call hist_in( urban_raing_t(:,:), 'URBAN_RAING_t', 'tendency of URBAN_RAING', 'K' )
271  call hist_in( urban_roff_t(:,:), 'URBAN_ROFF_t', 'tendency of URBAN_ROFF', 'K' )
272 
273  endif
274 
275  if ( statistics_checktotal ) then
276  call stat_total( total, urban_tr_t(:,:), 'URBAN_TR_t' )
277  call stat_total( total, urban_tb_t(:,:), 'URBAN_TB_t' )
278  call stat_total( total, urban_tg_t(:,:), 'URBAN_TG_t' )
279  call stat_total( total, urban_tc_t(:,:), 'URBAN_TC_t' )
280  call stat_total( total, urban_qc_t(:,:), 'URBAN_QC_t' )
281  call stat_total( total, urban_uc_t(:,:), 'URBAN_UC_t' )
282 
283  do k = uks, uke
284  write(sk,'(I2.2)') k
285 
286  call stat_total( total, urban_trl_t(k,:,:), 'URBAN_TRL_t'//sk )
287  call stat_total( total, urban_tbl_t(k,:,:), 'URBAN_TBL_t'//sk )
288  call stat_total( total, urban_tgl_t(k,:,:), 'URBAN_TGL_t'//sk )
289  enddo
290 
291  call stat_total( total, urban_rainr_t(:,:), 'URBAN_RAINR_t' )
292  call stat_total( total, urban_rainb_t(:,:), 'URBAN_RAINB_t' )
293  call stat_total( total, urban_raing_t(:,:), 'URBAN_RAING_t' )
294  call stat_total( total, urban_roff_t(:,:), 'URBAN_ROFF_t' )
295  endif
296 
297  return
298  end subroutine urban_phy_driver
299 
300 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:98
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:59
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)
subroutine, public urban_phy_setup(URBAN_TYPE)
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 x whole cells (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:132
module profiler
Definition: scale_prof.F90:10
integer, public ie
end point of inner domain: x, local
module ATMOSPHERE / Thermodynamics
real(rp), dimension(:,:), allocatable, public atmos_temp
integer, public const_i_sw
short-wave radiation index
Definition: scale_const.F90:99
module PRECISION
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: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_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 y whole cells (local, with HALO)