SCALE-RM
mod_atmos_phy_mp_driver.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
12 !-------------------------------------------------------------------------------
13 #include "inc_openmp.h"
15  !-----------------------------------------------------------------------------
16  !
17  !++ used modules
18  !
19  use scale_precision
20  use scale_stdio
21  use scale_prof
23  use scale_tracer
24  !-----------------------------------------------------------------------------
25  implicit none
26  private
27  !-----------------------------------------------------------------------------
28  !
29  !++ Public procedure
30  !
34  public :: atmos_phy_mp_driver
35 
36  !-----------------------------------------------------------------------------
37  !
38  !++ Public parameters & variables
39  !
40  !-----------------------------------------------------------------------------
41  !
42  !++ Private procedure
43  !
44  !-----------------------------------------------------------------------------
45  !
46  !++ Private parameters & variables
47  !
48  !-----------------------------------------------------------------------------
49 contains
50  !-----------------------------------------------------------------------------
53  use scale_atmos_phy_mp, only: &
55  use mod_atmos_admin, only: &
58  implicit none
59  !---------------------------------------------------------------------------
60 
61  if( io_l ) write(io_fid_log,*)
62  if( io_l ) write(io_fid_log,*) '++++++ Module[CONFIG] / Categ[ATMOS PHY_MP] / Origin[SCALE-RM]'
63 
64  if ( atmos_sw_phy_mp ) then
66  end if
67 
68  return
69  end subroutine atmos_phy_mp_driver_config
70 
71  !-----------------------------------------------------------------------------
73  subroutine atmos_phy_mp_driver_setup
74  use scale_atmos_phy_mp, only: &
76  use mod_atmos_admin, only: &
79  use mod_atmos_phy_mp_vars, only: &
80  sflx_rain => atmos_phy_mp_sflx_rain, &
81  sflx_snow => atmos_phy_mp_sflx_snow
82  implicit none
83  !---------------------------------------------------------------------------
84 
85  if( io_l ) write(io_fid_log,*)
86  if( io_l ) write(io_fid_log,*) '++++++ Module[DRIVER] / Categ[ATMOS PHY_MP] / Origin[SCALE-RM]'
87 
88  if ( atmos_sw_phy_mp ) then
89 
90  ! setup library component
92 
93  else
94 
95  if( io_l ) write(io_fid_log,*) '*** this component is never called.'
96  if( io_l ) write(io_fid_log,*) '*** SFLX_rain and SFLX_snow is set to zero.'
97  sflx_rain(:,:) = 0.0_rp
98  sflx_snow(:,:) = 0.0_rp
99 
100  endif
101 
102  return
103  end subroutine atmos_phy_mp_driver_setup
104 
105  !-----------------------------------------------------------------------------
107  subroutine atmos_phy_mp_driver_resume
108  use mod_atmos_admin, only: &
110  implicit none
111 
112  if ( atmos_sw_phy_mp ) then
113 
114  ! run once (only for the diagnostic value)
115  call prof_rapstart('ATM_Microphysics', 1)
116  call atmos_phy_mp_driver( update_flag = .true. )
117  call prof_rapend ('ATM_Microphysics', 1)
118 
119  end if
120 
121  return
122  end subroutine atmos_phy_mp_driver_resume
123 
124  !-----------------------------------------------------------------------------
126  subroutine atmos_phy_mp_driver( update_flag )
127  use scale_time, only: &
128  dt_mp => time_dtsec_atmos_phy_mp
129  use scale_rm_statistics, only: &
131  stat_total
132  use scale_history, only: &
133  hist_in
134  use scale_atmos_hydrometeor, only: &
135  n_hyd, &
136  i_hc, &
137  i_hr, &
138  i_hi, &
139  i_hs, &
140  i_hg, &
141  i_hh
142  use scale_atmos_phy_mp, only: &
143  atmos_phy_mp, &
145  qs_mp, &
146  qe_mp
147  use mod_atmos_vars, only: &
148  dens => dens_av, &
149  momz => momz_av, &
150  momx => momx_av, &
151  momy => momy_av, &
152  rhot => rhot_av, &
153  qtrc => qtrc_av, &
154  dens_t => dens_tp, &
155  momz_t => momz_tp, &
156  momx_t => momx_tp, &
157  momy_t => momy_tp, &
158  rhot_t => rhot_tp, &
159  rhoq_t => rhoq_tp, &
160  temp
161  use mod_atmos_phy_mp_vars, only: &
162  dens_t_mp => atmos_phy_mp_dens_t, &
163  momz_t_mp => atmos_phy_mp_momz_t, &
164  momx_t_mp => atmos_phy_mp_momx_t, &
165  momy_t_mp => atmos_phy_mp_momy_t, &
166  rhot_t_mp => atmos_phy_mp_rhot_t, &
167  rhoq_t_mp => atmos_phy_mp_rhoq_t, &
168  evaporate => atmos_phy_mp_evaporate, &
169  sflx_rain => atmos_phy_mp_sflx_rain, &
170  sflx_snow => atmos_phy_mp_sflx_snow
171  use mod_atmos_phy_ae_vars, only: &
172  ccn_t => atmos_phy_ae_ccn_t
173  implicit none
174 
175  logical, intent(in) :: update_flag
176 
177  real(RP) :: dens0(ka,ia,ja)
178  real(RP) :: momz0(ka,ia,ja)
179  real(RP) :: momx0(ka,ia,ja)
180  real(RP) :: momy0(ka,ia,ja)
181  real(RP) :: rhot0(ka,ia,ja)
182  real(RP) :: qtrc0(ka,ia,ja,qa)
183  real(RP) :: ccn(ka,ia,ja)
184 
185  real(RP) :: re (ka,ia,ja,n_hyd) ! effective radius
186  real(RP) :: precip(ia,ja)
187  real(RP) :: total ! dummy
188 
189  integer :: k, i, j, iq
190  !---------------------------------------------------------------------------
191 
192  if ( update_flag ) then
193 
194 !OCL XFILL
195  !$omp parallel do default(none) &
196  !$omp shared(JA,IA,KA,DENS0,MOMZ0,MOMX0,MOMY0,RHOT0,DENS,MOMZ,MOMX,MOMY,RHOT) &
197  !$omp private(i,j,k) OMP_SCHEDULE_ collapse(2)
198  do j = 1, ja
199  do i = 1, ia
200  do k = 1, ka
201  dens0(k,i,j) = dens(k,i,j) ! save
202  momz0(k,i,j) = momz(k,i,j) ! save
203  momx0(k,i,j) = momx(k,i,j) ! save
204  momy0(k,i,j) = momy(k,i,j) ! save
205  rhot0(k,i,j) = rhot(k,i,j) ! save
206  enddo
207  enddo
208  enddo
209 
210 !OCL XFILL
211  do iq = 1, qa
212  do j = 1, ja
213  do i = 1, ia
214  do k = 1, ka
215  qtrc0(k,i,j,iq) = qtrc(k,i,j,iq) ! save
216  enddo
217  enddo
218  enddo
219  enddo
220 
221  ccn(:,:,:) = ccn_t(:,:,:) * dt_mp
222 
223  call atmos_phy_mp( dens0(:,:,:), & ! [INOUT]
224  momz0(:,:,:), & ! [INOUT]
225  momx0(:,:,:), & ! [INOUT]
226  momy0(:,:,:), & ! [INOUT]
227  rhot0(:,:,:), & ! [INOUT]
228  qtrc0(:,:,:,:), & ! [INOUT]
229  ccn(:,:,:), & ! [IN]
230  evaporate(:,:,:), & ! [OUT]
231  sflx_rain(:,:), & ! [OUT]
232  sflx_snow(:,:) ) ! [OUT]
233 
234 !OCL XFILL
235  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
236  !$omp shared(JS,JE,IS,IE,KS,KE,DENS_t_MP,DENS0,DENS,MOMZ_t_MP,MOMZ0,MOMZ,MOMX_t_MP,MOMX0) &
237  !$omp shared(MOMX,MOMY_t_MP,MOMY0,MOMY,RHOT_t_MP,RHOT0,RHOT,dt_MP)
238  do j = js, je
239  do i = is, ie
240  do k = ks, ke
241  dens_t_mp(k,i,j) = ( dens0(k,i,j) - dens(k,i,j) ) / dt_mp
242  momz_t_mp(k,i,j) = ( momz0(k,i,j) - momz(k,i,j) ) / dt_mp
243  momx_t_mp(k,i,j) = ( momx0(k,i,j) - momx(k,i,j) ) / dt_mp
244  momy_t_mp(k,i,j) = ( momy0(k,i,j) - momy(k,i,j) ) / dt_mp
245  rhot_t_mp(k,i,j) = ( rhot0(k,i,j) - rhot(k,i,j) ) / dt_mp
246  enddo
247  enddo
248  enddo
249 
250 !OCL XFILL
251  do iq = qs_mp, qe_mp
252  !$omp parallel do default(none) &
253  !$omp shared(JS,JE,IS,IE,KS,KE,RHOQ_t_MP,iq,QTRC0,QTRC,DENS0,DENS,dt_MP) &
254  !$omp private(i,j,k) OMP_SCHEDULE_ collapse(2)
255  do j = js, je
256  do i = is, ie
257  do k = ks, ke
258  rhoq_t_mp(k,i,j,iq) = ( qtrc0(k,i,j,iq) * dens0(k,i,j) &
259  - qtrc(k,i,j,iq) * dens(k,i,j) ) / dt_mp
260  enddo
261  enddo
262  enddo
263  enddo
264 
265 !OCL XFILL
266  do j = js, je
267  do i = is, ie
268  precip(i,j) = sflx_rain(i,j) + sflx_snow(i,j)
269  end do
270  end do
271 
272  call hist_in( sflx_rain(:,:), 'RAIN_MP', 'surface rain rate by MP', 'kg/m2/s', nohalo=.true. )
273  call hist_in( sflx_snow(:,:), 'SNOW_MP', 'surface snow rate by MP', 'kg/m2/s', nohalo=.true. )
274  call hist_in( precip(:,:), 'PREC_MP', 'surface precipitation rate by MP', 'kg/m2/s', nohalo=.true. )
275  call hist_in( evaporate(:,:,:), 'EVAPORATE', 'evaporated cloud number', 'num/m3/s', nohalo=.true. )
276 
277  call hist_in( dens_t_mp(:,:,:), 'DENS_t_MP', 'tendency DENS in MP', 'kg/m3/s' , nohalo=.true. )
278  call hist_in( momz_t_mp(:,:,:), 'MOMZ_t_MP', 'tendency MOMZ in MP', 'kg/m2/s2' , nohalo=.true. )
279  call hist_in( momx_t_mp(:,:,:), 'MOMX_t_MP', 'tendency MOMX in MP', 'kg/m2/s2' , nohalo=.true. )
280  call hist_in( momy_t_mp(:,:,:), 'MOMY_t_MP', 'tendency MOMY in MP', 'kg/m2/s2' , nohalo=.true. )
281  call hist_in( rhot_t_mp(:,:,:), 'RHOT_t_MP', 'tendency RHOT in MP', 'K*kg/m3/s', nohalo=.true. )
282 
283  do iq = qs_mp, qe_mp
284  call hist_in( rhoq_t_mp(:,:,:,iq), trim(tracer_name(iq))//'_t_MP', &
285  'tendency rho*'//trim(tracer_name(iq))//'in MP', 'kg/m3/s', nohalo=.true. )
286  enddo
287 
288  call atmos_phy_mp_effectiveradius( re(:,:,:,:), & ! [OUT]
289  qtrc(:,:,:,:), & ! [IN]
290  dens(:,:,:), & ! [IN]
291  temp(:,:,:) ) ! [IN]
292 
293  call hist_in( re(:,:,:,i_hc), 'Re_QC', 'effective radius cloud water', 'cm', nohalo=.true. )
294  call hist_in( re(:,:,:,i_hr), 'Re_QR', 'effective radius rain', 'cm', nohalo=.true. )
295  call hist_in( re(:,:,:,i_hi), 'Re_QI', 'effective radius ice water', 'cm', nohalo=.true. )
296  call hist_in( re(:,:,:,i_hs), 'Re_QS', 'effective radius snow', 'cm', nohalo=.true. )
297  call hist_in( re(:,:,:,i_hg), 'Re_QG', 'effective radius graupel', 'cm', nohalo=.true. )
298  call hist_in( re(:,:,:,i_hh), 'Re_QH', 'effective radius hail', 'cm', nohalo=.true. )
299  endif
300 
301  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
302  !$omp shared(JS,JE,IS,IE,KS,KE,DENS_t,DENS_t_MP,MOMZ_t,MOMZ_t_MP,MOMX_t,MOMX_t_MP,MOMY_t) &
303  !$omp shared(MOMY_t_MP,RHOT_t,RHOT_t_MP)
304  do j = js, je
305  do i = is, ie
306  do k = ks, ke
307  dens_t(k,i,j) = dens_t(k,i,j) + dens_t_mp(k,i,j)
308  momz_t(k,i,j) = momz_t(k,i,j) + momz_t_mp(k,i,j)
309  momx_t(k,i,j) = momx_t(k,i,j) + momx_t_mp(k,i,j)
310  momy_t(k,i,j) = momy_t(k,i,j) + momy_t_mp(k,i,j)
311  rhot_t(k,i,j) = rhot_t(k,i,j) + rhot_t_mp(k,i,j)
312  enddo
313  enddo
314  enddo
315 
316  do iq = qs_mp, qe_mp
317  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ &
318  !$omp shared(JS,JE,IS,IE,KS,KE,RHOQ_t,iq,RHOQ_t_MP)
319  do j = js, je
320  do i = is, ie
321  do k = ks, ke
322  rhoq_t(k,i,j,iq) = rhoq_t(k,i,j,iq) + rhoq_t_mp(k,i,j,iq)
323  enddo
324  enddo
325  enddo
326  enddo
327 
328  if ( statistics_checktotal ) then
329  call stat_total( total, dens_t_mp(:,:,:), 'DENS_t_MP' )
330  call stat_total( total, momz_t_mp(:,:,:), 'MOMZ_t_MP' )
331  call stat_total( total, momx_t_mp(:,:,:), 'MOMX_t_MP' )
332  call stat_total( total, momy_t_mp(:,:,:), 'MOMY_t_MP' )
333  call stat_total( total, rhot_t_mp(:,:,:), 'RHOT_t_MP' )
334 
335  do iq = qs_mp, qe_mp
336  call stat_total( total, rhoq_t_mp(:,:,:,iq), trim(tracer_name(iq))//'_t_MP' )
337  enddo
338  endif
339 
340  return
341  end subroutine atmos_phy_mp_driver
342 
343 end module mod_atmos_phy_mp_driver
module ATMOS admin
real(rp), dimension(:,:,:), allocatable, public dens_tp
integer, public is
start point of inner domain: x, local
real(rp), dimension(:,:,:,:), allocatable, public atmos_phy_mp_rhoq_t
logical, public statistics_checktotal
calc&report variable totals to logfile?
real(rp), dimension(:,:,:), allocatable, target, public momz
procedure(mp), pointer, public atmos_phy_mp
real(rp), dimension(:,:,:), allocatable, public atmos_phy_mp_momy_t
integer, public je
end point of inner domain: y, local
real(rp), dimension(:,:,:), allocatable, target, public rhot
real(rp), dimension(:,:,:), allocatable, public momy_tp
real(dp), public time_dtsec_atmos_phy_mp
time interval of physics(microphysics) [sec]
Definition: scale_time.F90:41
module Atmosphere / Physics Cloud Microphysics
integer, parameter, public i_hs
snow
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:61
module ATMOSPHERE / Physics Cloud Microphysics
module ATMOSPHERIC Variables
real(rp), dimension(:,:,:,:), pointer, public qtrc_av
integer, parameter, public i_hr
liquid water rain
real(rp), dimension(:,:,:), allocatable, target, public momx
subroutine, public atmos_phy_mp_driver(update_flag)
Driver.
integer, parameter, public i_hi
ice water cloud
module STDIO
Definition: scale_stdio.F90:12
integer, public ke
end point of inner domain: z, local
real(rp), dimension(:,:,:), allocatable, public rhot_tp
integer, public qa
real(rp), dimension(:,:,:), allocatable, target, public dens
integer, parameter, public i_hh
hail
real(rp), dimension(:,:,:), allocatable, public atmos_phy_mp_momz_t
character(len=h_short), dimension(qa_max), public tracer_name
real(rp), dimension(:,:), allocatable, public atmos_phy_mp_sflx_rain
module Statistics
module grid index
real(rp), dimension(:,:,:), pointer, public momx_av
module TRACER
integer, public ia
of whole cells: x, local, with HALO
subroutine, public atmos_phy_mp_config(MP_TYPE)
Setup Cloud Microphysics.
real(rp), dimension(:,:,:), allocatable, public atmos_phy_mp_rhot_t
real(rp), dimension(:,:,:), allocatable, public temp
integer, public ka
of whole cells: z, local, with HALO
real(rp), dimension(:,:,:), allocatable, public atmos_phy_mp_dens_t
real(rp), dimension(:,:,:,:), allocatable, public rhoq_tp
integer, public js
start point of inner domain: y, local
module TIME
Definition: scale_time.F90:15
procedure(er), pointer, public atmos_phy_mp_effectiveradius
procedure(su), pointer, public atmos_phy_mp_setup
subroutine, public atmos_phy_mp_driver_setup
Setup.
real(rp), dimension(:,:,:), pointer, public dens_av
logical, public atmos_sw_phy_mp
integer, parameter, public i_hc
liquid water cloud
character(len=h_short), public atmos_phy_mp_type
integer, public ks
start point of inner domain: z, local
real(rp), dimension(:,:,:), allocatable, public momx_tp
real(rp), dimension(:,:,:), allocatable, target, public momy
real(rp), dimension(:,:,:), allocatable, public atmos_phy_mp_evaporate
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:156
real(rp), dimension(:,:), allocatable, public atmos_phy_mp_sflx_snow
real(rp), dimension(:,:,:), allocatable, public momz_tp
module profiler
Definition: scale_prof.F90:10
integer, public ie
end point of inner domain: x, local
module PRECISION
module HISTORY
real(rp), dimension(:,:,:), pointer, public momz_av
real(rp), dimension(:,:,:), pointer, public rhot_av
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, parameter, public n_hyd
subroutine, public atmos_phy_mp_driver_resume
resume
real(rp), dimension(:,:,:), pointer, public momy_av
module ATMOSPHERE / Physics Cloud Microphysics
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:204
module ATMOSPHERE / Physics Aerosol Microphysics
integer, parameter, public i_hg
graupel
subroutine, public atmos_phy_mp_driver_config
Config.
real(rp), dimension(:,:,:), allocatable, public atmos_phy_mp_momx_t
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc
real(rp), dimension(:,:,:), allocatable, public atmos_phy_ae_ccn_t
integer, public ja
of whole cells: y, local, with HALO