SCALE-RM
Data Types | Functions/Subroutines
scale_atmos_phy_cp_kf Module Reference

module atmosphere / physics / cumulus / Kain-Fritsch More...

Functions/Subroutines

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. More...
 
subroutine, public atmos_phy_cp_kf_finalize
 finalize More...
 
subroutine, public atmos_phy_cp_kf_putvar (in_delcape, in_deeplifetime, in_shallowlifetime)
 overwrite private variables More...
 
subroutine, public atmos_phy_cp_kf_getvar (out_delcape, out_deeplifetime, out_shallowlifetime)
 read private variables More...
 
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. More...
 

Detailed Description

module atmosphere / physics / cumulus / Kain-Fritsch

Description
Main module of Kain-Fritsch Cumulus Convection Parameterization
Author
Team SCALE
References
  • Kain J. S. and J. M. Fritsch, 1990: A one-dimensional entraining/detraining plume model and its application in convective parameterization. J. Atmos. Sci. 47:2784-2802.
  • Kain J. S., 2004: The Kain-Fritsch Convective Parameterization: An Update, J. Appl. Meteor., 43, 170-181.
  • Narita, M. and S. Ohmori, 2007: Improving Precipitation Forecasts by the Operational Nonhydrostatic Mesoscale Model with the Kain-Fritsch Convective Parameterization and Cloud Microphysics. 12th Conference on Mesoscale Processes, 6-9 August 2007, Waterville Valley conference & event center, Waterville Valley, NH, USA.
  • Ogura, Y., and H. R. Cho, 1973: Diagnostic determination of cumulus cloud populations from observed large-scale variables. J. Atmos. Sci., 30, 1276-1286.
This file was originally copied from WRF. The original file was published with the following notice.

WRF was developed at the National Center for Atmospheric Research (NCAR) which is operated by the University Corporation for Atmospheric Research (UCAR). NCAR and UCAR make no proprietary claims, either statutory or otherwise, to this version and release of WRF and consider WRF to be in the public domain for use by any person or entity for any purpose without any fee or charge. UCAR requests that any WRF user include this notice on any partial or full copies of WRF. WRF is provided on an "AS IS" basis and any warranties, either express or implied, including but not limited to implied warranties of non-infringement, originality, merchantability and fitness for a particular purpose, are disclaimed. In no event shall UCAR be liable for any damages, whatsoever, whether direct, indirect, consequential or special, that arise out of or in connection with the access, use or performance of WRF, including infringement actions.

WRF is a registered trademark of the University Corporation for Atmospheric Research (UCAR).

NAMELIST
  • PARAM_ATMOS_PHY_CP_KF
    nametypedefault valuecomment
    ATMOS_PHY_CP_KF_RATE real(RP) 0.03_RP ratio of cloud water and precipitation (Ogura and Cho 1973)
    ATMOS_PHY_CP_KF_TRIGGER_TYPE integer 1 trigger function type 1:KF2004 3:NO2007
    ATMOS_PHY_CP_KF_DLCAPE real(RP) 0.1_RP cape decleace rate
    ATMOS_PHY_CP_KF_DLIFETIME real(RP) 1800.0_RP minimum lifetime scale of deep convection [sec]
    ATMOS_PHY_CP_KF_SLIFETIME real(RP) 2400.0_RP lifetime of shallow convection [sec]
    ATMOS_PHY_CP_KF_DEPTH_USL real(RP) 300.0_RP depth of updraft source layer [hPa]
    ATMOS_PHY_CP_KF_PREC_TYPE integer 1 precipitation type 1:Ogura-Cho(1973) 2:Kessler
    ATMOS_PHY_CP_KF_DO_PREC logical .true. whether enable precipitation
    ATMOS_PHY_CP_KF_THRES real(RP) 1.E-3_RP autoconversion rate for Kessler
    ATMOS_PHY_CP_KF_LOG logical .false. output log?

History Output
namedescriptionunitvariable
KF_CONVFLAG CONVECTION FLAG R_convflag
KF_LIFETIME lifetime of KF scheme s lifetime
QC_t_KF_conv_DP QC tendency due to flux convergence (deep convection) in KF kg/m3/s QC_t_KF_conv_DP
QC_t_KF_conv_SH QC tendency due to flux convergence (shallow convection) in KF kg/m3/s QC_t_KF_conv_SH
QC_t_KF_updet_DP QC tendency due to detrainment of updraft (deep convection) in KF kg/m3/s QC_t_KF_updet_DP
QC_t_KF_updet_SH QC tendency due to detrainment of updraft (shallow convection) in KF kg/m3/s QC_t_KF_updet_SH
QI_t_KF_conv_DP QI tendency due to flux convergence (deep convection) in KF kg/m3/s QI_t_KF_conv_DP
QI_t_KF_conv_SH QI tendency due to flux convergence (shallow convection) in KF kg/m3/s QI_t_KF_conv_SH
QI_t_KF_updet_DP QI tendency due to detrainment of updraft (deep convection) in KF kg/m3/s QI_t_KF_updet_DP
QI_t_KF_updet_SH QI tendency due to detrainment of updraft (shallow convection) in KF kg/m3/s QI_t_KF_updet_SH
QR_t_KF_conv_DP QR tendency due to flux convergence (deep convection) in KF kg/m3/s QR_t_KF_conv_DP
QR_t_KF_conv_SH QR tendency due to flux convergence (shallow convection) in KF kg/m3/s QR_t_KF_conv_SH
QR_t_KF_rainfall_DP QR tendency due to rain fall (deep convection) in KF kg/m3/s QR_t_KF_rainfall_DP
QR_t_KF_rainfall_SH QR tendency due to rain fall (shallow convection) in KF kg/m3/s QR_t_KF_rainfall_SH
QS_t_KF_conv_DP QS tendency due to flux convergence (deep convection) in KF kg/m3/s QS_t_KF_conv_DP
QS_t_KF_conv_SH QS tendency due to flux convergence (shallow convection) in KF kg/m3/s QS_t_KF_conv_SH
QS_t_KF_snowfall_DP QS tendency due to snow fall (deep convection) in KF kg/m3/s QS_t_KF_snowfall_DP
QS_t_KF_snowfall_SH QS tendency due to snow fall (shallow convection) in KF kg/m3/s QS_t_KF_snowfall_SH
QV_t_KF_conv_DP QV tendency due to flux convergence (deep convection) in KF kg/m3/s QV_t_KF_conv_DP
QV_t_KF_conv_SH QV tendency due to flux convergence (shallow convection) in KF kg/m3/s QV_t_KF_conv_SH
QV_t_KF_downdet_DP QV tendency due to detrainment of downdraft (deep convection) in KF kg/m3/s QV_t_KF_downdet_DP
QV_t_KF_downent_DP QV tendency due to entrainment of downdraft (deep convection) in KF kg/m3/s QV_t_KF_downent_DP
QV_t_KF_negfix_DP QV tendency due to negative fixer (deep convection) in KF kg/m3/s QV_t_KF_negfix_DP
QV_t_KF_negfix_SH QV tendency due to negative fixer (shallow convection) in KF kg/m3/s QV_t_KF_negfix_SH
QV_t_KF_updet_DP QV tendency due to detrainment of updraft (deep convection) in KF kg/m3/s QV_t_KF_updet_DP
QV_t_KF_updet_SH QV tendency due to detrainment of updraft (shallow convection) in KF kg/m3/s QV_t_KF_updet_SH
QV_t_KF_upent_DP QV tendency due to entrainment of updraft (deep convection) in KF kg/m3/s QV_t_KF_upent_DP
QV_t_KF_upent_SH QV tendency due to entrainment of updraft (shallow convection) in KF kg/m3/s QV_t_KF_upent_SH

