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

module ATMOSPHERE / Physics Cumulus Parameterization More...

Functions/Subroutines

subroutine, public atmos_phy_cp_kf_setup (CP_TYPE)
 Setup. More...
 
subroutine, public atmos_phy_cp_kf (DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, DENS_t_CP, MOMZ_t_CP, MOMX_t_CP, MOMY_t_CP, RHOT_t_CP, RHOQ_t_CP, MFLX_cloudbase, SFLX_convrain, cloudtop, cloudbase, cldfrac_dp, cldfrac_sh, nca, w0avg)
 
subroutine kf_wmean (W0_avg, DENS, MOMZ)
 running mean vertical wind speed More...
 
subroutine cp_kf_param (RATE_in, TRIGGER_in, KF_DT_in, DELCAPE_in, DEEPLIFETIME_in, SHALLOWLIFETIME_in, DEPTH_USL_in, KF_prec_in, KF_threshold_in, KF_LOG_in, stepkf_in)
 
subroutine cp_kf_main (dens, MOMZ, MOMX, MOMY, RHOT, QTRC_in, w0avg, nca, DT_DENS, DT_RHOT, DT_RHOQ, rainrate_cp, cldfrac_sh, cldfrac_dp, timecp, cloudtop, zlcl, I_convflag)
 
subroutine kf_lutab
 

Variables

real(dp), save, public kf_dtsec = 300._RP
 

Detailed Description

module ATMOSPHERE / Physics Cumulus Parameterization

Description
Cumulus convection by Kain-Fritsch parameterization Reference: Kain and Fritsch(1990) Kain (2004) Narita and Ohmori (2007)
Author
Team SCALE
History
  • 2016-06-27 (S.Matsugishi) [new]
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
    PARAM_ATMOS_PHY_CP_KF_RATE real(RP) 0.03_RP ratio of cloud water and precipitation (Ogura and Cho 1973)
    PARAM_ATMOS_PHY_CP_KF_TRIGGER integer 1 trigger function type 1:KF2004 3:NO2007
    PARAM_ATMOS_PHY_CP_KF_DT real(DP) 5.0_DP KF convection check time interval [min]
    PARAM_ATMOS_PHY_CP_KF_DLCAPE real(RP) 0.1_RP cape decleace rate
    PARAM_ATMOS_PHY_CP_KF_DLIFETIME real(RP) 1800.0_RP minimum lifetime scale of deep convection [sec]
    PARAM_ATMOS_PHY_CP_KF_SLIFETIME real(RP) 2400.0_RP lifetime of shallow convection [sec]
    PARAM_ATMOS_PHY_CP_KF_DEPTH_USL real(RP) 300.0_RP depth of updraft source layer [hPa]
    PARAM_ATMOS_PHY_CP_KF_PREC integer 1 precipitation type 1:Ogura-Cho(1973) 2:Kessler
    PARAM_ATMOS_PHY_CP_KF_THRES real(RP) 1.E-3_RP autoconversion rate for Kessler
    PARAM_ATMOS_PHY_CP_KF_LOG logical .false. output log?
    PARAM_ATMOS_PHY_CP_KF_WADAPT logical .true.
    PARAM_ATMOS_PHY_CP_KF_W_TIME integer 16

