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_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 179 of file scale_atmos_phy_cp_kf.F90.

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

503  use scale_file_history, only: &
504  file_history_query, &
505  file_history_put, &
506  file_history_in
507  use scale_const, only: &
508  grav => const_grav, &
509  pre00 => const_pre00, &
510  epsvap => const_epsvap
511  use scale_atmos_hydrometeor, only: &
512  cp_vapor, &
513  cp_water, &
514  cp_ice, &
515  cv_vapor, &
516  cv_water, &
517  cv_ice, &
518  n_hyd, &
519  i_hc, &
520  i_hr, &
521  i_hi, &
522  i_hs
523  use scale_atmos_saturation ,only :&
524  saturation_psat_liq => atmos_saturation_psat_liq
525  implicit none
526  integer, intent(in) :: KA, KS, KE
527  integer, intent(in) :: IA, IS, IE
528  integer, intent(in) :: JA, JS, JE
529 
530  real(RP), intent(in) :: DENS (KA,IA,JA)
531  real(RP), intent(in) :: U (KA,IA,JA)
532  real(RP), intent(in) :: V (KA,IA,JA)
533  real(RP), intent(in) :: RHOT (KA,IA,JA)
534  real(RP), intent(in) :: TEMP (KA,IA,JA)
535  real(RP), intent(in) :: PRES (KA,IA,JA)
536  real(RP), intent(in) :: QDRY (KA,IA,JA)
537  real(RP), intent(in) :: QV_in(KA,IA,JA)
538  real(RP), intent(in) :: w0avg(KA,IA,JA)
539  real(RP), intent(in) :: FZ (0:KA,IA,JA)
540  real(DP), intent(in) :: KF_DTSECD
541 
542  real(RP), intent(inout) :: DENS_t (KA,IA,JA)
543  real(RP), intent(inout) :: RHOT_t (KA,IA,JA)
544  real(RP), intent(inout) :: RHOQV_t(KA,IA,JA)
545  real(RP), intent(inout) :: RHOQ_t (KA,IA,JA,N_HYD)
546  real(RP), intent(inout) :: nca (IA,JA)
547 
548  real(RP), intent(inout) :: SFLX_rain(IA,JA)
549  real(RP), intent(inout) :: SFLX_snow(IA,JA)
550  real(RP), intent(inout) :: SFLX_engi(IA,JA)
551  real(RP), intent(inout) :: cloudtop (IA,JA)
552  real(RP), intent(inout) :: cloudbase (IA,JA)
553  real(RP), intent(inout) :: cldfrac_dp (KA,IA,JA)
554  real(RP), intent(inout) :: cldfrac_sh (KA,IA,JA)
555 
556  integer :: k, i, j, iq
557  integer :: nic
558  integer :: k_lcl
559  integer :: k_top
560  integer :: k_ml
561  integer :: k_lc,k_let,k_pbl
562  integer :: k_lfs
563 
564  real(RP) :: qv (KA)
565  real(RP) :: PSAT
566  real(RP) :: QSAT (KA)
567  real(RP) :: rh (KA)
568  real(RP) :: deltap(KA)
569 
570  real(RP) :: cldfrac_KF(KA,2)
571  real(RP) :: dens_nw(KA)
572  real(RP) :: time_advec
573  real(RP) :: umf(KA)
574  real(RP) :: umflcl
575  real(RP) :: upent(KA)
576  real(RP) :: updet(KA)
577  ! updraft value
578  real(RP) :: temp_u(KA)
579  real(RP) :: tempv(KA)
580  real(RP) :: qv_u(KA)
581  real(RP) :: cape
582  ! water
583  real(RP) :: qc(KA)
584  real(RP) :: qi(KA)
585  real(RP) :: qvdet(KA)
586  real(RP) :: qcdet(KA)
587  real(RP) :: qidet(KA)
588  real(RP) :: flux_qr(KA)
589  real(RP) :: flux_qs(KA)
590  ! upward theta
591  real(RP) :: theta_eu(KA)
592  real(RP) :: theta_ee(KA)
593  real(RP) :: zmix
594  real(RP) :: presmix
595  real(RP) :: umfnewdold(KA)
596  real(RP) :: dpthmx
597  real(RP) :: ems(KA)
598  real(RP) :: emsd(KA)
599  ! output from downdraft
600  real(RP) :: wspd(3)
601  real(RP) :: dmf(KA)
602  real(RP) :: downent(KA)
603  real(RP) :: downdet(KA)
604  real(RP) :: theta_d(KA)
605  real(RP) :: qv_d(KA)
606  real(RP) :: rain_flux
607  real(RP) :: snow_flux
608  real(RP) :: prec_engi
609 
610  ! update variables
611  real(RP) :: theta_nw(KA)
612  real(RP) :: qv_nw(KA)
613  real(RP) :: qc_nw(KA)
614  real(RP) :: qi_nw(KA)
615  real(RP) :: qr_nw(KA)
616  real(RP) :: qs_nw(KA)
617 
618  real(RP) :: RHOD(KA)
619 
620  real(RP) :: dQV, dQC, dQI, dQR, dQS
621 
622  real(RP) :: R_convflag(IA,JA)
623 
624  real(RP) :: KF_DTSEC
625 
626  logical :: flag
627  integer :: n, m
628  ! ------
629 
630  log_progress(*) 'atmosphere / physics / cumulus / KF'
631  log_info("ATMOS_PHY_CP_kf_tendency",*) 'KF Convection Check '
632 
633  kf_dtsec = kf_dtsecd ! DP to RP
634 
635  call prof_rapstart('CP_kf', 3)
636 
637  do n = 1, hist_vmax
638  call file_history_query( hist_id(n,0), hist_flag )
639  if ( hist_flag ) exit
640  call file_history_query( hist_id(n,1), hist_flag )
641  if ( hist_flag ) exit
642  end do
643 
644  !$omp parallel do default(none) &
645  !$omp private(RHOD,tempv,qv_d,qc,qi,PSAT,QSAT,QV,rh, &
646  !$omp dens_nw,theta_nw,qv_nw,qc_nw,qi_nw,qr_nw,qs_nw, &
647  !$omp flux_qr,flux_qs,theta_eu,theta_ee,theta_d,umfnewdold,qvdet,qcdet,qidet, &
648  !$omp umf,umflcl,cape,presmix,upent,updet,temp_u,qv_u, &
649  !$omp dpthmx,zmix,rain_flux,snow_flux,prec_engi,time_advec, &
650  !$omp k_lcl,k_lc,k_lfs,k_pbl,k_top,k_let,k_ml, &
651  !$omp nic,deltap,cldfrac_KF,ems,emsd,wspd,dmf,downent,downdet, &
652  !$omp dQv,dQC,dQR,dQI,dQS) &
653  !$omp shared(DENS,RHOT,U,V,QDRY,TEMP,PRES,QV_in,FZ,Z,w0avg, &
654  !$omp DENS_t,RHOT_t,RHOQV_t,RHOQ_t, &
655  !$omp GRAV,PRE00,CP_VAPOR,CP_WATER,CP_ICE,CV_VAPOR,CV_WATER,CV_ICE,EPSvap, &
656  !$omp WARMRAIN,KF_DTSEC,SHALLOWLIFETIME, &
657  !$omp KA,KS,KE,IS,IE,JS,JE, &
658  !$omp nca,cloudtop,cloudbase,SFLX_rain,SFLX_snow,SFLX_engi, &
659  !$omp lifetime,deltaz,deltax,I_convflag, &
660  !$omp cldfrac_sh,cldfrac_dp,hist_flag,hist_work)
661  do j = js, je
662  do i = is, ie
663 
664  nca(i,j) = nca(i,j) - kf_dtsec
665 
666  ! check convection
667  if ( nca(i,j) .ge. 0.5_rp * kf_dtsec ) cycle
668 
669  do k = ks, ke
670  ! preparing a NON Hydriometeor condition to fit assumption in KF scheme
671  rhod(k) = dens(k,i,j) * qdry(k,i,j)
672  enddo
673 
674 
675  do k = ks, ke
676  ! temporary: WRF TYPE equations are used to maintain consistency with kf_main
677  !call SATURATION_psat_liq( PSAT, TEMP(k,i,j) )
678  !QSAT(k) = 0.622_RP * PSAT / ( PRES(k,i,j) - ( 1.0_RP-0.622_RP ) * PSAT )
679  psat = aliq*exp((bliq*temp(k,i,j)-cliq)/(temp(k,i,j)-dliq))
680  qsat(k) = epsvap * psat / ( pres(k,i,j) - psat )
681 
682  ! calculate water vaper and relative humidity
683  qv(k) = max( kf_eps, min( qsat(k), qv_in(k,i,j) / qdry(k,i,j) ) ) ! compare QSAT and QV, guess lower limit
684  rh(k) = qv(k) / qsat(k)
685  enddo
686 
687  ! calculate delta P by hydrostatic balance
688  ! deltap is the pressure interval between half levels(face levels) @ SCALE
689  do k = ks, ke
690  deltap(k) = rhod(k) * grav * ( fz(k,i,j) - fz(k-1,i,j) ) ! rho*g*dz
691  enddo
692 
693  call cp_kf_trigger ( &
694  ka, ks, ke, & ! [IN]
695  deltaz(:,i,j), z(:,i,j), & ! [IN]
696  qv(:), qsat(:), & ! [IN]
697  pres(:,i,j), & ! [IN]
698  deltap(:), deltax(i,j), & ! [IN]
699  temp(:,i,j), & ! [IN]
700  w0avg(:,i,j), & ! [IN]
701  i_convflag(i,j), & ! [OUT]
702  cloudtop(i,j), & ! [OUT]
703  temp_u(:), tempv(:), & ! [OUT]
704  qv_u(:), qc(:), qi(:), & ! [OUT]
705  qvdet(:), qcdet(:), qidet(:), & ! [OUT]
706  flux_qr(:), flux_qs(:), & ! [OUT]
707  theta_eu(:), theta_ee(:), & ! [OUT]
708  cape, & ! [OUT]
709  umf(:), umflcl, & ! [OUT]
710  upent(:), updet(:), & ! [OUT]
711  k_lcl, k_lc, k_pbl, & ! [OUT]
712  k_top, k_let, k_ml, & ! [OUT]
713  presmix, & ! [OUT]
714  dpthmx, & ! [OUT]
715  cloudbase(i,j), zmix, & ! [OUT]
716  umfnewdold(:) ) ! [OUT]
717 
718  if (i_convflag(i,j) /= 2) then ! convection allowed I_convflag=0 or 1
719 
720  ! calc ems(box weight[kg])
721  ems(k_top+1:ke) = 0._rp
722  emsd(k_top+1:ke) = 0._rp
723  do k = ks, k_top
724  ems(k) = deltap(k) * deltax(i,j)**2 / grav
725  emsd(k) = 1._rp/ems(k)
726  umfnewdold(k) = 1._rp/umfnewdold(k)
727  end do
728 
729  call cp_kf_downdraft ( &
730  ka, ks, ke, & ! [IN]
731  i_convflag(i,j), & ! [IN]
732  k_lcl, k_ml, k_top, k_pbl, k_let, k_lc, & ! [IN]
733  z(:,i,j), cloudbase(i,j), & ! [IN]
734  u(:,i,j), v(:,i,j), rh(:), qv(:), pres(:,i,j), & ! [IN]
735  deltap(:), deltax(i,j), & ! [IN]
736  ems(:), & ! [IN]
737  theta_ee(:), & ! [IN]
738  umf(:), temp_u(:), & ! [IN]
739  flux_qr(:), flux_qs(:), tempv(:), & ! [IN]
740  wspd(:), dmf(:), downent(:), downdet(:), & ! [OUT]
741  theta_d(:), qv_d(:), & ! [OUT]
742  rain_flux, snow_flux, prec_engi, & ! [OUT]
743  k_lfs ) ! [OUT]
744 
745  call cp_kf_compensational ( &
746  ka, ks, ke, & ! [IN]
747  k_top, k_lc, k_pbl, k_ml, k_lfs, & ! [IN]
748  deltaz(:,i,j), z(:,i,j), pres(:,i,j), deltap(:), deltax(i,j), & ! [IN]
749  temp(:,i,j), qv(:), ems(:), emsd(:), & ! [IN]
750  presmix, zmix, dpthmx, & ! [IN]
751  cape, & ! [IN]
752  temp_u(:), qvdet(:), umflcl, & ! [IN]
753  qc(:), qi(:), flux_qr(:), flux_qs(:), & ! [IN]
754  umfnewdold(:), & ! [IN]
755  wspd(:), & ! [IN]
756  qv_d(:), theta_d(:), & ! [IN]
757  prec_engi, & ! [IN]
758  kf_dtsec, & ! [IN]
759  rhod(:), i, j, & ! [IN]
760  i_convflag(i,j), k_lcl, & ! [INOUT]
761  umf(:), upent(:), updet(:), & ! [INOUT]
762  qcdet(:), qidet(:), dmf(:), downent(:), downdet(:), & ! [INOUT]
763  rain_flux, snow_flux, & ! [INOUT]
764  nic, & ! [OUT]
765  theta_nw(:), & ! [OUT]
766  qv_nw(:), qc_nw(:), qi_nw(:), qr_nw(:), qs_nw(:), & ! [OUT]
767  sflx_rain(i,j), sflx_snow(i,j), sflx_engi(i,j), & ! [OUT]
768  cldfrac_kf, lifetime(i,j), time_advec ) ! [OUT]
769 
770  end if
771 
772  if (i_convflag(i,j) == 2) then ! no convection
773 
774  sflx_rain(i,j) = 0.0_rp
775  sflx_snow(i,j) = 0.0_rp
776  sflx_engi(i,j) = 0.0_rp
777  cldfrac_kf(ks:ke,:) = 0.0_rp
778  lifetime(i,j) = 0.0_rp
779  cloudtop(i,j) = 0.0_rp
780  cloudbase(i,j) = 0.0_rp
781  nca(i,j) = -100.0_rp
782 
783  dens_t(ks:ke,i,j) = 0.0_rp
784  rhot_t(ks:ke,i,j) = 0.0_rp
785  rhoqv_t(ks:ke,i,j) = 0.0_rp
786  do iq = 1, n_hyd
787  rhoq_t(ks:ke,i,j,iq) = 0.0_rp
788  end do
789 
790  if ( hist_flag ) then
791  do m = 0, 1
792  do n = 1, hist_vmax
793  hist_work(:,i,j,n,m) = 0.0_rp
794  end do
795  end do
796  end if
797 
798  else
799 
800  ! compute tendencys
801 
802  ! check:
803  ! FEEDBACK TO RESOLVABLE SCALE TENDENCIES.
804  ! IF THE ADVECTIVE TIME PERIOD (time_advec) IS LESS THAN SPECIFIED MINIMUM
805  ! lifetime, ALLOW FEEDBACK TO OCCUR ONLY DURING TADVEC.
806  !
807  if (i_convflag(i,j) == 0) then ! deep
808  if (time_advec < lifetime(i,j)) nic=nint(time_advec/kf_dtsec)
809  nca(i,j) = nic * kf_dtsec ! convection feed back act this time span
810  elseif (i_convflag(i,j) == 1) then ! shallow
811  lifetime(i,j) = max(shallowlifetime, kf_dtsec)
812  nca(i,j) = kf_dtsec ! convection feed back act this time span
813  end if
814 
815  do k = ks, k_top
816  ! vapor
817  dqv = rhod(k) * ( qv_nw(k) - qv(k) )
818  rhoqv_t(k,i,j) = dqv / lifetime(i,j)
819  ! liquid water
820  dqc = qc_nw(k) * rhod(k)
821  dqr = qr_nw(k) * rhod(k)
822  rhoq_t(k,i,j,i_hc) = dqc / lifetime(i,j)
823  rhoq_t(k,i,j,i_hr) = dqr / lifetime(i,j)
824  ! ice water
825  if ( .not. warmrain ) then
826  dqi = qi_nw(k) * rhod(k)
827  dqs = qs_nw(k) * rhod(k)
828  else
829  dqi = 0.0_rp
830  dqs = 0.0_rp
831  end if
832  rhoq_t(k,i,j,i_hi) = dqi / lifetime(i,j)
833  rhoq_t(k,i,j,i_hs) = dqs / lifetime(i,j)
834 
835  dens_t(k,i,j) = rhoqv_t(k,i,j) &
836  + rhoq_t(k,i,j,i_hc) + rhoq_t(k,i,j,i_hr) &
837  + rhoq_t(k,i,j,i_hi) + rhoq_t(k,i,j,i_hs)
838  dens_nw(k) = dens(k,i,j) + dqv + dqc + dqr + dqi + dqs
839 
840  end do
841  do iq = i_hs+1, n_hyd
842  do k = ks, k_top
843  rhoq_t(k,i,j,iq) = 0.0_rp
844  end do
845  end do
846 
847  ! internal energy
848  do k = ks, k_top
849  rhot_t(k,i,j) = ( dens_nw(k) * theta_nw(k) - rhot(k,i,j) ) / lifetime(i,j)
850  end do
851 
852  do k=k_top+1, ke
853  dens_t(k,i,j) = 0.0_rp
854  rhot_t(k,i,j) = 0.0_rp
855  rhoqv_t(k,i,j) = 0.0_rp
856  do iq = 1, n_hyd
857  rhoq_t(k,i,j,iq) = 0.0_rp
858  end do
859  end do
860 
861  ! to keep conservation
862  ! if noconvection then nca is same value before call. nca only modifyed convectioned
863  end if
864 
865  cldfrac_sh(ks:ke,i,j) = cldfrac_kf(ks:ke,1)
866  cldfrac_dp(ks:ke,i,j) = cldfrac_kf(ks:ke,2)
867  end do
868  end do
869 
870  call prof_rapend('CP_kf', 3)
871 
872  ! history
873  call file_history_in( lifetime(:,:), 'KF_LIFETIME', 'lifetime of KF scheme', 's' )
874  r_convflag(:,:) = real(i_convflag(:,:),rp)
875  call file_history_in( r_convflag(:,:), 'KF_CONVFLAG', 'CONVECTION FLAG', '' )
876  do m = 0, 1
877  do n = 1, hist_vmax
878  call file_history_query( hist_id(n,m), flag )
879  if ( flag ) then
880  call file_history_put( hist_id(n,m), hist_work(:,:,:,n,m) )
881  end if
882  end do
883  end do
884 
885  return