Function/Subroutine Documentation

◆ atmos_phy_cp_kf_setup()

subroutine, public scale_atmos_phy_cp_kf::atmos_phy_cp_kf_setup ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension(ka,ia,ja), intent(in)  CZ,
real(rp), dimension(ia,ja), intent(in)  AREA,
logical, intent(in)  WARMRAIN_in 
)

Setup initial setup for Kain-Fritsch Cumulus Parameterization.

Parameters
[in]warmrain_intunning parameters, original parameter set is from KF2004 and NO2007

Definition at line 190 of file scale_atmos_phy_cp_kf.F90.

190  use scale_prc, only: &
191  prc_abort
192  use scale_file_history, only: &
194  implicit none
195  integer, intent(in) :: KA, KS, KE
196  integer, intent(in) :: IA, IS, IE
197  integer, intent(in) :: JA, JS, JE
198 
199  real(RP), intent(in) :: CZ(KA,IA,JA)
200  real(RP), intent(in) :: AREA(IA,JA)
201  logical, intent(in) :: WARMRAIN_in
202 
203 
204  integer :: ATMOS_PHY_CP_kf_trigger_type = 1
205  integer :: ATMOS_PHY_CP_kf_prec_type = 1
206  real(RP) :: ATMOS_PHY_CP_kf_rate = 0.03_rp
207  logical :: ATMOS_PHY_CP_kf_do_prec = .true.
208  real(RP) :: ATMOS_PHY_CP_kf_dlcape = 0.1_rp
209  real(RP) :: ATMOS_PHY_CP_kf_dlifetime = 1800.0_rp
210  real(RP) :: ATMOS_PHY_CP_kf_slifetime = 2400.0_rp
211  real(RP) :: ATMOS_PHY_CP_kf_DEPTH_USL = 300.0_rp
212  real(RP) :: ATMOS_PHY_CP_kf_thres = 1.e-3_rp
213  logical :: ATMOS_PHY_CP_kf_LOG = .false.
214 
215  namelist / param_atmos_phy_cp_kf / &
216  atmos_phy_cp_kf_rate, &
217  atmos_phy_cp_kf_trigger_type, &
218  atmos_phy_cp_kf_dlcape, &
219  atmos_phy_cp_kf_dlifetime, &
220  atmos_phy_cp_kf_slifetime, &
221  atmos_phy_cp_kf_depth_usl, &
222  atmos_phy_cp_kf_prec_type, &
223  atmos_phy_cp_kf_do_prec, &
224  atmos_phy_cp_kf_thres, &
225  atmos_phy_cp_kf_log
226 
227  integer :: k, i, j
228  integer :: n
229  integer :: ierr
230  !---------------------------------------------------------------------------
231 
232  log_newline
233  log_info("ATMOS_PHY_CP_kf_setup",*) 'Setup'
234  log_info("ATMOS_PHY_CP_kf_setup",*) 'Kain-Fritsch scheme'
235 
236  warmrain = warmrain_in
237 
238  !--- read namelist
239  rewind(io_fid_conf)
240  read(io_fid_conf,nml=param_atmos_phy_cp_kf,iostat=ierr)
241  if( ierr < 0 ) then !--- missing
242  log_info("ATMOS_PHY_CP_kf_setup",*) 'Not found namelist. Default used.'
243  elseif( ierr > 0 ) then !--- fatal error
244  log_error("ATMOS_PHY_CP_kf_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_CP_KF. Check!'
245  call prc_abort
246  endif
247  log_nml(param_atmos_phy_cp_kf)
248 
249  call cp_kf_lutab ! set up of KF look-up table
250 
251  call cp_kf_param( atmos_phy_cp_kf_prec_type, & ! [IN]
252  atmos_phy_cp_kf_rate, & ! [IN]
253  atmos_phy_cp_kf_do_prec, & ! [IN]
254  atmos_phy_cp_kf_dlcape, & ! [IN]
255  atmos_phy_cp_kf_dlifetime, & ! [IN]
256  atmos_phy_cp_kf_slifetime, & ! [IN]
257  atmos_phy_cp_kf_depth_usl, & ! [IN]
258  atmos_phy_cp_kf_thres, & ! [IN]
259  atmos_phy_cp_kf_log, & ! [IN]
260  atmos_phy_cp_kf_trigger_type ) ! [IN]
261 
262  ! output parameter lists
263  log_newline
264  log_info("ATMOS_PHY_CP_kf_setup",*) "Ogura-Cho condense material convert rate : ", atmos_phy_cp_kf_rate
265  log_info("ATMOS_PHY_CP_kf_setup",*) "Trigger function type, 1:KF2004 3:NO2007 : ", atmos_phy_cp_kf_trigger_type
266  log_info("ATMOS_PHY_CP_kf_setup",*) "CAPE decrease rate : ", atmos_phy_cp_kf_dlcape
267  log_info("ATMOS_PHY_CP_kf_setup",*) "Minimum lifetime scale of deep convection [sec] : ", atmos_phy_cp_kf_dlifetime
268  log_info("ATMOS_PHY_CP_kf_setup",*) "Lifetime of shallow convection [sec] : ", atmos_phy_cp_kf_slifetime
269  log_info("ATMOS_PHY_CP_kf_setup",*) "Updraft souce layer depth [hPa] : ", atmos_phy_cp_kf_depth_usl
270  log_info("ATMOS_PHY_CP_kf_setup",*) "Precipitation type 1:Ogura-Cho(1973) 2:Kessler : ", atmos_phy_cp_kf_prec_type
271  log_info("ATMOS_PHY_CP_kf_setup",*) "Kessler type precipitation's threshold : ", atmos_phy_cp_kf_thres
272  log_info("ATMOS_PHY_CP_kf_setup",*) "Enable sedimentation (precipitation) ? : ", atmos_phy_cp_kf_do_prec
273  log_info("ATMOS_PHY_CP_kf_setup",*) "Warm rain? : ", warmrain
274  log_info("ATMOS_PHY_CP_kf_setup",*) "Output log? : ", atmos_phy_cp_kf_log
275 
276  ! output variables
277  allocate( lifetime(ia,ja) )
278  allocate( i_convflag(ia,ja) )
279  lifetime(:,:) = 0.0_rp
280  i_convflag(:,:) = 2
281 
282  allocate( z(ka,ia,ja) )
283  z(:,:,:) = cz(:,:,:) ! because scale_atmos_phy_cp interface ,not use scale_grid
284 
285  allocate( deltaz(ka,ia,ja) )
286  ! deltaz is the interval of between model full levels(scalar point )
287  deltaz(:,:,:) = 0.0_rp
288  do j = js, je
289  do i = is, ie
290  do k = ks, ke-1
291  deltaz(k,i,j) = cz(k+1,i,j) - cz(k,i,j)
292  enddo
293  deltaz(ke,i,j) = 0.0_rp
294  enddo
295  enddo
296 
297  allocate( deltax(ia,ja) )
298  do j = js, je
299  do i = is, ie
300  deltax(i,j) = sqrt( area(i,j) )
301  end do
302  end do
303 
304 
305  ! history
306  hist_id(:,:) = -1
307  call file_history_reg( 'QV_t_KF_conv_DP', &
308  'QV tendency due to flux convergence (deep convection) in KF', &
309  'kg/m3/s', &
310  hist_id(i_hist_qv_fc,i_hist_d) )
311  call file_history_reg( 'QV_t_KF_conv_SH', &
312  'QV tendency due to flux convergence (shallow convection) in KF', &
313  'kg/m3/s', &
314  hist_id(i_hist_qv_fc,i_hist_s) )
315  call file_history_reg( 'QV_t_KF_upent_DP', &
316  'QV tendency due to entrainment of updraft (deep convection) in KF', &
317  'kg/m3/s', &
318  hist_id(i_hist_qv_ue,i_hist_d) )
319  call file_history_reg( 'QV_t_KF_upent_SH', &
320  'QV tendency due to entrainment of updraft (shallow convection) in KF', &
321  'kg/m3/s', &
322  hist_id(i_hist_qv_ue,i_hist_s) )
323  call file_history_reg( 'QV_t_KF_downent_DP', &
324  'QV tendency due to entrainment of downdraft (deep convection) in KF', &
325  'kg/m3/s', &
326  hist_id(i_hist_qv_de,i_hist_d) )
327  call file_history_reg( 'QV_t_KF_updet_DP', &
328  'QV tendency due to detrainment of updraft (deep convection) in KF', &
329  'kg/m3/s', &
330  hist_id(i_hist_qv_ud,i_hist_d) )
331  call file_history_reg( 'QV_t_KF_updet_SH', &
332  'QV tendency due to detrainment of updraft (shallow convection) in KF', &
333  'kg/m3/s', &
334  hist_id(i_hist_qv_ud,i_hist_s) )
335  call file_history_reg( 'QV_t_KF_downdet_DP', &
336  'QV tendency due to detrainment of downdraft (deep convection) in KF', &
337  'kg/m3/s', &
338  hist_id(i_hist_qv_dd,i_hist_d) )
339  call file_history_reg( 'QV_t_KF_negfix_DP', &
340  'QV tendency due to negative fixer (deep convection) in KF', &
341  'kg/m3/s', &
342  hist_id(i_hist_qv_nf,i_hist_d) )
343  call file_history_reg( 'QV_t_KF_negfix_SH', &
344  'QV tendency due to negative fixer (shallow convection) in KF', &
345  'kg/m3/s', &
346  hist_id(i_hist_qv_nf,i_hist_s) )
347 
348  call file_history_reg( 'QC_t_KF_conv_DP', &
349  'QC tendency due to flux convergence (deep convection) in KF', &
350  'kg/m3/s', &
351  hist_id(i_hist_qc_fc,i_hist_d) )
352  call file_history_reg( 'QC_t_KF_conv_SH', &
353  'QC tendency due to flux convergence (shallow convection) in KF', &
354  'kg/m3/s', &
355  hist_id(i_hist_qc_fc,i_hist_s) )
356  call file_history_reg( 'QC_t_KF_updet_DP', &
357  'QC tendency due to detrainment of updraft (deep convection) in KF', &
358  'kg/m3/s', &
359  hist_id(i_hist_qc_ud,i_hist_d) )
360  call file_history_reg( 'QC_t_KF_updet_SH', &
361  'QC tendency due to detrainment of updraft (shallow convection) in KF', &
362  'kg/m3/s', &
363  hist_id(i_hist_qc_ud,i_hist_s) )
364 
365  call file_history_reg( 'QI_t_KF_conv_DP', &
366  'QI tendency due to flux convergence (deep convection) in KF', &
367  'kg/m3/s', &
368  hist_id(i_hist_qi_fc,i_hist_d) )
369  call file_history_reg( 'QI_t_KF_conv_SH', &
370  'QI tendency due to flux convergence (shallow convection) in KF', &
371  'kg/m3/s', &
372  hist_id(i_hist_qi_fc,i_hist_s) )
373  call file_history_reg( 'QI_t_KF_updet_DP', &
374  'QI tendency due to detrainment of updraft (deep convection) in KF', &
375  'kg/m3/s', &
376  hist_id(i_hist_qi_ud,i_hist_d) )
377  call file_history_reg( 'QI_t_KF_updet_SH', &
378  'QI tendency due to detrainment of updraft (shallow convection) in KF', &
379  'kg/m3/s', &
380  hist_id(i_hist_qi_ud,i_hist_s) )
381 
382  call file_history_reg( 'QR_t_KF_conv_DP', &
383  'QR tendency due to flux convergence (deep convection) in KF', &
384  'kg/m3/s', &
385  hist_id(i_hist_qr_fc,i_hist_d) )
386  call file_history_reg( 'QR_t_KF_conv_SH', &
387  'QR tendency due to flux convergence (shallow convection) in KF', &
388  'kg/m3/s', &
389  hist_id(i_hist_qr_fc,i_hist_s) )
390  call file_history_reg( 'QR_t_KF_rainfall_DP', &
391  'QR tendency due to rain fall (deep convection) in KF', &
392  'kg/m3/s', &
393  hist_id(i_hist_qr_rf,i_hist_d) )
394  call file_history_reg( 'QR_t_KF_rainfall_SH', &
395  'QR tendency due to rain fall (shallow convection) in KF', &
396  'kg/m3/s', &
397  hist_id(i_hist_qr_rf,i_hist_s) )
398 
399  call file_history_reg( 'QS_t_KF_conv_DP', &
400  'QS tendency due to flux convergence (deep convection) in KF', &
401  'kg/m3/s', &
402  hist_id(i_hist_qs_fc,i_hist_d) )
403  call file_history_reg( 'QS_t_KF_conv_SH', &
404  'QS tendency due to flux convergence (shallow convection) in KF', &
405  'kg/m3/s', &
406  hist_id(i_hist_qs_fc,i_hist_s) )
407  call file_history_reg( 'QS_t_KF_snowfall_DP', &
408  'QS tendency due to snow fall (deep convection) in KF', &
409  'kg/m3/s', &
410  hist_id(i_hist_qs_sf,i_hist_d) )
411  call file_history_reg( 'QS_t_KF_snowfall_SH', &
412  'QS tendency due to snow fall (shallow convection) in KF', &
413  'kg/m3/s', &
414  hist_id(i_hist_qs_sf,i_hist_s) )
415 
416  do n = 1, hist_vmax
417  if ( hist_id(n,0) > 0 .or. hist_id(n,1) > 0 ) then
418  allocate( hist_work(ka,ia,ja,hist_vmax,0:1) )
419  hist_work(:,:,:,:,:) = 0.0_rp
420  exit
421  end if
422  end do
423 
424  !$acc update device( TTAB, QSTAB, THE0K, ALU, RDPR, RDTHK, PLUTOP, PBOT, GdCP, &
425  !$acc RATE, TRIGGER_type, DELCAPE, DEEPLIFETIME, SHALLOWLIFETIME, &
426  !$acc DEPTH_USL, WARMRAIN, KF_LOG, kf_threshold, prec_type, DO_prec, &
427  !$acc lifetime, I_convflag, deltaz, Z, deltax )
428 
429  return

References scale_file_history::file_history_reg(), scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ie, scale_io::io_fid_conf, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_tracer::k, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, and scale_prc::prc_abort().

Referenced by mod_atmos_phy_cp_driver::atmos_phy_cp_driver_setup().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_cp_kf_finalize()

subroutine, public scale_atmos_phy_cp_kf::atmos_phy_cp_kf_finalize

finalize

Definition at line 435 of file scale_atmos_phy_cp_kf.F90.

435 
436  deallocate( lifetime )
437  deallocate( i_convflag )
438  deallocate( z )
439  deallocate( deltaz )
440  deallocate( deltax )
441 
442  if ( allocated( hist_work ) ) deallocate( hist_work )
443 
444  return

Referenced by mod_atmos_phy_cp_driver::atmos_phy_cp_driver_finalize().

Here is the caller graph for this function:

◆ atmos_phy_cp_kf_putvar()

subroutine, public scale_atmos_phy_cp_kf::atmos_phy_cp_kf_putvar ( real(rp), intent(in), optional  in_delcape,
real(rp), intent(in), optional  in_deeplifetime,
real(rp), intent(in), optional  in_shallowlifetime 
)

overwrite private variables

Definition at line 453 of file scale_atmos_phy_cp_kf.F90.

453  implicit none
454 
455  real(RP), intent(in), optional :: in_delcape
456  real(RP), intent(in), optional :: in_deeplifetime
457  real(RP), intent(in), optional :: in_shallowlifetime
458  !---------------------------------------------------------------------------
459 
460  if( present(in_delcape ) ) delcape = in_delcape
461  if( present(in_deeplifetime ) ) deeplifetime = in_deeplifetime
462  if( present(in_shallowlifetime) ) shallowlifetime = in_shallowlifetime
463 
464  return

Referenced by mod_da_param_estimation::da_param_estimation_update().

Here is the caller graph for this function:

◆ atmos_phy_cp_kf_getvar()

subroutine, public scale_atmos_phy_cp_kf::atmos_phy_cp_kf_getvar ( real(rp), intent(out), optional  out_delcape,
real(rp), intent(out), optional  out_deeplifetime,
real(rp), intent(out), optional  out_shallowlifetime 
)

read private variables

Definition at line 473 of file scale_atmos_phy_cp_kf.F90.

473  implicit none
474 
475  real(RP), intent(out), optional :: out_delcape
476  real(RP), intent(out), optional :: out_deeplifetime
477  real(RP), intent(out), optional :: out_shallowlifetime
478  !---------------------------------------------------------------------------
479 
480  if( present(out_delcape ) ) out_delcape = delcape
481  if( present(out_deeplifetime ) ) out_deeplifetime = deeplifetime
482  if( present(out_shallowlifetime) ) out_shallowlifetime = shallowlifetime
483 
484  return

References scale_prc::prc_abort().

Referenced by mod_da_param_estimation::da_param_estimation_update().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_cp_kf_tendency()

subroutine, public scale_atmos_phy_cp_kf::atmos_phy_cp_kf_tendency ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension (ka,ia,ja), intent(in)  DENS,
real(rp), dimension (ka,ia,ja), intent(in)  U,
real(rp), dimension (ka,ia,ja), intent(in)  V,
real(rp), dimension (ka,ia,ja), intent(in)  RHOT,
real(rp), dimension (ka,ia,ja), intent(in)  TEMP,
real(rp), dimension (ka,ia,ja), intent(in)  PRES,
real(rp), dimension (ka,ia,ja), intent(in)  QDRY,
real(rp), dimension(ka,ia,ja), intent(in)  QV_in,
real(rp), dimension(ka,ia,ja), intent(in)  w0avg,
real(rp), dimension (0:ka,ia,ja), intent(in)  FZ,
real(dp), intent(in)  KF_DTSECD,
real(rp), dimension (ka,ia,ja), intent(inout)  DENS_t,
real(rp), dimension (ka,ia,ja), intent(inout)  RHOT_t,
real(rp), dimension(ka,ia,ja), intent(inout)  RHOQV_t,
real(rp), dimension (ka,ia,ja,n_hyd), intent(inout)  RHOQ_t,
real(rp), dimension (ia,ja), intent(inout)  nca,
real(rp), dimension(ia,ja), intent(inout)  SFLX_rain,
real(rp), dimension(ia,ja), intent(inout)  SFLX_snow,
real(rp), dimension(ia,ja), intent(inout)  SFLX_engi,
real(rp), dimension (ia,ja), intent(inout)  cloudtop,
real(rp), dimension (ia,ja), intent(inout)  cloudbase,
real(rp), dimension (ka,ia,ja), intent(inout)  cldfrac_dp,
real(rp), dimension (ka,ia,ja), intent(inout)  cldfrac_sh 
)

