SCALE-RM
mod_atmos_phy_cp_driver.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
11 #include "scalelib.h"
13  !-----------------------------------------------------------------------------
14  !
15  !++ used modules
16  !
17  use scale_precision
18  use scale_io
19  use scale_prof
21  use scale_tracer
22  !-----------------------------------------------------------------------------
23  implicit none
24  private
25  !-----------------------------------------------------------------------------
26  !
27  !++ Public procedure
28  !
32 
33  !-----------------------------------------------------------------------------
34  !
35  !++ Public parameters & variables
36  !
37  !-----------------------------------------------------------------------------
38  !
39  !++ Private procedure
40  !
41  !-----------------------------------------------------------------------------
42  !
43  !++ Private parameters & variables
44  !
45  !-----------------------------------------------------------------------------
46 contains
47  !-----------------------------------------------------------------------------
49  subroutine atmos_phy_cp_driver_setup
50  use scale_prc, only: &
51  prc_abort
52  use scale_atmos_phy_cp_common, only: &
54  use scale_atmos_phy_cp_kf, only: &
58  use mod_atmos_admin, only: &
61  use scale_time , only :&
63  use scale_atmos_grid_cartesc, only: &
64  dx
68  use scale_atmos_hydrometeor, only: &
70  implicit none
71 
72  logical :: warmrain
73  !---------------------------------------------------------------------------
74 
75  log_newline
76  log_info("ATMOS_PHY_CP_driver_setup",*) 'Setup'
77 
78  if ( atmos_sw_phy_cp ) then
79 
80  ! setup library component
82  select case ( atmos_phy_cp_type )
83  case ( 'KF' )
84  warmrain = ( .not. atmos_hydrometeor_ice_phase )
85  call atmos_phy_cp_kf_setup( ka, ks, ke, ia, is, ie, ja, js, je, &
86  cz, area, &
87  warmrain )
88  case ( 'KF-JMAPPLIB' )
90  dx, dt )
91  case default
92  log_error("ATMOS_PHY_CP_driver_setup",*) 'ATMOS_PHY_CP_TYPE (', trim(atmos_phy_cp_type), ') is invalid. Check!'
93  call prc_abort
94  end select
95 
96  else
97  log_info("ATMOS_PHY_CP_driver_setup",*) 'this component is never called.'
98  endif
99 
100  return
101  end subroutine atmos_phy_cp_driver_setup
102 
103  !-----------------------------------------------------------------------------
108  use mod_atmos_admin, only: &
111  implicit none
112  !---------------------------------------------------------------------------
113 
114  log_newline
115  log_info("ATMOS_PHY_CP_driver_finalize",*) 'Finalize'
116 
117  if ( atmos_sw_phy_cp ) then
118 
119  select case ( atmos_phy_cp_type )
120  case ( 'KF' )
122  end select
123 
124  endif
125 
126  return
127  end subroutine atmos_phy_cp_driver_finalize
128 
129  !-----------------------------------------------------------------------------
131  subroutine atmos_phy_cp_driver_calc_tendency( update_flag )
132  use scale_statistics, only: &
134  statistics_total
135  use scale_atmos_grid_cartesc_real, only: &
144  use scale_file_history, only: &
145  file_history_in
146  use scale_time , only :&
147  time_dtsec, &
149  use scale_atmos_grid_cartesc_real, only: &
152  use scale_atmos_hydrometeor, only: &
153  hyd_name, &
154  cv_water, &
155  cv_ice, &
156  lhf
157  use scale_atmos_phy_cp_kf, only: &
159  use scale_atmos_phy_cp_kf_jmapplib, only: &
161  use scale_atmos_phy_cp_common, only: &
163  use mod_atmos_phy_mp_vars, only: &
164  qs_mp, &
165  qe_mp
166  use mod_atmos_phy_mp_driver, only: &
168  use mod_atmos_admin, only: &
170  use mod_atmos_vars, only: &
171  dens => dens_av, &
172  momz => momz_av, &
173  momx => momx_av, &
174  momy => momy_av, &
175  rhot => rhot_av, &
176  qtrc => qtrc_av, &
177  dens_t => dens_tp, &
178  momz_t => momz_tp, &
179  momx_t => momx_tp, &
180  momy_t => momy_tp, &
181  rhot_t => rhot_tp, &
182  rhoq_t => rhoq_tp, &
183  u, &
184  v, &
185  w, &
186  temp, &
187  pott, &
188  pres, &
189  exner, &
190  qdry, &
191  qv, &
192  qc, &
193  qi
194  use scale_atmos_hydrometeor, only: &
195  n_hyd
196  use mod_atmos_phy_cp_vars, only: &
197  dens_t_cp => atmos_phy_cp_dens_t, &
198  rhot_t_cp => atmos_phy_cp_rhot_t, &
199  rhoqv_t_cp => atmos_phy_cp_rhoqv_t, &
200  rhohyd_t_cp => atmos_phy_cp_rhohyd_t, &
201  mflx_cloudbase => atmos_phy_cp_mflx_cloudbase, &
202  sflx_rain => atmos_phy_cp_sflx_rain, & ! convective rain [kg/m2/s]
203  sflx_snow => atmos_phy_cp_sflx_snow, & ! convective snow [kg/m2/s]
204  sflx_engi => atmos_phy_cp_sflx_engi, & ! internal energy flux [J/m2/s]
205  cloudtop => atmos_phy_cp_cloudtop, & ! cloud top height [m]
206  cloudbase => atmos_phy_cp_cloudbase, & ! cloud base height [m]
207  cldfrac_dp => atmos_phy_cp_cldfrac_dp, & ! cloud fraction (deep convection) (0-1)
208  cldfrac_sh => atmos_phy_cp_cldfrac_sh, & ! cloud fraction (shallow convection) (0-1)
209  w0mean => atmos_phy_cp_w0mean, & ! running mean vertical wind velocity [m/s]
210  kf_nca => atmos_phy_cp_kf_nca ! advection/cumulus convection timescale/dt for KF [step]
211  use mod_atmos_phy_bl_vars, only: &
212  pblh => atmos_phy_bl_zi, &
213  sflx_buoy => atmos_phy_bl_sflx_buoy
214  use mod_atmos_phy_sf_vars, only: &
215  us => atmos_phy_sf_ustar
216  implicit none
217 
218  logical, intent(in) :: update_flag
219 
220  real(rp) :: rhoq_t_cp(ka,ia,ja,qs_mp:qe_mp)
221 
222  real(rp) :: sflx_prec(ia,ja)
223 
224  integer :: k, i, j, iq
225  !---------------------------------------------------------------------------
226 
227  if ( atmos_phy_cp_type /= "KF-JMAPPLIB" ) then
228 
229  ! temporal running mean of vertical velocity
230 ! !$acc update host(W,w0mean)
231  call atmos_phy_cp_common_wmean( ka, ks, ke, ia, is, ie, ja, js, je, &
232  w(:,:,:), & ! [IN]
234  w0mean(:,:,:) ) ! [INOUT]
235 ! !$acc update device(w0mean)
236  end if
237 
238 
239  if ( update_flag ) then ! update
240 
241  !$acc data create(SFLX_prec)
242 
243  select case ( atmos_phy_cp_type )
244  case ( 'KF' )
245 ! !$acc update host(DENS,U,V,RHOT,TEMP,PRES,QDRY,QV,DENS_t_CP,RHOT_t_CP,RHOQV_t_CP,RHOHYD_t_CP,kf_nca,SFLX_rain,SFLX_snow,SFLX_ENGI,cloudtop,cloudbase,cldfrac_dp,cldfrac_sh)
246  call atmos_phy_cp_kf_tendency( ka, ks, ke, ia, is, ie, ja, js, je, &
247  dens(:,:,:), & ! [IN]
248  u(:,:,:), v(:,:,:), & ! [IN]
249  rhot(:,:,:), temp(:,:,:), pres(:,:,:), & ! [IN]
250  qdry(:,:,:), qv(:,:,:), & ! [IN]
251  w0mean(:,:,:), & ! [IN]
252  fz(:,:,:), & ! [IN]
253  time_dtsec_atmos_phy_cp, & ! [IN]
254  dens_t_cp(:,:,:), & ! [INOUT]
255  rhot_t_cp(:,:,:), & ! [INOUT]
256  rhoqv_t_cp(:,:,:), rhohyd_t_cp(:,:,:,:), & ! [INOUT]
257  kf_nca(:,:), & ! [INOUT]
258  sflx_rain(:,:), sflx_snow(:,:), & ! [INOUT]
259  sflx_engi(:,:), & ! [INOUT]
260  cloudtop(:,:), cloudbase(:,:), & ! [INOUT]
261  cldfrac_dp(:,:,:), cldfrac_sh(:,:,:) ) ! [INOUT]
262 ! !$acc update device(DENS_t_CP,RHOT_t_CP,RHOQV_t_CP,RHOHYD_t_CP,kf_nca,SFLX_rain,SFLX_snow,SFLX_ENGI,cloudtop,cloudbase,cldfrac_dp,cldfrac_sh)
263 
264  !$omp parallel do
265  !$acc kernels
266  do j = js, je
267  do i = is, ie
268  sflx_prec(i,j) = sflx_rain(i,j) + sflx_snow(i,j)
269  end do
270  end do
271  !$acc end kernels
272 
273  case ( 'KF-JMAPPLIB' )
274 
275  !$omp parallel do
276  !$acc kernels
277  do j = js, je
278  do i = is, ie
279  sflx_prec(i,j) = sflx_rain(i,j) + sflx_snow(i,j)
280  end do
281  end do
282  !$acc end kernels
283 
284  !$acc update host(DENS,U,V,W,TEMP,POTT,PRES,EXNER,QDRY,QV,QC,QI,us,PBLH,SFLX_BUOY,RHOT_t_CP,RHOQV_t_CP,RHOHYD_t_CP,w0mean,kf_nca,SFLX_rain,SFLX_snow,SFLX_prec,cloudtop,cloudbase)
286  dens(:,:,:), & ! [IN]
287  u(:,:,:), v(:,:,:), w(:,:,:), & ! [IN]
288  temp(:,:,:), pott(:,:,:), & ! [IN]
289  pres(:,:,:), exner(:,:,:), & ! [IN]
290  qdry(:,:,:), qv(:,:,:), & ! [IN]
291  qc(:,:,:), qi(:,:,:), & ! [IN]
292  us(:,:), pblh(:,:), sflx_buoy(:,:), & ! [IN]
293  cz(:,:,:), fz(:,:,:), & ! [IN]
294  dens_t_cp(:,:,:), & ! [OUT]
295  rhot_t_cp(:,:,:), & ! [INOUT]
296  rhoqv_t_cp(:,:,:), rhohyd_t_cp(:,:,:,:), & ! [INOUT]
297  w0mean(:,:,:), kf_nca(:,:), & ! [INOUT]
298  sflx_rain(:,:), sflx_snow(:,:), & ! [INOUT]
299  sflx_prec(:,:), & ! [INOUT]
300  cloudtop(:,:), cloudbase(:,:) ) ! [INOUT]
301  !$acc update device(DENS_t_CP,RHOT_t_CP,RHOQV_t_CP,RHOHYD_t_CP,w0mean,kf_nca,SFLX_rain,SFLX_snow,SFLX_prec,cloudtop,cloudbase)
302  !$omp parallel do
303  !$acc kernels
304  do j = js, je
305  do i = is, ie
306  ! assume that temperature of precipitation is the same as that at the lowest layer
307  sflx_engi(i,j) = sflx_rain(i,j) * cv_water * temp(ks,i,j) &
308  + sflx_snow(i,j) * ( cv_ice * temp(ks,i,j) - lhf )
309  end do
310  end do
311  !$acc end kernels
312 
313  end select
314 
315 !OCL XFILL
316  !$omp parallel do
317  !$acc kernels
318  do j = js, je
319  do i = is, ie
320  mflx_cloudbase(i,j) = 0.0_rp
321  enddo
322  enddo
323  !$acc end kernels
324 
325  ! diagnose tendency of number concentration
326 
327  call file_history_in( w0mean(:,:,:), 'w0mean', 'running mean vertical wind velocity', 'kg/m2/s', fill_halo=.true. )
328  call file_history_in( mflx_cloudbase(:,:), 'CBMFX', 'cloud base mass flux', 'kg/m2/s', fill_halo=.true. )
329  call file_history_in( sflx_rain(:,:), 'RAIN_CP', 'surface rain rate by CP', 'kg/m2/s', fill_halo=.true. )
330  call file_history_in( sflx_snow(:,:), 'SNOW_CP', 'surface snow rate by CP', 'kg/m2/s', fill_halo=.true. )
331  call file_history_in( sflx_prec(:,:), 'PREC_CP', 'surface precipitation rate by CP', 'kg/m2/s', fill_halo=.true. )
332  call file_history_in( cloudtop(:,:), 'CUMHGT', 'CP cloud top height', 'm', fill_halo=.true. )
333  call file_history_in( cloudbase(:,:), 'CUBASE', 'CP cloud base height', 'm', fill_halo=.true. )
334  call file_history_in( cldfrac_dp(:,:,:), 'CUMFRC_DP', 'CP cloud fraction (deep)', '1', fill_halo=.true. )
335  call file_history_in( cldfrac_sh(:,:,:), 'CUMFRC_SH', 'CP cloud fraction (shallow)', '1', fill_halo=.true. )
336  call file_history_in( kf_nca(:,:), 'kf_nca', 'advection or cumulus convection timescale for KF', 's', fill_halo=.true. )
337 
338  call file_history_in( dens_t_cp(:,:,:), 'DENS_t_CP', 'tendency DENS in CP', 'kg/m3/s' , fill_halo=.true. )
339  call file_history_in( rhot_t_cp(:,:,:), 'RHOT_t_CP', 'tendency RHOT in CP', 'K*kg/m3/s', fill_halo=.true. )
340 
341  call file_history_in( rhoqv_t_cp(:,:,:), 'QV_t_CP', 'tendency rho*QV in CP', 'kg/m3/s', fill_halo=.true. )
342  do iq = 1, n_hyd
343  call file_history_in( rhohyd_t_cp(:,:,:,iq), trim(hyd_name(iq))//'_t_CP', &
344  'tendency rho*'//trim(hyd_name(iq))//' in CP', 'kg/m3/s', fill_halo=.true. )
345  enddo
346 
347  !$acc end data
348 
349  endif ! update
350 
351  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
352  !$acc kernels
353  do j = js, je
354  do i = is, ie
355  do k = ks, ke
356  dens_t(k,i,j) = dens_t(k,i,j) + dens_t_cp(k,i,j)
357  rhot_t(k,i,j) = rhot_t(k,i,j) + rhot_t_cp(k,i,j)
358  enddo
359  enddo
360  enddo
361  !$acc end kernels
362 
363  !$acc data create(RHOQ_t_CP)
364 
365  call atmos_phy_mp_driver_qhyd2qtrc( ka, ks, ke, ia, is, ie, ja, js, je, &
366  rhoqv_t_cp(:,:,:), rhohyd_t_cp(:,:,:,:), & ! [IN]
367  rhoq_t_cp(:,:,:,qs_mp:qe_mp) ) ! [OUT]
368 
369  !$omp parallel do private(iq,i,j,k) OMP_SCHEDULE_ collapse(3)
370  !$acc kernels
371  do iq = qs_mp, qe_mp
372  do j = js, je
373  do i = is, ie
374  do k = ks, ke
375  rhoq_t(k,i,j,iq) = rhoq_t(k,i,j,iq) + rhoq_t_cp(k,i,j,iq)
376  enddo
377  enddo
378  enddo
379  enddo
380  !$acc end kernels
381 
382  if ( statistics_checktotal ) then
383  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
384  dens_t_cp(:,:,:), 'DENS_t_CP', &
385  atmos_grid_cartesc_real_vol(:,:,:), &
386  atmos_grid_cartesc_real_totvol )
387  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
388  rhot_t_cp(:,:,:), 'RHOT_t_CP', &
389  atmos_grid_cartesc_real_vol(:,:,:), &
390  atmos_grid_cartesc_real_totvol )
391 
392  do iq = qs_mp, qe_mp
393  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
394  rhoq_t_cp(:,:,:,iq), trim(tracer_name(iq))//'_t_CP', &
395  atmos_grid_cartesc_real_vol(:,:,:), &
396  atmos_grid_cartesc_real_totvol )
397  enddo
398  endif
399 
400  !$acc end data
401 
402  return
403  end subroutine atmos_phy_cp_driver_calc_tendency
404 
405 end module mod_atmos_phy_cp_driver
mod_atmos_vars::momz_av
real(rp), dimension(:,:,:), pointer, public momz_av
Definition: mod_atmos_vars.F90:91
mod_atmos_phy_cp_vars::atmos_phy_cp_mflx_cloudbase
real(rp), dimension(:,:), allocatable, public atmos_phy_cp_mflx_cloudbase
Definition: mod_atmos_phy_cp_vars.F90:73
scale_statistics
module Statistics
Definition: scale_statistics.F90:11
mod_atmos_phy_mp_vars
module Atmosphere / Physics Cloud Microphysics
Definition: mod_atmos_phy_mp_vars.F90:12
mod_atmos_phy_cp_vars::atmos_phy_cp_rhoqv_t
real(rp), dimension(:,:,:), allocatable, public atmos_phy_cp_rhoqv_t
Definition: mod_atmos_phy_cp_vars.F90:61
scale_atmos_grid_cartesc_index::ke
integer, public ke
end point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:52
mod_atmos_phy_cp_vars::atmos_phy_cp_rhohyd_t
real(rp), dimension(:,:,:,:), allocatable, public atmos_phy_cp_rhohyd_t
Definition: mod_atmos_phy_cp_vars.F90:62
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_atmos_phy_cp_kf::atmos_phy_cp_kf_setup
subroutine, public atmos_phy_cp_kf_setup(KA, KS, KE, IA, IS, IE, JA, JS, JE, CZ, AREA, WARMRAIN_in)
Setup initial setup for Kain-Fritsch Cumulus Parameterization.
Definition: scale_atmos_phy_cp_kf.F90:190
mod_atmos_vars::momx_av
real(rp), dimension(:,:,:), pointer, public momx_av
Definition: mod_atmos_vars.F90:92
scale_atmos_phy_cp_kf_jmapplib
module atmosphere / physics / cumulus / kf-jmapplib
Definition: scale_atmos_phy_cp_kf_jmapplib.F90:13
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvolzxv
real(rp), public atmos_grid_cartesc_real_totvolzxv
total volume (zxv, local) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:91
mod_atmos_vars::rhoq_tp
real(rp), dimension(:,:,:,:), allocatable, public rhoq_tp
Definition: mod_atmos_vars.F90:121
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_cz
geopotential height [m] (zxy)
Definition: scale_atmos_grid_cartesC_real.F90:39
mod_atmos_vars::pott
real(rp), dimension(:,:,:), allocatable, target, public pott
Definition: mod_atmos_vars.F90:133
scale_atmos_hydrometeor::hyd_name
character(len=h_short), dimension(n_hyd), parameter, public hyd_name
Definition: scale_atmos_hydrometeor.F90:104
mod_atmos_phy_cp_vars::atmos_phy_cp_rhot_t
real(rp), dimension(:,:,:), allocatable, public atmos_phy_cp_rhot_t
Definition: mod_atmos_phy_cp_vars.F90:60
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_volwxy
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_volwxy
control volume (wxy) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:85
mod_atmos_vars::qtrc_av
real(rp), dimension(:,:,:,:), pointer, public qtrc_av
Definition: mod_atmos_vars.F90:95
scale_atmos_phy_cp_common
module atmosphere / physics / cumulus / Common
Definition: scale_atmos_phy_cp_common.F90:13
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_atmos_grid_cartesc_index::ka
integer, public ka
Definition: scale_atmos_grid_cartesC_index.F90:47
mod_atmos_phy_bl_vars::atmos_phy_bl_zi
real(rp), dimension(:,:), allocatable, public atmos_phy_bl_zi
Definition: mod_atmos_phy_bl_vars.F90:67
mod_atmos_admin
module ATMOS admin
Definition: mod_atmos_admin.F90:11
mod_atmos_vars::qdry
real(rp), dimension(:,:,:), allocatable, target, public qdry
Definition: mod_atmos_vars.F90:140
mod_atmos_phy_mp_vars::qs_mp
integer, public qs_mp
Definition: mod_atmos_phy_mp_vars.F90:79
mod_atmos_phy_sf_vars::atmos_phy_sf_ustar
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_ustar
Definition: mod_atmos_phy_sf_vars.F90:90
mod_atmos_phy_bl_vars::atmos_phy_bl_sflx_buoy
real(rp), dimension(:,:), allocatable, public atmos_phy_bl_sflx_buoy
Definition: mod_atmos_phy_bl_vars.F90:68
scale_atmos_phy_cp_kf::atmos_phy_cp_kf_tendency
subroutine, public atmos_phy_cp_kf_tendency(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, U, V, RHOT, TEMP, PRES, QDRY, QV_in, w0avg, FZ, KF_DTSECD, DENS_t, RHOT_t, RHOQV_t, RHOQ_t, nca, SFLX_rain, SFLX_snow, SFLX_engi, cloudtop, cloudbase, cldfrac_dp, cldfrac_sh)
ATMOS_PHY_CP_kf calculate Kain-Fritsch Cumulus Parameterization.
Definition: scale_atmos_phy_cp_kf.F90:578
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvolwxy
real(rp), public atmos_grid_cartesc_real_totvolwxy
total volume (wxy, local) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:89
mod_atmos_vars::rhot_av
real(rp), dimension(:,:,:), pointer, public rhot_av
Definition: mod_atmos_vars.F90:94
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
mod_atmos_phy_sf_vars
module ATMOSPHERIC Surface Variables
Definition: mod_atmos_phy_sf_vars.F90:12
scale_atmos_phy_cp_kf
module atmosphere / physics / cumulus / Kain-Fritsch
Definition: scale_atmos_phy_cp_kf.F90:45
mod_atmos_vars::rhot
real(rp), dimension(:,:,:), allocatable, target, public rhot
Definition: mod_atmos_vars.F90:80
scale_atmos_grid_cartesc_real
module Atmosphere GRID CartesC Real(real space)
Definition: scale_atmos_grid_cartesC_real.F90:11
mod_atmos_vars::qtrc
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc
Definition: mod_atmos_vars.F90:81
scale_file_history
module file_history
Definition: scale_file_history.F90:15
scale_atmos_grid_cartesc::dx
real(rp), public dx
Definition: scale_atmos_grid_cartesC.F90:39
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_volzuy
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_volzuy
control volume (zuy) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:86
mod_atmos_phy_mp_driver
module atmosphere / physics / cloud microphysics
Definition: mod_atmos_phy_mp_driver.F90:12
scale_prc
module PROCESS
Definition: scale_prc.F90:11
mod_atmos_phy_cp_vars::atmos_phy_cp_sflx_rain
real(rp), dimension(:,:), allocatable, public atmos_phy_cp_sflx_rain
Definition: mod_atmos_phy_cp_vars.F90:63
scale_time::time_dtsec_atmos_phy_cp
real(dp), public time_dtsec_atmos_phy_cp
time interval of physics(cumulus ) [sec]
Definition: scale_time.F90:37
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_atmos_grid_cartesc_index::ie
integer, public ie
end point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:54
scale_io
module STDIO
Definition: scale_io.F90:10
scale_atmos_phy_cp_kf_jmapplib::atmos_phy_cp_kf_jmapplib_tendency
subroutine, public atmos_phy_cp_kf_jmapplib_tendency(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, U, V, W, TEMP, POTT, PRES, EXNER, QDRY, QV, QC, QI, us, PBLH, SFLX_BUOY, CZ, FZ, DENS_t, RHOT_t, RHOQV_t, RHOQ_t, w0avg, nca, SFLX_rain, SFLX_snow, SFLX_prec, cloudtop, cloudbase)
ATMOS_PHY_CP_KF_JMAPPLIB_tendency calculate tendency by the virtical eddy viscosity.
Definition: scale_atmos_phy_cp_kf_jmapplib.F90:217
mod_atmos_phy_bl_vars
module atmosphere / physics / PBL
Definition: mod_atmos_phy_bl_vars.F90:12
mod_atmos_vars::dens
real(rp), dimension(:,:,:), allocatable, target, public dens
Definition: mod_atmos_vars.F90:76
scale_tracer::k
real(rp), public k
Definition: scale_tracer.F90:45
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_volzxv
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_volzxv
control volume (zxv) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:87
scale_atmos_grid_cartesc_index
module atmosphere / grid / cartesC index
Definition: scale_atmos_grid_cartesC_index.F90:12
scale_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_vol
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_vol
control volume (zxy) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:84
mod_atmos_vars::momz
real(rp), dimension(:,:,:), allocatable, target, public momz
Definition: mod_atmos_vars.F90:77
mod_atmos_vars::momy_av
real(rp), dimension(:,:,:), pointer, public momy_av
Definition: mod_atmos_vars.F90:93
mod_atmos_vars::v
real(rp), dimension(:,:,:), allocatable, target, public v
Definition: mod_atmos_vars.F90:131
scale_tracer::tracer_name
character(len=h_short), dimension(qa_max), public tracer_name
Definition: scale_tracer.F90:39
mod_atmos_vars::w
real(rp), dimension(:,:,:), allocatable, target, public w
Definition: mod_atmos_vars.F90:129
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_area
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_area
horizontal area ( xy, normal z) [m2]
Definition: scale_atmos_grid_cartesC_real.F90:66
mod_atmos_vars::momz_tp
real(rp), dimension(:,:,:), allocatable, public momz_tp
Definition: mod_atmos_vars.F90:116
mod_atmos_vars::momx
real(rp), dimension(:,:,:), allocatable, target, public momx
Definition: mod_atmos_vars.F90:78
mod_atmos_vars::exner
real(rp), dimension(:,:,:), allocatable, target, public exner
Definition: mod_atmos_vars.F90:136
scale_atmos_grid_cartesc_index::is
integer, public is
start point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:53
mod_atmos_admin::atmos_sw_phy_cp
logical, public atmos_sw_phy_cp
Definition: mod_atmos_admin.F90:59
scale_atmos_phy_cp_kf::atmos_phy_cp_kf_finalize
subroutine, public atmos_phy_cp_kf_finalize
finalize
Definition: scale_atmos_phy_cp_kf.F90:435
mod_atmos_vars::temp
real(rp), dimension(:,:,:), allocatable, target, public temp
Definition: mod_atmos_vars.F90:134
mod_atmos_vars::dens_tp
real(rp), dimension(:,:,:), allocatable, public dens_tp
Definition: mod_atmos_vars.F90:115
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvolzuy
real(rp), public atmos_grid_cartesc_real_totvolzuy
total volume (zuy, local) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:90
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
mod_atmos_vars::momy
real(rp), dimension(:,:,:), allocatable, target, public momy
Definition: mod_atmos_vars.F90:79
scale_time
module TIME
Definition: scale_time.F90:11
mod_atmos_vars::qv
real(rp), dimension(:,:,:), allocatable, pointer, target, public qv
Definition: mod_atmos_vars.F90:97
mod_atmos_phy_cp_driver
module ATMOSPHERE / Physics Cumulus
Definition: mod_atmos_phy_cp_driver.F90:12
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
mod_atmos_admin::atmos_phy_cp_type
character(len=h_short), public atmos_phy_cp_type
Definition: mod_atmos_admin.F90:43
mod_atmos_vars::pres
real(rp), dimension(:,:,:), allocatable, target, public pres
Definition: mod_atmos_vars.F90:135
mod_atmos_vars::dens_av
real(rp), dimension(:,:,:), pointer, public dens_av
Definition: mod_atmos_vars.F90:90
mod_atmos_vars::u
real(rp), dimension(:,:,:), allocatable, target, public u
Definition: mod_atmos_vars.F90:130
scale_atmos_grid_cartesc_index::ks
integer, public ks
start point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:51
scale_atmos_hydrometeor::atmos_hydrometeor_ice_phase
logical, public atmos_hydrometeor_ice_phase
Definition: scale_atmos_hydrometeor.F90:90
scale_atmos_phy_cp_common::atmos_phy_cp_common_wmean
subroutine, public atmos_phy_cp_common_wmean(KA, KS, KE, IA, IS, IE, JA, JS, JE, W, TIME_DTSEC, CP_DTSEC, W0_mean)
ATMOS_PHY_CP_wmean running mean vertical wind velocity comment for W0 imported from WRF.
Definition: scale_atmos_phy_cp_common.F90:102
mod_atmos_phy_cp_driver::atmos_phy_cp_driver_setup
subroutine, public atmos_phy_cp_driver_setup
Setup.
Definition: mod_atmos_phy_cp_driver.F90:50
scale_statistics::statistics_checktotal
logical, public statistics_checktotal
calc&report variable totals to logfile?
Definition: scale_statistics.F90:109
scale_time::time_dtsec
real(dp), public time_dtsec
time interval of model [sec]
Definition: scale_time.F90:33
mod_atmos_vars
module ATMOSPHERIC Variables
Definition: mod_atmos_vars.F90:12
mod_atmos_vars::qi
real(rp), dimension(:,:,:), pointer, public qi
Definition: mod_atmos_vars.F90:100
mod_atmos_vars::rhot_tp
real(rp), dimension(:,:,:), allocatable, public rhot_tp
Definition: mod_atmos_vars.F90:119
scale_atmos_grid_cartesc_index::js
integer, public js
start point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:55
mod_atmos_phy_cp_driver::atmos_phy_cp_driver_finalize
subroutine, public atmos_phy_cp_driver_finalize
finalize
Definition: mod_atmos_phy_cp_driver.F90:106
mod_atmos_phy_cp_vars::atmos_phy_cp_dens_t
real(rp), dimension(:,:,:), allocatable, public atmos_phy_cp_dens_t
Definition: mod_atmos_phy_cp_vars.F90:58
mod_atmos_vars::momy_tp
real(rp), dimension(:,:,:), allocatable, public momy_tp
Definition: mod_atmos_vars.F90:125
mod_atmos_phy_cp_driver::atmos_phy_cp_driver_calc_tendency
subroutine, public atmos_phy_cp_driver_calc_tendency(update_flag)
Driver.
Definition: mod_atmos_phy_cp_driver.F90:132
scale_atmos_hydrometeor::lhf
real(rp), public lhf
latent heat of fusion for use [J/kg]
Definition: scale_atmos_hydrometeor.F90:146
mod_atmos_vars::momx_tp
real(rp), dimension(:,:,:), allocatable, public momx_tp
Definition: mod_atmos_vars.F90:124
scale_atmos_phy_cp_kf_jmapplib::atmos_phy_cp_kf_jmapplib_setup
subroutine, public atmos_phy_cp_kf_jmapplib_setup(KA, KS, KE, IA, JA, dx, dt)
ATMOS_PHY_CP_KF_JMAPPLIB_setup Setup.
Definition: scale_atmos_phy_cp_kf_jmapplib.F90:62
mod_atmos_vars::qc
real(rp), dimension(:,:,:), pointer, public qc
Definition: mod_atmos_vars.F90:98
mod_atmos_phy_mp_vars::qe_mp
integer, public qe_mp
Definition: mod_atmos_phy_mp_vars.F90:80
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_fz
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_fz
geopotential height [m] (wxy)
Definition: scale_atmos_grid_cartesC_real.F90:43
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvol
real(rp), public atmos_grid_cartesc_real_totvol
total volume (zxy, local) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:88
scale_atmos_grid_cartesc
module atmosphere / grid / cartesC
Definition: scale_atmos_grid_cartesC.F90:12
scale_atmos_phy_cp_common::atmos_phy_cp_common_setup
subroutine, public atmos_phy_cp_common_setup
Setup.
Definition: scale_atmos_phy_cp_common.F90:52
mod_atmos_phy_mp_driver::atmos_phy_mp_driver_qhyd2qtrc
subroutine, public atmos_phy_mp_driver_qhyd2qtrc(KA, KS, KE, IA, IS, IE, JA, JS, JE, QV, QHYD, QTRC, QNUM)
Definition: mod_atmos_phy_mp_driver.F90:1553
mod_atmos_phy_cp_vars
module Atmosphere / Physics Cumulus
Definition: mod_atmos_phy_cp_vars.F90:12
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56
scale_atmos_hydrometeor::n_hyd
integer, parameter, public n_hyd
Definition: scale_atmos_hydrometeor.F90:95
scale_atmos_hydrometeor::cv_water
real(rp), public cv_water
CV for water [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:151
scale_atmos_hydrometeor::cv_ice
real(rp), public cv_ice
CV for ice [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:153