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  !
33  public :: atmos_phy_mp_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 atmos_phy_mp_driver_setup
52  use scale_atmos_phy_mp, only: &
54  use mod_atmos_admin, only: &
57  use mod_atmos_phy_mp_vars, only: &
58  sflx_rain => atmos_phy_mp_sflx_rain, &
59  sflx_snow => atmos_phy_mp_sflx_snow
60  implicit none
61  !---------------------------------------------------------------------------
62 
63  if( io_l ) write(io_fid_log,*)
64  if( io_l ) write(io_fid_log,*) '++++++ Module[DRIVER] / Categ[ATMOS PHY_MP] / Origin[SCALE-RM]'
65 
66  if ( atmos_sw_phy_mp ) then
67 
68  ! setup library component
70 
71  else
72 
73  if( io_l ) write(io_fid_log,*) '*** this component is never called.'
74  if( io_l ) write(io_fid_log,*) '*** SFLX_rain and SFLX_snow is set to zero.'
75  sflx_rain(:,:) = 0.0_rp
76  sflx_snow(:,:) = 0.0_rp
77 
78  endif
79 
80  return
81  end subroutine atmos_phy_mp_driver_setup
82 
83  !-----------------------------------------------------------------------------
86  use mod_atmos_admin, only: &
88  implicit none
89 
90  if ( atmos_sw_phy_mp ) then
91 
92  ! run once (only for the diagnostic value)
93  call prof_rapstart('ATM_Microphysics', 1)
94  call atmos_phy_mp_driver( update_flag = .true. )
95  call prof_rapend ('ATM_Microphysics', 1)
96 
97  end if
98 
99  return
100  end subroutine atmos_phy_mp_driver_resume
101 
102  !-----------------------------------------------------------------------------
104  subroutine atmos_phy_mp_driver( update_flag )
105  use scale_time, only: &
106  dt_mp => time_dtsec_atmos_phy_mp
107  use scale_rm_statistics, only: &
109  stat_total
110  use scale_history, only: &
111  hist_in
112  use scale_atmos_phy_mp, only: &
114  use mod_atmos_vars, only: &
115  dens, &
116  momz, &
117  momx, &
118  momy, &
119  rhot, &
120  qtrc, &
121  dens_t => dens_tp, &
122  momz_t => momz_tp, &
123  momx_t => momx_tp, &
124  momy_t => momy_tp, &
125  rhot_t => rhot_tp, &
126  rhoq_t => rhoq_tp
127  use mod_atmos_phy_mp_vars, only: &
128  dens_t_mp => atmos_phy_mp_dens_t, &
129  momz_t_mp => atmos_phy_mp_momz_t, &
130  momx_t_mp => atmos_phy_mp_momx_t, &
131  momy_t_mp => atmos_phy_mp_momy_t, &
132  rhot_t_mp => atmos_phy_mp_rhot_t, &
133  rhoq_t_mp => atmos_phy_mp_rhoq_t, &
134  evaporate => atmos_phy_mp_evaporate, &
135  sflx_rain => atmos_phy_mp_sflx_rain, &
136  sflx_snow => atmos_phy_mp_sflx_snow
137  use mod_atmos_phy_ae_vars, only: &
138  ccn_t => atmos_phy_ae_ccn_t
139  implicit none
140 
141  logical, intent(in) :: update_flag
142 
143  real(RP) :: DENS0(ka,ia,ja)
144  real(RP) :: MOMZ0(ka,ia,ja)
145  real(RP) :: MOMX0(ka,ia,ja)
146  real(RP) :: MOMY0(ka,ia,ja)
147  real(RP) :: RHOT0(ka,ia,ja)
148  real(RP) :: QTRC0(ka,ia,ja,qa)
149  real(RP) :: CCN(ka,ia,ja)
150 
151  real(RP) :: precip(ia,ja)
152  real(RP) :: total ! dummy
153 
154  integer :: k, i, j, iq
155  !---------------------------------------------------------------------------
156 
157  if ( update_flag ) then
158 
159 !OCL XFILL
160  do j = 1, ja
161  do i = 1, ia
162  do k = 1, ka
163  dens0(k,i,j) = dens(k,i,j) ! save
164  momz0(k,i,j) = momz(k,i,j) ! save
165  momx0(k,i,j) = momx(k,i,j) ! save
166  momy0(k,i,j) = momy(k,i,j) ! save
167  rhot0(k,i,j) = rhot(k,i,j) ! save
168  enddo
169  enddo
170  enddo
171 
172 !OCL XFILL
173  do iq = 1, qa
174  do j = 1, ja
175  do i = 1, ia
176  do k = 1, ka
177  qtrc0(k,i,j,iq) = qtrc(k,i,j,iq) ! save
178  enddo
179  enddo
180  enddo
181  enddo
182 
183  ccn(:,:,:) = ccn_t(:,:,:) * dt_mp
184 
185  call atmos_phy_mp( dens0(:,:,:), & ! [INOUT]
186  momz0(:,:,:), & ! [INOUT]
187  momx0(:,:,:), & ! [INOUT]
188  momy0(:,:,:), & ! [INOUT]
189  rhot0(:,:,:), & ! [INOUT]
190  qtrc0(:,:,:,:), & ! [INOUT]
191  ccn(:,:,:), & ! [IN]
192  evaporate(:,:,:), & ! [OUT]
193  sflx_rain(:,:), & ! [OUT]
194  sflx_snow(:,:) ) ! [OUT]
195 
196 !OCL XFILL
197  do j = js, je
198  do i = is, ie
199  do k = ks, ke
200  dens_t_mp(k,i,j) = ( dens0(k,i,j) - dens(k,i,j) ) / dt_mp
201  momz_t_mp(k,i,j) = ( momz0(k,i,j) - momz(k,i,j) ) / dt_mp
202  momx_t_mp(k,i,j) = ( momx0(k,i,j) - momx(k,i,j) ) / dt_mp
203  momy_t_mp(k,i,j) = ( momy0(k,i,j) - momy(k,i,j) ) / dt_mp
204  rhot_t_mp(k,i,j) = ( rhot0(k,i,j) - rhot(k,i,j) ) / dt_mp
205  enddo
206  enddo
207  enddo
208 
209 !OCL XFILL
210  do iq = 1, qa
211  do j = js, je
212  do i = is, ie
213  do k = ks, ke
214  rhoq_t_mp(k,i,j,iq) = ( qtrc0(k,i,j,iq) * dens0(k,i,j) &
215  - qtrc(k,i,j,iq) * dens(k,i,j) ) / dt_mp
216  enddo
217  enddo
218  enddo
219  enddo
220 
221 !OCL XFILL
222  do j = js, je
223  do i = is, ie
224  precip(i,j) = sflx_rain(i,j) + sflx_snow(i,j)
225  end do
226  end do
227 
228  call hist_in( sflx_rain(:,:), 'RAIN', 'surface rain rate', 'kg/m2/s', nohalo=.true. )
229  call hist_in( sflx_snow(:,:), 'SNOW', 'surface snow rate', 'kg/m2/s', nohalo=.true. )
230  call hist_in( precip(:,:), 'PREC', 'surface precipitation rate', 'kg/m2/s', nohalo=.true. )
231  call hist_in( evaporate(:,:,:), 'EVAPORATE', 'evaporated cloud number', 'num/m3/s', nohalo=.true. )
232 
233  call hist_in( dens_t_mp(:,:,:), 'DENS_t_MP', 'tendency DENS in MP', 'kg/m3/s' , nohalo=.true. )
234  call hist_in( momz_t_mp(:,:,:), 'MOMZ_t_MP', 'tendency MOMZ in MP', 'kg/m2/s2' , nohalo=.true. )
235  call hist_in( momx_t_mp(:,:,:), 'MOMX_t_MP', 'tendency MOMX in MP', 'kg/m2/s2' , nohalo=.true. )
236  call hist_in( momy_t_mp(:,:,:), 'MOMY_t_MP', 'tendency MOMY in MP', 'kg/m2/s2' , nohalo=.true. )
237  call hist_in( rhot_t_mp(:,:,:), 'RHOT_t_MP', 'tendency RHOT in MP', 'K*kg/m3/s', nohalo=.true. )
238 
239  do iq = 1, qa
240  call hist_in( rhoq_t_mp(:,:,:,iq), trim(aq_name(iq))//'_t_MP', &
241  'tendency rho*'//trim(aq_name(iq))//'in MP', 'kg/m3/s', nohalo=.true. )
242  enddo
243 
244  endif
245 
246  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
247  do j = js, je
248  do i = is, ie
249  do k = ks, ke
250  dens_t(k,i,j) = dens_t(k,i,j) + dens_t_mp(k,i,j)
251  momz_t(k,i,j) = momz_t(k,i,j) + momz_t_mp(k,i,j)
252  momx_t(k,i,j) = momx_t(k,i,j) + momx_t_mp(k,i,j)
253  momy_t(k,i,j) = momy_t(k,i,j) + momy_t_mp(k,i,j)
254  rhot_t(k,i,j) = rhot_t(k,i,j) + rhot_t_mp(k,i,j)
255  enddo
256  enddo
257  enddo
258 
259  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(3)
260  do iq = 1, qa
261  do j = js, je
262  do i = is, ie
263  do k = ks, ke
264  rhoq_t(k,i,j,iq) = rhoq_t(k,i,j,iq) + rhoq_t_mp(k,i,j,iq)
265  enddo
266  enddo
267  enddo
268  enddo
269 
270  if ( statistics_checktotal ) then
271  call stat_total( total, dens_t_mp(:,:,:), 'DENS_t_MP' )
272  call stat_total( total, momz_t_mp(:,:,:), 'MOMZ_t_MP' )
273  call stat_total( total, momx_t_mp(:,:,:), 'MOMX_t_MP' )
274  call stat_total( total, momy_t_mp(:,:,:), 'MOMY_t_MP' )
275  call stat_total( total, rhot_t_mp(:,:,:), 'RHOT_t_MP' )
276 
277  do iq = 1, qa
278  call stat_total( total, rhoq_t_mp(:,:,:,iq), trim(aq_name(iq))//'_t_MP' )
279  enddo
280  endif
281 
282  return
283  end subroutine atmos_phy_mp_driver
284 
285 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
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
module ATMOSPHERE / Physics Cloud Microphysics
module ATMOSPHERIC Variables
real(rp), dimension(:,:,:), allocatable, target, public momx
subroutine, public atmos_phy_mp_driver(update_flag)
Driver.
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
real(rp), dimension(:,:,:), allocatable, public atmos_phy_mp_momz_t
real(rp), dimension(:,:), allocatable, public atmos_phy_mp_sflx_rain
module Statistics
subroutine, public atmos_phy_mp_setup(MP_TYPE)
Setup Cloud Microphysics.
module grid index
module TRACER
integer, public ia
of x whole cells (local, with HALO)
real(rp), dimension(:,:,:), allocatable, public atmos_phy_mp_rhot_t
integer, public ka
of z whole cells (local, with HALO)
real(rp), dimension(:,:,:), allocatable, public atmos_phy_mp_dens_t
real(rp), dimension(:,:,:,:), allocatable, public rhoq_tp
character(len=h_short), dimension(:), allocatable, public aq_name
integer, public js
start point of inner domain: y, local
module TIME
Definition: scale_time.F90:15
subroutine, public atmos_phy_mp_driver_setup
Setup.
logical, public atmos_sw_phy_mp
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:132
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
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
subroutine, public atmos_phy_mp_driver_resume
resume
module ATMOSPHERE / Physics Cloud Microphysics
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:178
module ATMOSPHERE / Physics Aerosol Microphysics
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 y whole cells (local, with HALO)