ATMOS_PHY_CP_kf calculate Kain-Fritsch Cumulus Parameterization.

Parameters
[in]densDensity [kg/m3]
[in]uvelocity u [m/s]
[in]vvelocity v [m/s]
[in]rhotDENS * POTT [k*kg/m3]
[in]temptemperature [K]
[in]prespressure of dry air [Pa]
[in]qdrydry air [1]
[in]qv_inspecific humidity [kg/kg]
[in]w0avgrunning mean of vertical velocity [m/s]
[in,out]dens_ttendency DENS [kg/m3/s]
[in,out]rhot_ttendency RHOT [K*kg/m3/s]
[in,out]rhoqv_ttendency rho*QV [kg/kg/s]
[in,out]rhoq_ttendency rho*QTRC [kg/kg/s]
[in,out]ncaconvection active time [sec]
[in,out]sflx_rainconvective rain rate [kg/m2/s]
[in,out]sflx_snowconvective snow rate [kg/m2/s]
[in,out]sflx_engiprecipitation internal energy [J/m2/s]
[in,out]cloudtopcloud top height [m]
[in,out]cloudbasecloud base height [m]
[in,out]cldfrac_dpcloud fraction (deep convection)
[in,out]cldfrac_shcloud fraction (shallow convection)

Definition at line 578 of file scale_atmos_phy_cp_kf.F90.