History Output
namedescriptionunitvariable
KF_LIFETIME lifetime of KF scheme s lifetime
RP) KF_CONVFLAG CONVECTION FLAG real(I_convflag

Function/Subroutine Documentation

◆ atmos_phy_cp_kf_setup()

subroutine, public scale_atmos_phy_cp_kf::atmos_phy_cp_kf_setup ( character(len=*), intent(in)  CP_TYPE)

Setup.

Definition at line 147 of file scale_atmos_phy_cp_kf.F90.

References cp_kf_param(), scale_grid::dx, scale_grid::dy, scale_atmos_hydrometeor::i_qc, scale_atmos_hydrometeor::i_qi, scale_atmos_hydrometeor::i_qr, scale_atmos_hydrometeor::i_qs, scale_atmos_hydrometeor::i_qv, scale_grid_index::ia, scale_grid_index::ie, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_fid_nml, scale_stdio::io_l, scale_stdio::io_nml, scale_grid_index::is, scale_grid_index::ja, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ka, scale_grid_index::ke, kf_lutab(), scale_grid_index::ks, scale_process::prc_mpistop(), scale_grid_real::real_cz, scale_time::time_dtsec, and scale_time::time_dtsec_atmos_phy_cp.

Referenced by scale_atmos_phy_cp::atmos_phy_cp_setup().

147  use scale_process, only: &
149  use scale_time , only :&
150  time_dtsec, &
152  use scale_grid_real, only: &
153  cz => real_cz
154  use scale_grid,only: &
155  dx => dx, &
156  dy => dy
157  use scale_atmos_hydrometeor, only: &
158  i_qv, &
159  i_qc, &
160  i_qr, &
161  i_qi, &
162  i_qs
163  implicit none
164 
165  character(len=*), intent(in) :: CP_TYPE
166 
167  ! tunning parameters, original parameter set is from KF2004 and NO2007
168  real(RP) :: PARAM_ATMOS_PHY_CP_kf_rate = 0.03_rp ! ratio of cloud water and precipitation (Ogura and Cho 1973)
169  integer :: PARAM_ATMOS_PHY_CP_kf_trigger = 1 ! trigger function type 1:KF2004 3:NO2007
170  real(DP) :: PARAM_ATMOS_PHY_CP_kf_dt = 5.0_dp ! KF convection check time interval [min]
171  real(RP) :: PARAM_ATMOS_PHY_CP_kf_dlcape = 0.1_rp ! cape decleace rate
172  real(RP) :: PARAM_ATMOS_PHY_CP_kf_dlifetime = 1800.0_rp ! minimum lifetime scale of deep convection [sec]
173  real(RP) :: PARAM_ATMOS_PHY_CP_kf_slifetime = 2400.0_rp ! lifetime of shallow convection [sec]
174  real(RP) :: PARAM_ATMOS_PHY_CP_kf_DEPTH_USL = 300.0_rp ! depth of updraft source layer [hPa]
175  integer :: PARAM_ATMOS_PHY_CP_kf_prec = 1 ! precipitation type 1:Ogura-Cho(1973) 2:Kessler
176  real(RP) :: PARAM_ATMOS_PHY_CP_kf_thres = 1.e-3_rp ! autoconversion rate for Kessler
177  logical :: PARAM_ATMOS_PHY_CP_kf_LOG = .false. ! output log?
178 
179  namelist / param_atmos_phy_cp_kf / &
180  param_atmos_phy_cp_kf_rate, &
181  param_atmos_phy_cp_kf_trigger, &
182  param_atmos_phy_cp_kf_dt, &
183  param_atmos_phy_cp_kf_dlcape, &
184  param_atmos_phy_cp_kf_dlifetime, &
185  param_atmos_phy_cp_kf_slifetime, &
186  param_atmos_phy_cp_kf_depth_usl, &
187  param_atmos_phy_cp_kf_prec, &
188  param_atmos_phy_cp_kf_thres, &
189  param_atmos_phy_cp_kf_log, &
190  param_atmos_phy_cp_kf_wadapt, &
191  param_atmos_phy_cp_kf_w_time
192 
193  integer :: k, i, j
194  integer :: QS
195  integer :: ierr
196  !---------------------------------------------------------------------------
197 
198  if( io_l ) write(io_fid_log,*)
199  if( io_l ) write(io_fid_log,*) '++++++ Module[CUMULUS] / Categ[ATMOS PHYSICS] / Origin[SCALElib]'
200  if( io_l ) write(io_fid_log,*) '*** Kain-Fritsch scheme'
201 
202  if ( cp_type /= 'KF' ) then
203  write(*,*) 'xxx ATMOS_PHY_CP_TYPE is not KF. Check!'
204  call prc_mpistop
205  endif
206 
207  if ( i_qv < 1 ) then
208  write(*,*) 'xxx QV is not registered'
209  call prc_mpistop
210  end if
211  if ( i_qc < 1 ) then
212  write(*,*) 'xxx QC is not registered'
213  call prc_mpistop
214  end if
215  if ( i_qr < 1 ) then
216  write(*,*) 'xxx QR is not registered'
217  call prc_mpistop
218  end if
219 
220  if ( i_qi < 1 ) then
221  warmrain = .true.
222  else
223  warmrain = .false.
224  if ( i_qs < 1 ) then
225  write(*,*) 'xxx QI is registerd, but QS is not registered'
226  call prc_mpistop
227  end if
228  endif
229 
230  if ( abs(time_dtsec_atmos_phy_cp-time_dtsec) > 0.0_dp ) then
231  write(*,*) 'xxx TIME_DTSEC_ATMOS_PHY_CP should be same as TIME_DTSEC for KF scheme. STOP.'
232  call prc_mpistop
233  endif
234 
235  !--- read namelist
236  rewind(io_fid_conf)
237  read(io_fid_conf,nml=param_atmos_phy_cp_kf,iostat=ierr)
238  if( ierr < 0 ) then !--- missing
239  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
240  elseif( ierr > 0 ) then !--- fatal error
241  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_CP_KF. Check!'
242  call prc_mpistop
243  endif
244  if( io_nml ) write(io_fid_nml,nml=param_atmos_phy_cp_kf)
245 
246  call kf_lutab ! set kf look up table
247 
248  ! set kf convection check step
249  time_dstep_kf = nint( param_atmos_phy_cp_kf_dt * 60.0_dp / time_dtsec_atmos_phy_cp )
250  time_dstep_kf = max(time_dstep_kf,1) ! kf time interval step
251  time_res_kf = -1 ! initialize to keep consistent for below step
252  time_dokf = .true. ! initialize
253 
254  call cp_kf_param( param_atmos_phy_cp_kf_rate, &
255  param_atmos_phy_cp_kf_trigger, &
256  param_atmos_phy_cp_kf_dt, &
257  param_atmos_phy_cp_kf_dlcape, &
258  param_atmos_phy_cp_kf_dlifetime, &
259  param_atmos_phy_cp_kf_slifetime, &
260  param_atmos_phy_cp_kf_depth_usl, &
261  param_atmos_phy_cp_kf_prec, &
262  param_atmos_phy_cp_kf_thres, &
263  param_atmos_phy_cp_kf_log , &
264  time_dstep_kf )
265 
266  ! output parameter lists
267  if( io_l ) write(io_fid_log,*)
268  if( io_l ) write(io_fid_log,*) "*** Interval for check [step] : ", time_dstep_kf
269  if( io_l ) write(io_fid_log,*) "*** Ogura-Cho condense material convert rate : ", param_atmos_phy_cp_kf_rate
270  if( io_l ) write(io_fid_log,*) "*** Trigger function type, 1:KF2004 3:NO2007 : ", param_atmos_phy_cp_kf_trigger
271  if( io_l ) write(io_fid_log,*) "*** CAPE decrease rate : ", param_atmos_phy_cp_kf_dlcape
272  if( io_l ) write(io_fid_log,*) "*** Minimum lifetime scale of deep convection [sec] : ", param_atmos_phy_cp_kf_dlifetime
273  if( io_l ) write(io_fid_log,*) "*** Lifetime of shallow convection [sec] : ", param_atmos_phy_cp_kf_slifetime
274  if( io_l ) write(io_fid_log,*) "*** Updraft souce layer depth [hPa] : ", param_atmos_phy_cp_kf_depth_usl
275  if( io_l ) write(io_fid_log,*) "*** Precipitation type 1:Ogura-Cho(1973) 2:Kessler : ", param_atmos_phy_cp_kf_prec
276  if( io_l ) write(io_fid_log,*) "*** Kessler type precipitation's threshold : ", param_atmos_phy_cp_kf_thres
277  if( io_l ) write(io_fid_log,*) "*** Warm rain? : ", warmrain
278  if( io_l ) write(io_fid_log,*) "*** Use running mean of w in adaptive timestep? : ", param_atmos_phy_cp_kf_wadapt
279  if( io_l ) write(io_fid_log,*) "*** Fixed time scale for running mean of w : ", param_atmos_phy_cp_kf_w_time
280  if( io_l ) write(io_fid_log,*) "*** Output log? : ", param_atmos_phy_cp_kf_log
281 
282  ! output variables
283  allocate( lifetime(ia,ja) )
284  allocate( i_convflag(ia,ja) )
285  lifetime(:,:) = 0.0_rp
286  i_convflag(:,:) = 2
287 
288  allocate( z(ka,ia,ja) )
289  z(:,:,:) = cz(:,:,:) ! becouse scale_atmos_phy_cp interface ,not use scale_grid
290 
291  allocate( deltaz(ka,ia,ja) )
292  ! deltaz is the interval of between model full levels(scalar point )
293  deltaz(:,:,:) = 0.0_rp
294  do j = js, je
295  do i = is, ie
296  do k = ks, ke
297  deltaz(k,i,j) = cz(k+1,i,j) - cz(k,i,j)
298  enddo
299  enddo
300  enddo
301  deltaz(ke,:,:) = 0.0_rp
302 
303  deltax = sqrt( dx*dy )
304 
305  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
real(rp), public dy
length in the main region [m]: y
subroutine, public prc_mpistop
Abort MPI.
real(rp), public dx
length in the main region [m]: x
integer, public ke
end point of inner domain: z, local
real(rp), dimension(:,:,:), allocatable, public real_cz
geopotential height [m] (cell center)
module GRID (real space)
real(dp), public time_dtsec
time interval of model [sec]
Definition: scale_time.F90:36
integer, public js
start point of inner domain: y, local
module TIME
Definition: scale_time.F90:15
module PROCESS
integer, public ks
start point of inner domain: z, local
real(dp), public time_dtsec_atmos_phy_cp
time interval of physics(cumulus ) [sec]
Definition: scale_time.F90:40
module GRID (cartesian)
integer, public ie
end point of inner domain: x, local
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_cp_kf()

subroutine, public scale_atmos_phy_cp_kf::atmos_phy_cp_kf ( real(rp), dimension (ka,ia,ja), intent(in)  DENS,
real(rp), dimension (ka,ia,ja), intent(in)  MOMZ,
real(rp), dimension (ka,ia,ja), intent(in)  MOMX,
real(rp), dimension (ka,ia,ja), intent(in)  MOMY,
real(rp), dimension (ka,ia,ja), intent(in)  RHOT,
real(rp), dimension (ka,ia,ja,qa), intent(in)  QTRC,
real(rp), dimension (ka,ia,ja), intent(inout)  DENS_t_CP,
real(rp), dimension (ka,ia,ja), intent(inout)  MOMZ_t_CP,
real(rp), dimension (ka,ia,ja), intent(inout)  MOMX_t_CP,
real(rp), dimension (ka,ia,ja), intent(inout)  MOMY_t_CP,
real(rp), dimension (ka,ia,ja), intent(inout)  RHOT_t_CP,
real(rp), dimension (ka,ia,ja,qa_mp), intent(inout)  RHOQ_t_CP,
real(rp), dimension(ia,ja), intent(inout)  MFLX_cloudbase,
real(rp), dimension (ia,ja), intent(inout)  SFLX_convrain,
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,
real(rp), dimension (ia,ja), intent(inout)  nca,
real(rp), dimension (ka,ia,ja), intent(inout)  w0avg 
)

Definition at line 330 of file scale_atmos_phy_cp_kf.F90.

References cp_kf_main(), scale_stdio::io_fid_log, scale_stdio::io_l, kf_wmean(), scale_prof::prof_rapend(), scale_prof::prof_rapstart(), and scale_atmos_phy_mp::qa_mp.

Referenced by scale_atmos_phy_cp::atmos_phy_cp_setup().

330  use scale_grid_index
331  use scale_history, only: &
332  hist_in
333  use scale_atmos_phy_mp, only: &
334  qa_mp
335  implicit none
336 
337  real(RP), intent(in) :: DENS (KA,IA,JA)
338  real(RP), intent(in) :: MOMX (KA,IA,JA)
339  real(RP), intent(in) :: MOMY (KA,IA,JA)
340  real(RP), intent(in) :: MOMZ (KA,IA,JA)
341  real(RP), intent(in) :: RHOT (KA,IA,JA)
342  real(RP), intent(in) :: QTRC (KA,IA,JA,QA)
343  real(RP), intent(inout) :: DENS_t_CP (KA,IA,JA)
344  real(RP), intent(inout) :: MOMZ_t_CP (KA,IA,JA) ! not used
345  real(RP), intent(inout) :: MOMX_t_CP (KA,IA,JA) ! not used
346  real(RP), intent(inout) :: MOMY_t_CP (KA,IA,JA) ! not used
347  real(RP), intent(inout) :: RHOT_t_CP (KA,IA,JA)
348  real(RP), intent(inout) :: RHOQ_t_CP (KA,IA,JA,QA_MP)
349  real(RP), intent(inout) :: MFLX_cloudbase(IA,JA) ! not used
350  real(RP), intent(inout) :: SFLX_convrain (IA,JA) ! convective rain rate [kg/m2/s]
351  real(RP), intent(inout) :: cloudtop (IA,JA) ! cloud top height [m]
352  real(RP), intent(inout) :: cloudbase (IA,JA) ! cloud base height [m]
353  real(RP), intent(inout) :: cldfrac_dp (KA,IA,JA) ! cloud fraction (deep convection)
354  real(RP), intent(inout) :: cldfrac_sh (KA,IA,JA) ! cloud fraction (shallow convection)
355  real(RP), intent(inout) :: nca (IA,JA) ! convection active time [sec]
356  real(RP), intent(inout) :: w0avg (KA,IA,JA) ! running mean of vertical velocity [m/s]
357 
358  real(RP) :: cldfrac(KA,2) ! 1 shallow , 2 deep
359 
360  !---------------------------------------------------------------------------
361 
362  if( io_l ) write(io_fid_log,*) '*** Atmos physics step: Cumulus Parameterization(KF)'
363 
364  call kf_wmean( w0avg(:,:,:), & ! [OUT]
365  dens(:,:,:), & ! [IN]
366  momz(:,:,:) ) ! [IN]
367 
368  time_dokf = .false.
369  time_res_kf = time_res_kf + 1
370  if ( time_res_kf == time_dstep_kf ) then
371  time_dokf = .true.
372  time_res_kf = 0
373  endif
374 
375  if ( time_dokf ) then
376  if( io_l ) write(io_fid_log,*) '*** KF Convection Check '
377 
378  call prof_rapstart('CP_kf', 3)
379 
380  ! calc cumulus convection
381  call cp_kf_main( dens(:,:,:), & ! [IN]
382  momz(:,:,:), & ! [IN]
383  momx(:,:,:), & ! [IN]
384  momy(:,:,:), & ! [IN]
385  rhot(:,:,:), & ! [IN]
386  qtrc(:,:,:,:), & ! [IN]
387  w0avg(:,:,:), & ! [IN]
388  nca(:,:), & ! [INOUT]
389  dens_t_cp(:,:,:), & ! [INOUT]
390  rhot_t_cp(:,:,:), & ! [INOUT]
391  rhoq_t_cp(:,:,:,:), & ! [INOUT]
392  sflx_convrain(:,:), & ! [INOUT]
393  cldfrac_sh(:,:,:), & ! [INOUT]
394  cldfrac_dp(:,:,:), & ! [INOUT]
395  lifetime(:,:), & ! [INOUT]
396  cloudtop(:,:), & ! [INOUT]
397  cloudbase(:,:), & ! [INOUT]
398  i_convflag(:,:) ) ! [INOUT]
399 
400  call prof_rapend('CP_kf', 3)
401  endif
402 
403  call hist_in( lifetime(:,:), 'KF_LIFETIME', 'lifetime of KF scheme', 's' )
404  call hist_in( real(I_convflag(:,:),RP), 'KF_CONVFLAG', 'CONVECTION FLAG', '' )
405 
406  return
module ATMOSPHERE / Physics Cloud Microphysics
module grid index
module HISTORY
Here is the call graph for this function:
Here is the caller graph for this function:

◆ kf_wmean()

subroutine scale_atmos_phy_cp_kf::kf_wmean ( real(rp), dimension(ka,ia,ja), intent(inout)  W0_avg,
real(rp), dimension (ka,ia,ja), intent(in)  DENS,
real(rp), dimension (ka,ia,ja), intent(in)  MOMZ 
)

running mean vertical wind speed

Definition at line 421 of file scale_atmos_phy_cp_kf.F90.

References scale_grid_index::ie, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ke, kf_dtsec, scale_grid_index::ks, and scale_time::time_dtsec.

Referenced by atmos_phy_cp_kf().

421  use scale_time , only :&
422  time_dtsec
423  implicit none
424 
425  real(RP), intent(inout) :: W0_avg(KA,IA,JA)
426  real(RP), intent(in) :: DENS (KA,IA,JA)
427  real(RP), intent(in) :: MOMZ (KA,IA,JA)
428 
429  real(RP) :: W0
430  real(RP) :: fact_W0_avg, fact_W0
431 
432  integer :: k, i, j
433  !---------------------------------------------------------------------------
434 
435  if ( param_atmos_phy_cp_kf_wadapt ) then
436  fact_w0_avg = 2.0_rp * max(kf_dtsec,time_dtsec) - time_dtsec
437  fact_w0 = time_dtsec
438  else ! w_time is tuning parameter
439  fact_w0_avg = real(PARAM_ATMOS_PHY_CP_kf_w_time,RP)
440  fact_w0 = 1.0_rp
441  endif
442 
443  do j = js, je
444  do i = is, ie
445  do k = ks, ke
446  w0 = 0.5_rp * ( momz(k,i,j) + momz(k-1,i,j) ) / dens(k,i,j)
447 
448  w0_avg(k,i,j) = ( w0_avg(k,i,j) * fact_w0_avg &
449  + w0 * fact_w0 ) / ( fact_w0_avg + fact_w0 )
450  enddo
451  enddo
452  enddo
453 
454  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
integer, public ke
end point of inner domain: z, local
real(dp), public time_dtsec
time interval of model [sec]
Definition: scale_time.F90:36
integer, public js
start point of inner domain: y, local
module TIME
Definition: scale_time.F90:15
integer, public ks
start point of inner domain: z, local
integer, public ie
end point of inner domain: x, local
Here is the caller graph for this function:

◆ cp_kf_param()

subroutine scale_atmos_phy_cp_kf::cp_kf_param ( real(rp), intent(in)  RATE_in,
integer, intent(inout)  TRIGGER_in,
real(dp), intent(in)  KF_DT_in,
real(rp), intent(in)  DELCAPE_in,
real(rp), intent(in)  DEEPLIFETIME_in,
real(rp), intent(in)  SHALLOWLIFETIME_in,
real(rp), intent(in)  DEPTH_USL_in,
integer, intent(in)  KF_prec_in,
real(rp), intent(in)  KF_threshold_in,
logical, intent(in)  KF_LOG_in,
integer, intent(in)  stepkf_in 
)

Definition at line 471 of file scale_atmos_phy_cp_kf.F90.

References scale_stdio::io_fid_log, scale_stdio::io_l, kf_dtsec, and scale_process::prc_mpistop().

Referenced by atmos_phy_cp_kf_setup().

471  use scale_process, only: &
473  implicit none
474  real(RP),intent(in) :: RATE_in
475  integer, intent(inout) :: TRIGGER_in
476  real(DP),intent(in) :: KF_DT_in
477  real(RP),intent(in) :: DELCAPE_in
478  real(RP),intent(in) :: DEEPLIFETIME_in
479  real(RP),intent(in) :: SHALLOWLIFETIME_in
480  real(RP),intent(in) :: DEPTH_USL_in
481  integer, intent(in) :: KF_prec_in
482  real(RP),intent(in) :: KF_threshold_in
483  logical, intent(in) :: KF_LOG_in
484  integer, intent(in) :: stepkf_in
485  !
486  rate = rate_in
487  ! TRIGGER must be 1 or 3
488  if (trigger_in /= 1 .and. trigger_in /=3) then
489  if( io_l ) write(io_fid_log,*) "TRIGGER must be 1 or 3 but now :",trigger_in
490  if( io_l ) write(io_fid_log,*) "CHAGNGE ",trigger_in," to 3"
491  trigger_in = 3
492  end if
493  trigger = trigger_in
494  kf_dtsec = kf_dt_in * 60.0_dp ! convert from min to sec
495  delcape = delcape_in
496  deeplifetime = deeplifetime_in
497  shallowlifetime = shallowlifetime_in
498  depth_usl = depth_usl_in
499  kf_prec = kf_prec_in
500  kf_threshold = kf_threshold_in
501  kf_log = kf_log_in
502  stepkf = stepkf_in
503  if (kf_prec == 1) then
504  p_precipitation => precipitation_oc1973 ! Ogura and Cho (1973)
505  elseif( kf_prec == 2) then
506  p_precipitation => precipitation_kessler ! Kessler type
507  else
508  write(*,*) 'xxx ERROR at KF namelist'
509  write(*,*) 'KF_prec must be 1 or 2'
510  call prc_mpistop
511  end if
512  return
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
Here is the call graph for this function:
Here is the caller graph for this function:

◆ cp_kf_main()

subroutine scale_atmos_phy_cp_kf::cp_kf_main ( real(rp), dimension(ka,ia,ja), intent(in)  dens,
real(rp), dimension(ka,ia,ja), intent(in)  MOMZ,
real(rp), dimension(ka,ia,ja), intent(in)  MOMX,
real(rp), dimension(ka,ia,ja), intent(in)  MOMY,
real(rp), dimension(ka,ia,ja), intent(in)  RHOT,
real(rp), dimension(ka,ia,ja,qa), intent(in)  QTRC_in,
real(rp), dimension(ka,ia,ja), intent(in)  w0avg,
real(rp), dimension(ia,ja), intent(inout)  nca,
real(rp), dimension(ka,ia,ja), intent(inout)  DT_DENS,
real(rp), dimension(ka,ia,ja), intent(inout)  DT_RHOT,
real(rp), dimension(ka,ia,ja,qs_mp:qe_mp), intent(inout)  DT_RHOQ,
real(rp), dimension(ia,ja), intent(inout)  rainrate_cp,
real(rp), dimension(ka,ia,ja), intent(inout)  cldfrac_sh,
real(rp), dimension(ka,ia,ja), intent(inout)  cldfrac_dp,
real(rp), dimension(ia,ja), intent(inout)  timecp,
real(rp), dimension(ia,ja), intent(inout)  cloudtop,
real(rp), dimension(ia,ja), intent(inout)  zlcl,
integer, dimension(ia,ja), intent(inout)  I_convflag 
)

I_convflag = 0 ==> deep convection = 1 ==> shallow convection = 2 ==> NONE !!

Definition at line 536 of file scale_atmos_phy_cp_kf.F90.

References scale_const::const_cpdry, scale_const::const_emelt, scale_const::const_grav, scale_const::const_pre00, scale_const::const_rdry, scale_const::const_tem00, scale_atmos_hydrometeor::i_qc, scale_atmos_hydrometeor::i_qi, scale_atmos_hydrometeor::i_qr, scale_atmos_hydrometeor::i_qs, scale_atmos_hydrometeor::i_qv, scale_grid_index::ie, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::is, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ke, kf_dtsec, scale_grid_index::ks, scale_process::prc_mpistop(), scale_atmos_phy_mp::qa_mp, scale_atmos_phy_mp::qe_mp, scale_atmos_phy_mp::qs_mp, scale_grid_real::real_fz, scale_time::time_dtsec_atmos_phy_cp, scale_tracer::tracer_cv, scale_tracer::tracer_mass, and scale_tracer::tracer_r.

Referenced by atmos_phy_cp_kf().

536  use scale_precision
537  use scale_grid_index
538  use scale_tracer
539  use scale_const, only: &
540  grav => const_grav, &
541  r => const_rdry
542  use scale_atmos_phy_mp, only: &
543  qa_mp, &
544  qs_mp, &
545  qe_mp
546  use scale_atmos_hydrometeor, only: &
547  i_qv, &
548  i_qc, &
549  i_qr, &
550  i_qi, &
551  i_qs
552  use scale_time , only :&
554  use scale_grid_real, only: &
555  fz => real_fz
556  use scale_atmos_thermodyn, only: &
557  thermodyn_temp_pres => atmos_thermodyn_temp_pres, &
558  thermodyn_rhoe => atmos_thermodyn_rhoe, &
559  thermodyn_temp_pres_e => atmos_thermodyn_temp_pres_e, &
560  thermodyn_qd => atmos_thermodyn_qd, &
561  thermodyn_pott => atmos_thermodyn_pott
562  use scale_atmos_saturation ,only :&
563  saturation_psat_liq => atmos_saturation_psat_liq
564  implicit none
565  ![IN]
566  real(RP),intent(in) :: DENS(KA,IA,JA) ! density [kg/m**3]
567  real(RP),intent(in) :: MOMZ(KA,IA,JA) ! momentum
568  real(RP),intent(in) :: MOMX(KA,IA,JA) ! momentum
569  real(RP),intent(in) :: MOMY(KA,IA,JA) ! momentum
570  real(RP),intent(in) :: RHOT(KA,IA,JA) ! density*PT
571  real(RP),intent(in) :: QTRC_in(KA,IA,JA,QA) ! raito of water elements
572  real(RP),intent(in) :: w0avg(KA,IA,JA) ! running mean vertical wind velocity [m/s]
573  ! [INOUT]
574  real(RP),intent(inout) :: nca(IA,JA) ! num of step convection active [step]
575  real(RP),intent(inout) :: DT_DENS(KA,IA,JA) ! dens/dt
576  real(RP),intent(inout) :: DT_RHOT(KA,IA,JA) ! dens*PT/dt
577  real(RP),intent(inout) :: DT_RHOQ(KA,IA,JA,QS_MP:QE_MP) ! dens*q/dt
578  real(RP),intent(inout) :: rainrate_cp(IA,JA) ! rain PPTFLX(prcp_flux)/deltax**2*dt ! convective rain
579  real(RP),intent(inout) :: cldfrac_sh(KA,IA,JA) ! cloud fraction
580  real(RP),intent(inout) :: cldfrac_dp(KA,IA,JA) ! cloud fraction
581  real(RP),intent(inout) :: timecp(IA,JA) ! timescale of cumulus parameterization
582  real(RP),intent(inout) :: cloudtop(IA,JA) ! cloud top height
583  real(RP),intent(inout) :: zlcl(IA,JA) ! hight of lcl cloud bottom height[m]
584  integer, intent(inout) :: I_convflag(IA,JA) ! convection type
588 
589  ! [Internal Work]
590  real(RP) :: u (KA) ! x-direction wind velocity [m/s]
591  real(RP) :: v (KA) ! y-direction wind velocity [m/s]
592  real(RP) :: temp (KA) ! temperature [K]
593  real(RP) :: pres (KA) ! pressure of dry air [Pa]
594  real(RP) :: qv (KA) ! water vapor mixing ratio[kg/kg] original in KF scheme
595  real(RP) :: QDRY (KA) ! ratio of dry air
596  !real(RP) :: qes (KA) ! saturate vapor [kg/kg]
597  real(RP) :: PSAT (KA) ! saturation vaper pressure
598  real(RP) :: QSAT (KA) ! saturate water vaper mixing ratio [kg/kg]
599  real(RP) :: rh (KA) ! saturate vapor [%]
600  real(RP) :: deltap(KA) ! delta Pressure [Pa]
601 
602  real(RP) :: q_hyd(KA,QA_MP-1) ! water mixing ratio [kg/kg]
603  real(RP) :: dens_nw(KA) ! density [kg/m**3]
604  integer :: nic
605  real(RP) :: time_advec ! advection timescale
606  real(RP) :: umf(KA) ! Updraft Mass Flux
607  real(RP) :: umflcl ! Updraft Mass Flux
608  real(RP) :: upent(KA) ! Updraft Mass flux ENTrainment
609  real(RP) :: updet(KA) ! Updraft Mass flux DETrainment
610  ! updraft value
611  real(RP) :: temp_u(KA) ! updraft temperature [K]
612  real(RP) :: tempv(KA) ! vertual temperature [K]
613  real(RP) :: qv_u(KA) ! updraft qv
614  real(RP) :: cape ! cape
615  ! water
616  real(RP) :: qc(KA) ! cloud water
617  real(RP) :: qi(KA) ! cloud ice
618  real(RP) :: qvdet(KA) ! vapor detrainment
619  real(RP) :: qcdet(KA) ! liquit detrainment
620  real(RP) :: qidet(KA) ! ice detrainment
621  real(RP) :: totalprcp ! total precipitation
622  real(RP) :: flux_qr(KA) ! rain flux
623  real(RP) :: flux_qs(KA) ! snow flux
624  ! upward theta
625  real(RP) :: theta_eu(KA) ! updraft equivalent PT
626  real(RP) :: theta_ee(KA) ! environment equivalent PT
627  ! convection type
628  ! index valuables
629  integer :: k_lcl ! index of LCL layer
630  integer :: k_top ! index of cloud top hight
631  integer :: k_ml ! index of melt layer (temp < tem00)
632  integer :: k_lc,k_let,k_pbl ! indexs
633  real(RP) :: zmix ! usl(up source layer) layer depth [m]
634  real(RP) :: presmix ! usl layer depth [Pa]
635  real(RP) :: umfnewdold(KA) ! umfnew/umfold
636  real(RP) :: dpthmx ! max depth of pressure
637  ! move internalwark
638  real(RP) :: ems(KA) ! box weight[kg]
639  real(RP) :: emsd(KA) ! 1/ems
640  ! [OUT]@downdraft
641  real(RP) :: wspd(3) ! horizontal wind speed 1 k_lcl, 2 k_z5,3 k_top
642  real(RP) :: dmf(KA) ! downdraft massflux
643  real(RP) :: downent(KA) ! downdraft entrainment
644  real(RP) :: downdet(KA) ! downdraft detrainment
645  real(RP) :: theta_d(KA) ! downdraft theta
646  real(RP) :: qv_d(KA) ! water vapor of downdraft
647  real(RP) :: prcp_flux ! precpitation flux
648  real(RP) :: tder ! temperature change from evap
649  real(RP) :: CPR ! all precipitation before consider cloud bottom evaporation
650  integer :: k_lfs ! LFS layer index
651  ![OUT] @ compensasional subsidence
652  real(RP) :: temp_g(KA) ! temporarly temperature -> after iteration then new variable
653  real(RP) :: qv_nw(KA) ! new vapor mixing ratio [kg/kg]
654  real(RP) :: qc_nw(KA) ! new cloud water mixing ratio [kg/kg]
655  real(RP) :: qi_nw(KA) ! new cloud ice mixing ratio [kg/kg]
656  real(RP) :: qr_nw(KA) ! new rain water mixing ratio [kg/kg]
657  real(RP) :: qs_nw(KA) ! new snow water mixing ratio [kg/kg]
658  ! new variable
659  real(RP) :: qtrc_nw(KA,QA) ! qv,qc,qr,qi,qs (qg not change)
660  real(RP) :: pott_nw(KA) ! new PT
661  !
662  real(RP) :: cldfrac_KF(KA,2) ! cloud fraction
663  !
664  real(RP) :: RHOD(KA) ! dry density
665  real(RP) :: QV_resd(KA) ! residual vapor
666 
667  ! do loop variable
668  integer :: k, i, j, iq, iqa, ii
669 
670  i_convflag(:,:) = 2
671 
672  do j = js, je
673  do i = is, ie
674 
675  nca(i,j) = nca(i,j) - real(TIME_DSTEP_KF,RP) * dt
676 
677  ! check convection
678  if ( nca(i,j) .ge. 0.5_dp * dt ) cycle
679 
680  do k = ks, ke
681  ! preparing a NON Hydriometeor condition to fit assumption in KF scheme
682  call thermodyn_temp_pres( temp(k), & ! [OUT]
683  pres(k), & ! [OUT] !dummy
684  dens(k,i,j), & ! [IN]
685  rhot(k,i,j), & ! [IN]
686  qtrc_in(k,i,j,:), & ! [IN]
687  tracer_cv(:), & ! [IN]
688  tracer_r(:), & ! [IN]
689  tracer_mass(:) ) ! [IN]
690 
691  call thermodyn_qd( qdry(k), qtrc_in(k,i,j,:), tracer_mass(:) )
692  rhod(k) = dens(k,i,j) * qdry(k)
693 
694  ! calculate u(x-directin velocity ), v(y-direction velocity)
695  u(k) = 0.5_rp * ( momx(k,i,j) + momx(k,i-1,j) ) / dens(k,i,j)
696  v(k) = 0.5_rp * ( momy(k,i,j) + momy(k,i,j-1) ) / dens(k,i,j)
697  enddo
698 
699  ! initialize variables
700  cloudtop(i,j) = 0.0_rp
701  zlcl(i,j) = 0.0_rp
702  rainrate_cp(i,j)= 0.0_rp
703  timecp(i,j) = 0.0_rp
704  cldfrac_kf(:,:) = 0.0_rp
705  tempv(:) = 0.0_rp
706  qc(:) = 0.0_rp
707  qi(:) = 0.0_rp
708  flux_qr(:) = 0.0_rp
709  flux_qs(:) = 0.0_rp
710  theta_eu(:) = 0.0_rp
711  theta_ee(:) = 0.0_rp
712  theta_d(:) = 0.0_rp
713  umfnewdold(:) = 0.0_rp
714  wspd(:) = 0.0_rp
715  temp_g(:) = 0.0_rp
716  qv_nw(:) = 0.0_rp
717  qc_nw(:) = 0.0_rp
718  qi_nw(:) = 0.0_rp
719  qr_nw(:) = 0.0_rp
720  qs_nw(:) = 0.0_rp
721  pott_nw(:) = 0.0_rp
722  totalprcp = 0.0_rp
723  umflcl = 0.0_rp
724  cape = 0.0_rp
725  presmix = 0.0_rp
726  dpthmx = 0.0_rp
727  zmix = 0.0_rp
728  prcp_flux = 0.0_rp
729  cpr = 0.0_rp
730  tder = 0.0_rp
731  time_advec = 0.0_rp
732  k_lcl = 0
733  k_lc = 0
734  k_lfs = 0
735  k_pbl = 0
736  k_top = 0
737  k_let = 0
738  k_ml = 0
739  nic = 0
740 
741  do k = ks, ke
742  ! temporary: WRF TYPE equations are used to maintain consistency with kf_main
743  !call SATURATION_psat_liq( PSAT(k), TEMP(k) )
744  !QSAT(k) = 0.622_RP * PSAT(k) / ( PRES(k) - ( 1.0_RP-0.622_RP ) * PSAT(k) )
745  psat(k) = aliq*exp((bliq*temp(k)-cliq)/(temp(k)-dliq))
746  qsat(k) = 0.622_rp * psat(k) / ( pres(k) - psat(k) )
747 
748  ! calculate water vaper and relative humidity
749 ! QV (k) = max( 0.000001_RP, min( QSAT(k), QV(k) ) ) ! conpare QSAT and QV, guess lower limit
750  qv(k) = max( kf_eps, min( qsat(k), qtrc_in(k,i,j,i_qv) / qdry(k) ) ) ! conpare QSAT and QV, guess lower limit
751  rh(k) = qv(k) / qsat(k)
752  enddo
753 
754  ! calculate delta P by hydrostatic balance
755  ! deltap is the pressure interval between half levels(face levels) @ SCALE
756  do k = ks, ke
757  deltap(k) = rhod(k) * grav * ( fz(k+1,i,j) - fz(k,i,j) ) ! rho*g*dz
758  enddo
759 
760  cldfrac_kf(ks:ke,:) = 0.0_rp
761  rainrate_cp(i,j) = 0.0_rp
762  timecp(i,j) = 0.0_rp
763  cloudtop(i,j) = 0.0_rp
764  zlcl(i,j) = 0.0_rp
765 
766  call cp_kf_trigger ( &
767  ! [IN]
768  deltaz(:,i,j), z(:,i,j), & ! deltaz and height[m]
769  qv(:), & ! water vapor mixing ratio [kg/kg]
770  qsat(:), & ! saturation water vapor mixing ratio
771  pres(:), & ! pressure [Pa]
772  deltap(:), & ! interval of pressure [Pa]
773  temp(:), & ! temperature [K]
774  w0avg(:,i,j), & ! average w
775  ! [OUT]
776  i_convflag(i,j), & ! convection flag
777  cloudtop(i,j), & ! cloud top height [m]
778  temp_u(:), & ! updraft temperature
779  tempv(:), & ! virtual temperature
780  ! water values
781  qv_u(:), & ! updraft water vapor mixing ratio
782  qc(:), & ! cloud water mixing ratio
783  qi(:), & ! ice mixing ratio
784  qvdet(:), & ! updraft detrain vapor
785  qcdet(:), & ! updraft detrain water
786  qidet(:), & ! updraft detrain ice
787  flux_qr(:), & ! rain flux
788  flux_qs(:), & ! snow flux
789  totalprcp, & ! total precipitation(rain and snow ) fall
790  ! thermodynamics values
791  theta_eu(:), & ! updraft equivalent PT
792  theta_ee(:), & ! environment equivalent PT
793  cape, & ! cape
794  ! massflux values
795  umf(:), & ! updraft mass flux
796  umflcl, & ! updraft mass flux@lcl
797  upent(:), & ! updraft entrainment
798  updet(:), & ! updraft detrainment
799  ! layer nambers
800  k_lcl, & ! index of LCL layer
801  k_lc, & ! index of USL layer botom
802  k_pbl, & ! index of USL layer top
803  k_top, & ! index of cloud top
804  k_let, & ! index of below detrain layer
805  k_ml, & ! index of temprerure < tem00
806  ! pbl layer values
807  presmix, & ! USL layer pressure [Pa]
808  dpthmx, & ! USL layer depth ( [Pa])
809  zlcl(i,j), & ! lcl hight
810  zmix, & ! USL layer depth
811  umfnewdold(:) )
812  !------------------------------------------------------------------------
813  ! calc ems(box weight[kg])
814  !------------------------------------------------------------------------
815 
816  if(i_convflag(i,j) /= 2) then
817  ems(k_top+1:ke) = 0._rp
818  emsd(k_top+1:ke) = 0._rp
819  do k = ks, k_top
820  ems(k) = deltap(k)*deltax*deltax/grav
821  emsd(k) = 1._rp/ems(k)
822  umfnewdold(k) = 1._rp/umfnewdold(k)
823  end do
824  end if
825  !------------------------------------------------------------------------
826  call cp_kf_downdraft ( &
827  ! [IN]
828  deltaz(:,i,j), z(:,i,j) ,& ! deltaz and height [m]
829  u(:) ,& ! u-wind m/s
830  v(:) ,& ! v-wind m/s
831  zlcl(i,j) ,& ! lcl height [m]
832  rh(:) ,& ! relative humidity
833  deltap(:) ,& ! interval of pressure [Pa]
834  pres(:) ,& ! pressure [Pa]
835  qv(:) ,& ! water vapor mixing ratio [kg/kg]
836  ems(:) ,& ! ems(box weight[kg])
837  emsd(:) ,& ! 1/ems
838  theta_ee(:) ,& ! environment theta_E
839  umf(:) ,& ! updraft mass flux
840  totalprcp ,& ! total precipitation
841  flux_qs(:) ,& ! snow flux
842  tempv(:) ,& ! vertual temperature [K]
843  i_convflag(i,j) ,& ! convection type
844  ! layer index
845  k_lcl ,& ! index of LCL layer
846  k_ml ,& ! index of melt layer
847  k_top ,& ! index of cloud top
848  k_pbl ,& ! index of USL layer top
849  k_let ,& ! index of below detrain layer
850  k_lc ,& ! index of USL layer botom
851  ! [out]
852  wspd(:) ,& ! wind speed 1 k_lcl, 2 k_z5,3 k_top
853  dmf(:) ,& ! downward mass flux
854  downent(:) ,& ! downdraft entrainment
855  downdet(:) ,& ! downdraft detrainment
856  theta_d(:) ,& ! downdraft theta
857  qv_d(:) ,& ! downdraft water vaper mixng ratio
858  prcp_flux ,& ! precipitateion
859  k_lfs ,& ! index of LFS layer
860  cpr ,& ! all precipitation before consider cloud bottom evaporation
861  tder)
862  !------------------------------------------------------------------------
863  call cp_kf_compensational ( &
864  ![IN]
865  deltaz(:,i,j),z(:,i,j) ,& ! deltaz and height [m]
866  pres(:), & ! pressure [Pa]
867  deltap(:), & ! deltap [Pa]
868  temp(:), & ! temperature [K]
869  qv(:), & ! water vapor mixing ratio
870  ! form kf_trigger
871  ems(:), & ! box weight[kg]
872  emsd(:), & ! 1/ems
873  presmix, & ! usl layer pressure depth
874  zmix, & ! usl layer pressure depth
875  dpthmx, & ! usl layer depth
876  cape, & ! CAPE
877  temp_u(:), & ! updraft temperature
878  qvdet(:), & ! detrainment water vaper mixing ratio
879  umflcl, & ! umf @LCL
880  umf(:), & ! upward mass flux (updraft)
881  upent(:), & ! updraft entrainment
882  updet(:), & ! updraft detrainment
883  qcdet(:), & ! updraft detrainment qc
884  qidet(:), & ! updraft detrainment qi
885  qc(:), & ! cloud water mixing ratio
886  qi(:), & ! cloud ice mixing ratio
887  flux_qr(:), & ! rain flux
888  flux_qs(:), & ! snow flux
889  umfnewdold(:), & ! umfnew/umfold ratio
890  ! from kf_downdraft
891  wspd(:), & ! wind speed 1. 2 is used
892  dmf(:), & ! downward mass flux (downdraft)
893  downent(:), & ! downdraft entrainment
894  downdet(:), & ! downdraft detrainment
895  qv_d(:), & ! downdraft detrainment qv
896  theta_d(:), & ! potential temperature @downdraft
897  prcp_flux, & ! precipitation
898  tder, & ! evapolation
899  cpr, & ! all precipitation before consider cloud bottom evaporation
900  i_convflag(i,j), & ! intent inout convective flag
901  ! layer index
902  ! from kf_triger
903  k_top, & ! index of cloud top hight
904  k_lcl, & ! index of LCL layer
905  k_lc , & ! index of LCL layer botoom
906  k_pbl, & ! index of USL layer top
907  k_ml, & ! index of melt layer (temp < tem00)
908  ! from kf_downdraft
909  k_lfs, & ! index of LFS layer
910  ![OUT] new valuavles after timestep
911  temp_g(:), & ! update temperature [K]
912  qv_nw(:), & ! update water vapor mixing ratio
913  qc_nw(:), & ! update cloud water mixing ratio
914  qi_nw(:), & ! update cloud ice mixing ratio
915  qr_nw(:), & ! update rain water mixing ratio
916  qs_nw(:), & ! update snow water mixing ratio
917  rainrate_cp(i,j), & ! rain PPTFLX(prcp_flux)/deltax**2*dt ! convective rain
918  nic, & ! (timescale of cumulus parameterization)/dt integer
919  cldfrac_kf, & ! cloud fraction
920  timecp(i,j), & ! timescale of cumulus parameterization
921  time_advec) ! advection timescale
922 
923  !------------------------------------------------------------------------
924  ! compute tendencys
925  !------------------------------------------------------------------------
926 
927  if(i_convflag(i,j) == 2) then ! no convection
928  rainrate_cp(i,j) = 0.0_rp
929  cldfrac_kf(ks:ke,:) = 0.0_rp
930  timecp(i,j) = 0.0_rp
931  cloudtop(i,j) = 0.0_rp
932  zlcl(i,j) = 0.0_rp
933  nca(i,j) = -100.0_rp
934 
935  dt_dens(ks:ke,i,j) = 0.0_rp
936  dt_rhot(ks:ke,i,j) = 0.0_rp
937  do iq = qs_mp, qe_mp
938  dt_rhoq(ks:ke,i,j,iq) = 0.0_rp
939  end do
940 
941  else ! convection allowed I_convflag=0 or 1
942  ! chek
943  !
944  !...FEEDBACK TO RESOLVABLE SCALE TENDENCIES.
945  !
946  !...IF THE ADVECTIVE TIME PERIOD (time_advec) IS LESS THAN SPECIFIED MINIMUM
947  !...TIMECP, ALLOW FEEDBACK TO OCCUR ONLY DURING TADVEC...
948  !
949  if (i_convflag(i,j) == 0) then ! deep
950  if (time_advec < timecp(i,j)) nic=nint(time_advec/dt)
951  nca(i,j) = real(nic,RP)*dt ! convection feed back act this time span
952  elseif (i_convflag(i,j) == 1) then ! shallow
953  timecp(i,j) = shallowlifetime
954  nca(i,j) = kf_dtsec ! convection feed back act this time span
955  end if
956 
957  dt_rhoq(:,i,j,:) = 0.0_rp
958  do k=ks, k_top
959  dt_rhoq(k,i,j,i_qv) = ( rhod(k) * ( qv_nw(k) - qv(k) ) ) / timecp(i,j)
960  dt_rhoq(k,i,j,i_qc) = qc_nw(k) * rhod(k) / timecp(i,j)
961  dt_rhoq(k,i,j,i_qr) = qr_nw(k) * rhod(k) / timecp(i,j)
962  dt_dens(k,i,j) = dt_rhoq(k,i,j,i_qv) + dt_rhoq(k,i,j,i_qc) + dt_rhoq(k,i,j,i_qr)
963  if ( i_qi>0 ) then
964  dt_rhoq(k,i,j,i_qi) = qi_nw(k) * rhod(k) / timecp(i,j)
965  dt_dens(k,i,j) = dt_dens(k,i,j) + dt_rhoq(k,i,j,i_qi)
966  end if
967  if ( i_qs>0 ) then
968  dt_rhoq(k,i,j,i_qs) = qs_nw(k) * rhod(k) / timecp(i,j)
969  dt_dens(k,i,j) = dt_dens(k,i,j) + dt_rhoq(k,i,j,i_qs)
970  end if
971 
972  dens_nw(k) = dens(k,i,j) + dt_dens(k,i,j) * timecp(i,j)
973 
974  do iq = qs_mp, qe_mp
975  qtrc_nw(k,iq) = ( qtrc_in(k,i,j,iq) * dens(k,i,j) + dt_rhoq(k,i,j,iq) * timecp(i,j) ) / dens_nw(k)
976  end do
977  do iq = 1, qs_mp - 1
978  qtrc_nw(k,iq) = qtrc_in(k,i,j,iq) * dens(k,i,j) / dens_nw(k)
979  end do
980  do iq = qe_mp + 1, qa
981  qtrc_nw(k,iq) = qtrc_in(k,i,j,iq) * dens(k,i,j) / dens_nw(k)
982  end do
983  end do
984 
985  ! treatment for not evaluated layers
986  do k=k_top+1, ke
987  do iq = 1, qa
988  qtrc_nw(k,iq) = qtrc_in(k,i,j,iq)
989  end do
990  dens_nw(k) = dens(k,i,j)
991  enddo
992 
993  ! calc new potential temperature
994  call thermodyn_pott(&
995  pott_nw(:), &
996  temp_g(:), pres(:), qtrc_nw(:,:), &
997  tracer_cv(:), tracer_r(:), tracer_mass(:) )
998  ! update rhot
999  do k = ks, k_top
1000  dt_rhot(k,i,j) = ( dens_nw(k)*pott_nw(k) - rhot(k,i,j) ) / timecp(i,j)
1001  end do
1002 
1003  do k=k_top+1, ke
1004  dt_dens(k,i,j) = 0.0_rp
1005  dt_rhot(k,i,j) = 0.0_rp
1006  do iq = qs_mp, qe_mp
1007  dt_rhoq(k,i,j,iq) = 0.0_rp
1008  end do
1009  end do
1010 
1011  ! to keep conservation
1012  ! if noconvection then nca is same value before call. nca only modifyed convectioned
1013  end if
1014 
1015  cldfrac_sh(ks:ke,i,j) = cldfrac_kf(ks:ke,1)
1016  cldfrac_dp(ks:ke,i,j) = cldfrac_kf(ks:ke,2)
1017 
1018  end do
1019  end do
1020 
1021  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
module ATMOSPHERE / Saturation adjustment
real(rp), dimension(qa_max), public tracer_r
module ATMOSPHERE / Physics Cloud Microphysics
integer, public ke
end point of inner domain: z, local
real(rp), dimension(:,:,:), allocatable, public real_fz
geopotential height [m] (cell face )
real(rp), dimension(qa_max), public tracer_cv
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
module grid index
module TRACER
module GRID (real space)
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:48
integer, public js
start point of inner domain: y, local
module TIME
Definition: scale_time.F90:15
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
real(dp), public time_dtsec_atmos_phy_cp
time interval of physics(cumulus ) [sec]
Definition: scale_time.F90:40
integer, public ie
end point of inner domain: x, local
module ATMOSPHERE / Thermodynamics
module PRECISION
real(rp), dimension(qa_max), public tracer_mass
Here is the call graph for this function:
Here is the caller graph for this function:

◆ kf_lutab()

subroutine scale_atmos_phy_cp_kf::kf_lutab ( )

Definition at line 3642 of file scale_atmos_phy_cp_kf.F90.

References scale_const::const_cpdry, scale_const::const_grav, and scale_const::const_pre00.

Referenced by atmos_phy_cp_kf_setup().

3642  !
3643  ! This subroutine is a lookup table.
3644  ! Given a series of series of saturation equivalent potential
3645  ! temperatures, the temperature is calculated.
3646  !
3647  !--------------------------------------------------------------------
3648  use scale_const, only :&
3649  cp => const_cpdry , &
3650  pre00 => const_pre00, &
3651  grav => const_grav
3652  IMPLICIT NONE
3653  !--------------------------------------------------------------------
3654  ! Lookup table variables
3655  ! INTEGER, SAVE, PARAMETER :: KFNT=250,KFNP=220
3656  ! REAL, SAVE, DIMENSION(1:KFNT,1:KFNP) :: TTAB,QSTAB
3657  ! REAL, SAVE, DIMENSION(1:KFNP) :: THE0K
3658  ! REAL, SAVE, DIMENSION(1:200) :: ALU
3659  ! REAL, SAVE :: RDPR,RDTHK,PLUTOP
3660  ! End of Lookup table variables
3661  !! use scale_const, only : &
3662  !! SVPT0 -> CONST_TEM00,&
3663  integer :: KP,IT,ITCNT,I
3664  real(RP) :: DTH=1._rp,tmin=150._rp,toler=0.001_rp
3665  real(RP) :: PBOT,DPR, &
3666  TEMP,P,ES,QS,PI,THES,TGUES,THGUES,F0,T1,T0,THGS,F1,DT, &
3667  ASTRT,AINC,A1,THTGS
3668  ! REAL :: ALIQ,BLIQ,CLIQ,DLIQ,SVP1,SVP2,SVP3,SVPT0
3669  !! to module variables
3670  ! real(RP) :: ALIQ,BLIQ,CLIQ,DLIQ
3671  !! real(RP), INTENT(IN) :: SVP1,SVP2,SVP3,SVPT0
3672  !
3673  ! equivalent potential temperature increment
3674  ! data dth/1._RP/
3675  ! minimum starting temp
3676  ! data tmin/150._RP/
3677  ! tolerance for accuracy of temperature
3678  ! data toler/0.001_RP/
3679  ! top pressure (pascals)
3680  plutop=5000.0_rp
3681  ! bottom pressure (pascals)
3682  pbot=110000.0_rp
3683  !! to module variable
3684  !! ALIQ = SVP1*1000.
3685  !! BLIQ = SVP2
3686  !! CLIQ = SVP2*SVPT0
3687  !! DLIQ = SVP3
3688 
3689  !
3690  ! compute parameters
3691  !
3692  ! 1._over_(sat. equiv. theta increment)
3693  rdthk=1._rp/dth
3694  ! pressure increment
3695  !
3696  dpr=(pbot-plutop)/REAL(KFNP-1)
3697  ! dpr=(pbot-plutop)/REAL(kfnp-1)
3698  ! 1._over_(pressure increment)
3699  rdpr=1._rp/dpr
3700  ! compute the spread of thes
3701  ! thespd=dth*(kfnt-1)
3702  !
3703  ! calculate the starting sat. equiv. theta
3704  !
3705  temp=tmin
3706  p=plutop-dpr
3707  do kp=1,kfnp
3708  p=p+dpr
3709  es=aliq*exp((bliq*temp-cliq)/(temp-dliq))
3710  qs=0.622_rp*es/(p-es)
3711  pi=(1.e5_rp/p)**(0.2854_rp*(1.-0.28_rp*qs))
3712  the0k(kp)=temp*pi*exp((3374.6525_rp/temp-2.5403_rp)*qs* &
3713  (1._rp+0.81_rp*qs))
3714  enddo
3715  !
3716  ! compute temperatures for each sat. equiv. potential temp.
3717  !
3718  p=plutop-dpr
3719  do kp=1,kfnp
3720  thes=the0k(kp)-dth
3721  p=p+dpr
3722  do it=1,kfnt
3723  ! define sat. equiv. pot. temp.
3724  thes=thes+dth
3725  ! iterate to find temperature
3726  ! find initial guess
3727  if(it.eq.1) then
3728  tgues=tmin
3729  else
3730  tgues=ttab(it-1,kp)
3731  endif
3732  es=aliq*exp((bliq*tgues-cliq)/(tgues-dliq))
3733  qs=0.622_rp*es/(p-es)
3734  pi=(1.e5_rp/p)**(0.2854_rp*(1._rp-0.28_rp*qs))
3735  thgues=tgues*pi*exp((3374.6525_rp/tgues-2.5403_rp)*qs* &
3736  (1._rp + 0.81_rp*qs))
3737  f0=thgues-thes
3738  t1=tgues-0.5_rp*f0
3739  t0=tgues
3740  itcnt=0
3741  ! iteration loop
3742  do itcnt=1,11
3743  es=aliq*exp((bliq*t1-cliq)/(t1-dliq))
3744  qs=0.622_rp*es/(p-es)
3745  pi=(1.e5_rp/p)**(0.2854_rp*(1._rp-0.28_rp*qs))
3746  thtgs=t1*pi*exp((3374.6525_rp/t1-2.5403_rp)*qs*(1._rp + 0.81_rp*qs))
3747  f1=thtgs-thes
3748  if(abs(f1).lt.toler)then
3749  exit
3750  endif
3751  ! itcnt=itcnt+1
3752  dt=f1*(t1-t0)/(f1-f0)
3753  t0=t1
3754  f0=f1
3755  t1=t1-dt
3756  enddo
3757  ttab(it,kp)=t1
3758  qstab(it,kp)=qs
3759  enddo
3760  enddo
3761  !
3762  ! lookup table for tlog(emix/aliq)
3763  !
3764  ! set up intial variable for lookup tables
3765  !
3766  astrt=1.e-3_rp
3767  ainc=0.075_rp
3768  !
3769  a1=astrt-ainc
3770  do i=1,200
3771  a1=a1+ainc
3772  alu(i)=log(a1)
3773  enddo
3774  !GdCP is g/cp add for SCALE
3775  gdcp = - grav/cp ! inital set
3776  return
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:58
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:90
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:48
module CONSTANT
Definition: scale_const.F90:14
Here is the caller graph for this function:

Variable Documentation

◆ kf_dtsec

real(dp), save, public scale_atmos_phy_cp_kf::kf_dtsec = 300._RP

Definition at line 63 of file scale_atmos_phy_cp_kf.F90.

Referenced by cp_kf_main(), cp_kf_param(), and kf_wmean().

63  real(DP), public, save :: KF_DTSEC = 300._rp ! kf time scale [sec]