References scale_const::const_cpdry, scale_const::const_cpvap, scale_const::const_emelt, 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:46
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_atmos_hydrometeor::i_hr
integer, parameter, public i_hr
liquid water rain
Definition: scale_atmos_hydrometeor.F90:82
scale_atmos_hydrometeor::cp_water
real(rp), public cp_water
CP for water [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:133
scale_atmos_hydrometeor::i_hs
integer, parameter, public i_hs
snow
Definition: scale_atmos_hydrometeor.F90:84
scale_const::const_epsvap
real(rp), public const_epsvap
Rdry / Rvap.
Definition: scale_const.F90:69
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:130
scale_atmos_hydrometeor::i_hi
integer, parameter, public i_hi
ice water cloud
Definition: scale_atmos_hydrometeor.F90:83
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:81
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:650
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:88
scale_atmos_hydrometeor::n_hyd
integer, parameter, public n_hyd
Definition: scale_atmos_hydrometeor.F90:79
scale_atmos_hydrometeor::cp_vapor
real(rp), public cp_vapor
CP for vapor [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:131
scale_atmos_hydrometeor::cv_water
real(rp), public cv_water
CV for water [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:132
scale_atmos_hydrometeor::cv_ice
real(rp), public cv_ice
CV for ice [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:134
scale_atmos_hydrometeor::cp_ice
real(rp), public cp_ice
CP for ice [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:135