578  use scale_prc, only: &
579  prc_abort
580  use scale_file_history, only: &
581  file_history_query, &
582  file_history_put, &
583  file_history_in
584  use scale_const, only: &
585  grav => const_grav, &
586  pre00 => const_pre00, &
587  epsvap => const_epsvap
588  use scale_atmos_hydrometeor, only: &
589  cp_vapor, &
590  cp_water, &
591  cp_ice, &
592  cv_vapor, &
593  cv_water, &
594  cv_ice, &
595  n_hyd, &
596  i_hc, &
597  i_hr, &
598  i_hi, &
599  i_hs
600  use scale_atmos_saturation ,only :&
601  saturation_psat_liq => atmos_saturation_psat_liq
602  implicit none
603  integer, intent(in) :: KA, KS, KE
604  integer, intent(in) :: IA, IS, IE
605  integer, intent(in) :: JA, JS, JE
606 
607  real(RP), intent(in) :: DENS (KA,IA,JA)
608  real(RP), intent(in) :: U (KA,IA,JA)
609  real(RP), intent(in) :: V (KA,IA,JA)
610  real(RP), intent(in) :: RHOT (KA,IA,JA)
611  real(RP), intent(in) :: TEMP (KA,IA,JA)
612  real(RP), intent(in) :: PRES (KA,IA,JA)
613  real(RP), intent(in) :: QDRY (KA,IA,JA)
614  real(RP), intent(in) :: QV_in(KA,IA,JA)
615  real(RP), intent(in) :: w0avg(KA,IA,JA)
616  real(RP), intent(in) :: FZ (0:KA,IA,JA)
617  real(DP), intent(in) :: KF_DTSECD
618 
619  real(RP), intent(inout) :: DENS_t (KA,IA,JA)
620  real(RP), intent(inout) :: RHOT_t (KA,IA,JA)
621  real(RP), intent(inout) :: RHOQV_t(KA,IA,JA)
622  real(RP), intent(inout) :: RHOQ_t (KA,IA,JA,N_HYD)
623  real(RP), intent(inout) :: nca (IA,JA)
624 
625  real(RP), intent(inout) :: SFLX_rain(IA,JA)
626  real(RP), intent(inout) :: SFLX_snow(IA,JA)
627  real(RP), intent(inout) :: SFLX_engi(IA,JA)
628  real(RP), intent(inout) :: cloudtop (IA,JA)
629  real(RP), intent(inout) :: cloudbase (IA,JA)
630  real(RP), intent(inout) :: cldfrac_dp (KA,IA,JA)
631  real(RP), intent(inout) :: cldfrac_sh (KA,IA,JA)
632 
633  integer :: k, i, j, iq
634  integer :: nic
635  integer :: k_lcl
636  integer :: k_top
637  integer :: k_ml
638  integer :: k_lc,k_let,k_pbl
639  integer :: k_lfs
640 
641  real(RP) :: qv (KA)
642  real(RP) :: PSAT
643  real(RP) :: QSAT (KA)
644  real(RP) :: rh (KA)
645  real(RP) :: deltap(KA)
646 
647  real(RP) :: cldfrac_KF(KA,2)
648  real(RP) :: dens_nw(KA)
649  real(RP) :: time_advec
650  real(RP) :: umf(KA)
651  real(RP) :: umflcl
652  real(RP) :: upent(KA)
653  real(RP) :: updet(KA)
654  ! updraft value
655  real(RP) :: temp_u(KA)
656  real(RP) :: tempv(KA)
657  real(RP) :: qv_u(KA)
658  real(RP) :: cape
659  ! water
660  real(RP) :: qc(KA)
661  real(RP) :: qi(KA)
662  real(RP) :: qvdet(KA)
663  real(RP) :: qcdet(KA)
664  real(RP) :: qidet(KA)
665  real(RP) :: flux_qr(KA)
666  real(RP) :: flux_qs(KA)
667  ! upward theta
668  real(RP) :: theta_eu(KA)
669  real(RP) :: theta_ee(KA)
670  real(RP) :: zmix
671  real(RP) :: presmix
672  real(RP) :: umfnewdold(KA)
673  real(RP) :: dpthmx
674  real(RP) :: ems(KA)
675  real(RP) :: emsd(KA)
676  ! output from downdraft
677  real(RP) :: wspd(3)
678  real(RP) :: dmf(KA)
679  real(RP) :: downent(KA)
680  real(RP) :: downdet(KA)
681  real(RP) :: theta_d(KA)
682  real(RP) :: qv_d(KA)
683  real(RP) :: rain_flux
684  real(RP) :: snow_flux
685  real(RP) :: prec_engi
686 
687  ! update variables
688  real(RP) :: theta_nw(KA)
689  real(RP) :: qv_nw(KA)
690  real(RP) :: qc_nw(KA)
691  real(RP) :: qi_nw(KA)
692  real(RP) :: qr_nw(KA)
693  real(RP) :: qs_nw(KA)
694 
695  real(RP) :: RHOD(KA)
696 
697  real(RP) :: dQV, dQC, dQI, dQR, dQS
698 
699  real(RP) :: R_convflag(IA,JA)
700 
701  real(RP) :: KF_DTSEC
702 
703  logical :: error
704 
705  logical :: flag, logical
706  integer :: n, m
707 
708 #ifdef _OPENACC
709  real(RP) :: work(KA,24)
710  real(RP) :: work2(0:KA,7)
711  integer :: iwork(KA,1)
712 #endif
713 
714  ! ------
715 
716  log_progress(*) 'atmosphere / physics / cumulus / KF'
717  log_info("ATMOS_PHY_CP_kf_tendency",*) 'KF Convection Check '
718 
719  kf_dtsec = kf_dtsecd ! DP to RP
720 
721  call prof_rapstart('CP_kf', 3)
722 
723  do n = 1, hist_vmax
724  call file_history_query( hist_id(n,0), hist_flag )
725  if ( hist_flag ) exit
726  call file_history_query( hist_id(n,1), hist_flag )
727  if ( hist_flag ) exit
728  end do
729  !$acc update device( hist_flag )
730 
731  error = .false.
732 
733  !$omp parallel do default(none) schedule(dynamic) collapse(2) reduction(.or.:error) &
734  !$omp private(RHOD,tempv,qv_d,qc,qi,PSAT,QSAT,QV,rh, &
735  !$omp dens_nw,theta_nw,qv_nw,qc_nw,qi_nw,qr_nw,qs_nw, &
736  !$omp flux_qr,flux_qs,theta_eu,theta_ee,theta_d,umfnewdold,qvdet,qcdet,qidet, &
737  !$omp umf,umflcl,cape,presmix,upent,updet,temp_u,qv_u, &
738  !$omp dpthmx,zmix,rain_flux,snow_flux,prec_engi,time_advec, &
739  !$omp k_lcl,k_lc,k_lfs,k_pbl,k_top,k_let,k_ml, &
740  !$omp nic,deltap,cldfrac_KF,ems,emsd,wspd,dmf,downent,downdet, &
741  !$omp dQv,dQC,dQR,dQI,dQS) &
742  !$omp shared(DENS,RHOT,U,V,QDRY,TEMP,PRES,QV_in,FZ,Z,w0avg, &
743  !$omp DENS_t,RHOT_t,RHOQV_t,RHOQ_t, &
744  !$omp GRAV,PRE00,CP_VAPOR,CP_WATER,CP_ICE,CV_VAPOR,CV_WATER,CV_ICE,EPSvap, &
745  !$omp WARMRAIN,KF_DTSEC,SHALLOWLIFETIME, &
746  !$omp KA,KS,KE,IS,IE,JS,JE, &
747  !$omp nca,cloudtop,cloudbase,SFLX_rain,SFLX_snow,SFLX_engi, &
748  !$omp lifetime,deltaz,deltax,I_convflag, &
749  !$omp cldfrac_sh,cldfrac_dp,hist_flag,hist_work)
750  !$acc kernels
751  !$acc loop collapse(2) reduction(.or.:error) independent &
752  !$acc private(RHOD,tempv,qv_d,qc,qi,PSAT,QSAT,QV,rh, &
753  !$acc dens_nw,theta_nw,qv_nw,qc_nw,qi_nw,qr_nw,qs_nw, &
754  !$acc flux_qr,flux_qs,theta_eu,theta_ee,theta_d,umfnewdold,qvdet,qcdet,qidet, &
755  !$acc umf,umflcl,cape,presmix,upent,updet,temp_u,qv_u, &
756  !$acc dpthmx,zmix,rain_flux,snow_flux,prec_engi,time_advec, &
757  !$acc k_lcl,k_lc,k_lfs,k_pbl,k_top,k_let,k_ml, &
758  !$acc nic,deltap,cldfrac_KF,ems,emsd,wspd,dmf,downent,downdet, &
759  !$acc dQv,dQC,dQR,dQI,dQS, &
760  !$acc work,work2,iwork)
761  do j = js, je
762  do i = is, ie
763 
764  nca(i,j) = nca(i,j) - kf_dtsec
765 
766  ! check convection
767  if ( nca(i,j) .lt. 0.5_rp * kf_dtsec ) then
768 
769  do k = ks, ke
770  ! preparing a NON Hydriometeor condition to fit assumption in KF scheme
771  rhod(k) = dens(k,i,j) * qdry(k,i,j)
772  enddo
773 
774  do k = ks, ke
775  ! temporary: WRF TYPE equations are used to maintain consistency with kf_main
776  !call SATURATION_psat_liq( PSAT, TEMP(k,i,j) )
777  !QSAT(k) = 0.622_RP * PSAT / ( PRES(k,i,j) - ( 1.0_RP-0.622_RP ) * PSAT )
778  psat = aliq*exp((bliq*temp(k,i,j)-cliq)/(temp(k,i,j)-dliq))
779  qsat(k) = epsvap * psat / ( pres(k,i,j) - psat )
780 
781  ! calculate water vaper and relative humidity
782  qv(k) = max( kf_eps, min( qsat(k), qv_in(k,i,j) / qdry(k,i,j) ) ) ! compare QSAT and QV, guess lower limit
783  rh(k) = qv(k) / qsat(k)
784  enddo
785 
786  ! calculate delta P by hydrostatic balance
787  ! deltap is the pressure interval between half levels(face levels) @ SCALE
788  do k = ks, ke
789  deltap(k) = rhod(k) * grav * ( fz(k,i,j) - fz(k-1,i,j) ) ! rho*g*dz
790  enddo
791 
792  call cp_kf_trigger ( &
793  ka, ks, ke, & ! [IN]
794  deltaz(:,i,j), z(:,i,j), & ! [IN]
795  qv(:), qsat(:), & ! [IN]
796  pres(:,i,j), & ! [IN]
797  deltap(:), deltax(i,j), & ! [IN]
798  temp(:,i,j), & ! [IN]
799  w0avg(:,i,j), & ! [IN]
800 #ifdef _OPENACC
801  work(:,:), iwork(:,:), & ! [WORK]
802 #endif
803  i_convflag(i,j), & ! [OUT]
804  cloudtop(i,j), & ! [OUT]
805  temp_u(:), tempv(:), & ! [OUT]
806  qv_u(:), qc(:), qi(:), & ! [OUT]
807  qvdet(:), qcdet(:), qidet(:), & ! [OUT]
808  flux_qr(:), flux_qs(:), & ! [OUT]
809  theta_eu(:), theta_ee(:), & ! [OUT]
810  cape, & ! [OUT]
811  umf(:), umflcl, & ! [OUT]
812  upent(:), updet(:), & ! [OUT]
813  k_lcl, k_lc, k_pbl, & ! [OUT]
814  k_top, k_let, k_ml, & ! [OUT]
815  presmix, & ! [OUT]
816  dpthmx, & ! [OUT]
817  cloudbase(i,j), zmix, & ! [OUT]
818  umfnewdold(:), & ! [OUT]
819  error ) ! [OUT]
820 
821  if (i_convflag(i,j) /= 2) then ! convection allowed I_convflag=0 or 1
822 
823  ! calc ems(box weight[kg])
824  ems(k_top+1:ke) = 0._rp
825  emsd(k_top+1:ke) = 0._rp
826  do k = ks, k_top
827  ems(k) = deltap(k) * deltax(i,j)**2 / grav
828  emsd(k) = 1._rp/ems(k)
829  umfnewdold(k) = 1._rp/umfnewdold(k)
830  end do
831 
832  call cp_kf_downdraft ( &
833  ka, ks, ke, & ! [IN]
834  i_convflag(i,j), & ! [IN]
835  k_lcl, k_ml, k_top, k_pbl, k_let, k_lc, & ! [IN]
836  z(:,i,j), cloudbase(i,j), & ! [IN]
837  u(:,i,j), v(:,i,j), rh(:), qv(:), pres(:,i,j), & ! [IN]
838  deltap(:), deltax(i,j), & ! [IN]
839  ems(:), & ! [IN]
840  theta_ee(:), & ! [IN]
841  umf(:), temp_u(:), & ! [IN]
842  flux_qr(:), flux_qs(:), tempv(:), & ! [IN]
843 #ifdef _OPENACC
844  work(:,:), & ! [WORK]
845 #endif
846  wspd(:), dmf(:), downent(:), downdet(:), & ! [OUT]
847  theta_d(:), qv_d(:), & ! [OUT]
848  rain_flux, snow_flux, prec_engi, & ! [OUT]
849  k_lfs ) ! [OUT]
850 
851  call cp_kf_compensational ( &
852  ka, ks, ke, & ! [IN]
853  k_top, k_lc, k_pbl, k_ml, k_lfs, & ! [IN]
854  deltaz(:,i,j), z(:,i,j), pres(:,i,j), deltap(:), deltax(i,j), & ! [IN]
855  temp(:,i,j), qv(:), ems(:), emsd(:), & ! [IN]
856  presmix, zmix, dpthmx, & ! [IN]
857  cape, & ! [IN]
858  temp_u(:), qvdet(:), umflcl, & ! [IN]
859  qc(:), qi(:), flux_qr(:), flux_qs(:), & ! [IN]
860  umfnewdold(:), & ! [IN]
861  wspd(:), & ! [IN]
862  qv_d(:), theta_d(:), & ! [IN]
863  prec_engi, & ! [IN]
864  kf_dtsec, & ! [IN]
865  rhod(:), i, j, & ! [IN]
866 #ifdef _OPENACC
867  work(:,:), work2(:,:), & ! [WORK]
868 #endif
869  i_convflag(i,j), k_lcl, & ! [INOUT]
870  umf(:), upent(:), updet(:), & ! [INOUT]
871  qcdet(:), qidet(:), dmf(:), downent(:), downdet(:), & ! [INOUT]
872  rain_flux, snow_flux, & ! [INOUT]
873  nic, & ! [OUT]
874  theta_nw(:), & ! [OUT]
875  qv_nw(:), qc_nw(:), qi_nw(:), qr_nw(:), qs_nw(:), & ! [OUT]
876  sflx_rain(i,j), sflx_snow(i,j), sflx_engi(i,j), & ! [OUT]
877  cldfrac_kf, lifetime(i,j), time_advec, & ! [OUT]
878  error ) ! [OUT]
879  end if
880 
881  if (i_convflag(i,j) == 2) then ! no convection
882 
883  sflx_rain(i,j) = 0.0_rp
884  sflx_snow(i,j) = 0.0_rp
885  sflx_engi(i,j) = 0.0_rp
886  cldfrac_kf(ks:ke,:) = 0.0_rp
887  lifetime(i,j) = 0.0_rp
888  cloudtop(i,j) = 0.0_rp
889  cloudbase(i,j) = 0.0_rp
890  nca(i,j) = -100.0_rp
891 
892  dens_t(ks:ke,i,j) = 0.0_rp
893  rhot_t(ks:ke,i,j) = 0.0_rp
894  rhoqv_t(ks:ke,i,j) = 0.0_rp
895  do iq = 1, n_hyd
896  rhoq_t(ks:ke,i,j,iq) = 0.0_rp
897  end do
898 
899  if ( hist_flag ) then
900  do m = 0, 1
901  do n = 1, hist_vmax
902  hist_work(:,i,j,n,m) = 0.0_rp
903  end do
904  end do
905  end if
906 
907  else
908 
909  ! compute tendencys
910 
911  ! check:
912  ! FEEDBACK TO RESOLVABLE SCALE TENDENCIES.
913  ! IF THE ADVECTIVE TIME PERIOD (time_advec) IS LESS THAN SPECIFIED MINIMUM
914  ! lifetime, ALLOW FEEDBACK TO OCCUR ONLY DURING TADVEC.
915  !
916  if (i_convflag(i,j) == 0) then ! deep
917  if (time_advec < lifetime(i,j)) nic=nint(time_advec/kf_dtsec)
918  nca(i,j) = nic * kf_dtsec ! convection feed back act this time span
919  elseif (i_convflag(i,j) == 1) then ! shallow
920  lifetime(i,j) = max(shallowlifetime, kf_dtsec)
921  nca(i,j) = kf_dtsec ! convection feed back act this time span
922  end if
923 
924  do k = ks, k_top
925  ! vapor
926  dqv = rhod(k) * ( qv_nw(k) - qv(k) )
927  rhoqv_t(k,i,j) = dqv / lifetime(i,j)
928  ! liquid water
929  dqc = qc_nw(k) * rhod(k)
930  dqr = qr_nw(k) * rhod(k)
931  rhoq_t(k,i,j,i_hc) = dqc / lifetime(i,j)
932  rhoq_t(k,i,j,i_hr) = dqr / lifetime(i,j)
933  ! ice water
934  if ( .not. warmrain ) then
935  dqi = qi_nw(k) * rhod(k)
936  dqs = qs_nw(k) * rhod(k)
937  else
938  dqi = 0.0_rp
939  dqs = 0.0_rp
940  end if
941  rhoq_t(k,i,j,i_hi) = dqi / lifetime(i,j)
942  rhoq_t(k,i,j,i_hs) = dqs / lifetime(i,j)
943 
944  dens_t(k,i,j) = rhoqv_t(k,i,j) &
945  + rhoq_t(k,i,j,i_hc) + rhoq_t(k,i,j,i_hr) &
946  + rhoq_t(k,i,j,i_hi) + rhoq_t(k,i,j,i_hs)
947  dens_nw(k) = dens(k,i,j) + dqv + dqc + dqr + dqi + dqs
948 
949  end do
950 
951  do k = ks, k_top
952  do iq = i_hs+1, n_hyd
953  rhoq_t(k,i,j,iq) = 0.0_rp
954  end do
955  end do
956 
957  ! internal energy
958  do k = ks, k_top
959  rhot_t(k,i,j) = ( dens_nw(k) * theta_nw(k) - rhot(k,i,j) ) / lifetime(i,j)
960  end do
961 
962  do k=k_top+1, ke
963  dens_t(k,i,j) = 0.0_rp
964  rhot_t(k,i,j) = 0.0_rp
965  rhoqv_t(k,i,j) = 0.0_rp
966  do iq = 1, n_hyd
967  rhoq_t(k,i,j,iq) = 0.0_rp
968  end do
969  end do
970 
971  ! to keep conservation
972  ! if noconvection then nca is same value before call. nca only modifyed convectioned
973  end if
974 
975  cldfrac_sh(ks:ke,i,j) = cldfrac_kf(ks:ke,1)
976  cldfrac_dp(ks:ke,i,j) = cldfrac_kf(ks:ke,2)
977 
978  endif ! nca
979 
980  end do
981  end do
982  !$acc end kernels
983 
984 
985  if ( error ) call prc_abort
986 
987 
988  call prof_rapend('CP_kf', 3)
989 
990  ! history
991  call file_history_in( lifetime(:,:), 'KF_LIFETIME', 'lifetime of KF scheme', 's' )
992  r_convflag(:,:) = real(i_convflag(:,:),rp)
993  call file_history_in( r_convflag(:,:), 'KF_CONVFLAG', 'CONVECTION FLAG', '' )
994  do m = 0, 1
995  do n = 1, hist_vmax
996  call file_history_query( hist_id(n,m), flag )
997  if ( flag ) then
998  call file_history_put( hist_id(n,m), hist_work(:,:,:,n,m) )
999  end if
1000  end do
1001  end do
1002 
1003  return

References scale_const::const_cpdry, scale_const::const_cpvap, scale_const::const_emelt, scale_const::const_eps, scale_const::const_epstvap, scale_const::const_epsvap, scale_const::const_grav, scale_const::const_pre00, scale_const::const_rdry, scale_const::const_rvap, scale_const::const_tem00, scale_atmos_hydrometeor::cp_ice, scale_atmos_hydrometeor::cp_vapor, scale_atmos_hydrometeor::cp_water, scale_atmos_hydrometeor::cv_ice, scale_atmos_hydrometeor::cv_vapor, scale_atmos_hydrometeor::cv_water, scale_atmos_hydrometeor::i_hc, scale_atmos_hydrometeor::i_hi, scale_atmos_hydrometeor::i_hr, scale_atmos_hydrometeor::i_hs, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::js, scale_tracer::k, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, scale_atmos_hydrometeor::lhf, scale_atmos_hydrometeor::lhs, scale_atmos_hydrometeor::n_hyd, scale_prc::prc_abort(), scale_prof::prof_rapend(), scale_prof::prof_rapstart(), and scale_precision::rp.

Referenced by mod_atmos_phy_cp_driver::atmos_phy_cp_driver_calc_tendency().

Here is the call graph for this function:
Here is the caller graph for this function:
scale_const::const_grav
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:49
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_atmos_hydrometeor::i_hr
integer, parameter, public i_hr
liquid water rain
Definition: scale_atmos_hydrometeor.F90:98
scale_atmos_hydrometeor::cp_water
real(rp), public cp_water
CP for water [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:152
scale_atmos_hydrometeor::i_hs
integer, parameter, public i_hs
snow
Definition: scale_atmos_hydrometeor.F90:100
scale_const::const_epsvap
real(rp), public const_epsvap
Rdry / Rvap.
Definition: scale_const.F90:75
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
scale_file_history
module file_history
Definition: scale_file_history.F90:15
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_atmos_hydrometeor::cv_vapor
real(rp), public cv_vapor
CV for vapor [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:149
scale_atmos_hydrometeor::i_hi
integer, parameter, public i_hi
ice water cloud
Definition: scale_atmos_hydrometeor.F90:99
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_hydrometeor::i_hc
integer, parameter, public i_hc
liquid water cloud
Definition: scale_atmos_hydrometeor.F90:97
scale_file_history::file_history_reg
subroutine, public file_history_reg(name, desc, unit, itemid, standard_name, ndims, dim_type, cell_measures, fill_halo)
Register/Append variable to history file.
Definition: scale_file_history.F90:685
scale_atmos_saturation
module atmosphere / saturation
Definition: scale_atmos_saturation.F90:12
scale_const::const_pre00
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:97
scale_atmos_hydrometeor::n_hyd
integer, parameter, public n_hyd
Definition: scale_atmos_hydrometeor.F90:95
scale_atmos_hydrometeor::cp_vapor
real(rp), public cp_vapor
CP for vapor [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:150
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
scale_atmos_hydrometeor::cp_ice
real(rp), public cp_ice
CP for ice [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:154