SCALE-RM
Functions/Subroutines
scale_atmos_phy_cp_kf_jmapplib Module Reference

module atmosphere / physics / cumulus / kf-jmapplib More...

Functions/Subroutines

subroutine, public atmos_phy_cp_kf_jmapplib_setup (KA, KS, KE, IA, JA, dx, dt)
 ATMOS_PHY_CP_KF_JMAPPLIB_setup Setup. More...
 
subroutine, public atmos_phy_cp_kf_jmapplib_finalize
 
subroutine, public atmos_phy_cp_kf_jmapplib_tendency (KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, U, V, W, TEMP, POTT, PRES, EXNER, QDRY, QV, QC, QI, us, PBLH, SFLX_BUOY, CZ, FZ, DENS_t, RHOT_t, RHOQV_t, RHOQ_t, w0avg, nca, SFLX_rain, SFLX_snow, SFLX_prec, cloudtop, cloudbase)
 ATMOS_PHY_CP_KF_JMAPPLIB_tendency calculate tendency by the virtical eddy viscosity. More...
 

Detailed Description

module atmosphere / physics / cumulus / kf-jmapplib

Description
Cumulus parameterization Kain-Fritsch scheme implemented in JMA Physics Process library
Author
Team SCALE
NAMELIST
  • PARAM_ATMOS_PHY_CP_KF_JMAPPLIB
    nametypedefault valuecomment
    ATMOS_PHY_CP_KF_JMAPPLIB_TYPE character(len=6) "KF1701" > KF1701 for MSM; KF for LFM
    ATMOS_PHY_CP_KF_JMAPPLIB_RADIUS_TYPE integer 3 > 3 for KF1701
    ATMOS_PHY_CP_KF_JMAPPLIB_DTLCL_FACT real(RP) 0.5_RP > ignored for KF1701
    ATMOS_PHY_CP_KF_JMAPPLIB_DPMIN real(RP) 5000.0_RP > 5000 for KF1701; 1000 for KF
    ATMOS_PHY_CP_KF_JMAPPLIB_DLIFETIME real(RP) 600.0_RP > 600 for KF1701; 3600 for KF
    ATMOS_PHY_CP_KF_JMAPPLIB_SLIFETIME real(RP) 600.0_RP > 600 for KF1701; 3600 for KF
    ATMOS_PHY_CP_KF_JMAPPLIB_DETLQ_QR_RATE real(RP) 0.0_RP > 0.0 for KF1701; 0.0 for KF
    ATMOS_PHY_CP_KF_JMAPPLIB_DETIC_QS_RATE real(RP) 1.0_RP > 1.0 for KF1701; 0.0 for KF

History Output
No history output

Function/Subroutine Documentation

◆ atmos_phy_cp_kf_jmapplib_setup()

subroutine, public scale_atmos_phy_cp_kf_jmapplib::atmos_phy_cp_kf_jmapplib_setup ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  JA,
real(rp), intent(in)  dx,
real(dp), intent(in)  dt 
)

ATMOS_PHY_CP_KF_JMAPPLIB_setup Setup.

Definition at line 62 of file scale_atmos_phy_cp_kf_jmapplib.F90.

62  use scale_prc, only: &
63  prc_abort
64  use scale_const, only: &
65  pre00 => const_pre00
66 #ifdef JMAPPLIB
67  use conv_grid, only: &
68  conv_grid_ini
69  use conv_kf_const, only: &
70  conv_kf_const_ini
71  use conv_kf_parm, only: &
72  conv_kf_parm_ini, &
73  conv_kf_parm_get_stepcu
74  use conv_kf_parm1701, only: &
75  conv_kf_parm_ini1701, &
76  conv_kf_parm_get_stepcu1701
77  use conv_kf_lut, only: &
78  conv_kf_lut_ini
79 #endif
80  implicit none
81 
82  integer, intent(in) :: KA, KS, KE
83  integer, intent(in) :: IA, JA
84 
85  real(RP), intent(in) :: dx
86  real(DP), intent(in) :: dt
87 
88 #ifdef JMAPPLIB
89 
90  integer :: ATMOS_PHY_CP_KF_JMAPPLIB_radius_type = 3
94  real(RP) :: ATMOS_PHY_CP_KF_JMAPPLIB_dtlcl_fact = 0.5_rp
95  real(RP) :: ATMOS_PHY_CP_KF_JMAPPLIB_dpmin = 5000.0_rp
96  real(RP) :: ATMOS_PHY_CP_KF_JMAPPLIB_dlifetime = 600.0_rp
97  real(RP) :: ATMOS_PHY_CP_KF_JMAPPLIB_slifetime = 600.0_rp
98  real(RP) :: ATMOS_PHY_CP_KF_JMAPPLIB_detlq_qr_rate = 0.0_rp
99  real(RP) :: ATMOS_PHY_CP_KF_JMAPPLIB_detic_qs_rate = 1.0_rp
100 
101  namelist / param_atmos_phy_cp_kf_jmapplib / &
102  atmos_phy_cp_kf_jmapplib_type, &
103  atmos_phy_cp_kf_jmapplib_radius_type, &
104  atmos_phy_cp_kf_jmapplib_dtlcl_fact, &
105  atmos_phy_cp_kf_jmapplib_dpmin, &
106  atmos_phy_cp_kf_jmapplib_dlifetime, &
107  atmos_phy_cp_kf_jmapplib_slifetime, &
108  atmos_phy_cp_kf_jmapplib_detlq_qr_rate, &
109  atmos_phy_cp_kf_jmapplib_detic_qs_rate
110 
111  integer :: KMAX
112 
113  integer :: ierr
114  integer :: k
115  !---------------------------------------------------------------------------
116 
117  log_newline
118  log_info("ATMOS_PHY_CP_KF_JMAPPLIB_setup",*) 'Setup'
119  log_info("ATMOS_PHY_CP_KF_JMAPPLIB_setup",*) 'KF scheme implemented in the JMA Physics Process Library'
120 
121  !--- read namelist
122  rewind(io_fid_conf)
123  read(io_fid_conf,nml=param_atmos_phy_cp_kf_jmapplib,iostat=ierr)
124  if( ierr > 0 ) then !--- fatal error
125  log_error("ATMOS_PHY_CP_KF_JMAPPLIB_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_CP_KF_JMAPPLIB. Check!'
126  call prc_abort
127  endif
128  log_nml(param_atmos_phy_cp_kf_jmapplib)
129 
130 
131  kmax = ke - ks + 1
132  call conv_grid_ini( kmax )
133  call conv_kf_const_ini( dx, pre00, dt )
134 
135  select case ( atmos_phy_cp_kf_jmapplib_type )
136  case ( "KF" )
137 
138  call conv_kf_parm_ini( &
139  dt, &
140  kfrad_var_in = atmos_phy_cp_kf_jmapplib_radius_type, &
141  cudt_in = -1.0_rp, &
142  dtlcl_fct_in = atmos_phy_cp_kf_jmapplib_dtlcl_fact, &
143  dpmin_in = atmos_phy_cp_kf_jmapplib_dpmin, &
144  cu_lifetime_min_in = atmos_phy_cp_kf_jmapplib_dlifetime, &
145  shallow_lifetime_in = atmos_phy_cp_kf_jmapplib_slifetime, &
146  detlq_qr_r_in = atmos_phy_cp_kf_jmapplib_detlq_qr_rate, &
147  detic_qs_r_in = atmos_phy_cp_kf_jmapplib_detic_qs_rate )
148 
149 ! call conv_kf_parm_get_stepcu(stepcu) ! (out)
150 
151  case ( "KF1701" )
152 
153  call conv_kf_parm_ini1701( &
154  dt, &
155  kfrad_var_in = atmos_phy_cp_kf_jmapplib_radius_type, &
156  kfrad_in = 1000.0_rp, &
157  cudt_in = -1.0_rp, &
158  dpmin_in = atmos_phy_cp_kf_jmapplib_dpmin, &
159  cu_lifetime_min_in = atmos_phy_cp_kf_jmapplib_dlifetime, &
160  shallow_lifetime_in = atmos_phy_cp_kf_jmapplib_slifetime, &
161  detlq_qr_r_in = atmos_phy_cp_kf_jmapplib_detlq_qr_rate, &
162  detic_qs_r_in = atmos_phy_cp_kf_jmapplib_detic_qs_rate, &
163  tkemax_in = 1.0_rp )
164 
165 ! call conv_kf_parm_get_stepcu1701(stepcu) ! (out)
166 
167  case default
168  log_error("ATMOS_PHY_CP_KF_JMAPPLIB_setup",*) 'ATMOS_PHY_CP_KF_TYPE must be either "KF" or "KF1701"'
169  call prc_abort
170  end select
171 
172  call conv_kf_lut_ini
173 
174  allocate( ishall_counter(ia,ja), ideep_counter(ia,ja) )
175  ishall_counter(:,:) = 0.0_rp
176  ideep_counter(:,:) = 0.0_rp
177 
178 #else
179 
180  log_error("ATMOS_PHY_CP_KF_JMAPPLIB_setup",*) 'To use "KF-JMAPPLIB", compile SCALE with "SCALE_ENABLE_JMAPPLIB=T" option.'
181  call prc_abort
182 
183 #endif
184 
185  return

References scale_const::const_pre00, scale_io::io_fid_conf, 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_jmapplib_finalize()

subroutine, public scale_atmos_phy_cp_kf_jmapplib::atmos_phy_cp_kf_jmapplib_finalize

Definition at line 192 of file scale_atmos_phy_cp_kf_jmapplib.F90.

192  implicit none
193 
194 #ifdef JMAPPLIB
195  deallocate( ishall_counter, ideep_counter )
196 #endif
197 
198  return

◆ atmos_phy_cp_kf_jmapplib_tendency()

subroutine, public scale_atmos_phy_cp_kf_jmapplib::atmos_phy_cp_kf_jmapplib_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)  W,
real(rp), dimension (ka,ia,ja), intent(in)  TEMP,
real(rp), dimension (ka,ia,ja), intent(in)  POTT,
real(rp), dimension (ka,ia,ja), intent(in)  PRES,
real(rp), dimension(ka,ia,ja), intent(in)  EXNER,
real(rp), dimension (ka,ia,ja), intent(in)  QDRY,
real(rp), dimension (ka,ia,ja), intent(in)  QV,
real(rp), dimension (ka,ia,ja), intent(in)  QC,
real(rp), dimension (ka,ia,ja), intent(in)  QI,
real(rp), dimension (ia,ja), intent(in)  us,
real(rp), dimension (ia,ja), intent(in)  PBLH,
real(rp), dimension(ia,ja), intent(in)  SFLX_BUOY,
real(rp), dimension (ka,ia,ja), intent(in)  CZ,
real(rp), dimension (0:ka,ia,ja), intent(in)  FZ,
real(rp), dimension (ka,ia,ja), intent(out)  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 (ka,ia,ja), intent(inout)  w0avg,
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_prec,
real(rp), dimension (ia,ja), intent(inout)  cloudtop,
real(rp), dimension(ia,ja), intent(inout)  cloudbase 
)

ATMOS_PHY_CP_KF_JMAPPLIB_tendency calculate tendency by the virtical eddy viscosity.

Definition at line 217 of file scale_atmos_phy_cp_kf_jmapplib.F90.

217  use scale_atmos_hydrometeor, only: &
218  n_hyd, &
219  i_hc, &
220  i_hr, &
221  i_hi, &
222  i_hs, &
223  i_hg
224 #ifdef JMAPPLIB
225  use conv_kf_main, only: &
226  conv_kf_run
227  use conv_kf_main1701, only: &
228  conv_kf_run1701
229  use pbl_diag, only: &
230  pbl_diag_buoyancy_excess
231 #endif
232  implicit none
233 
234  integer, intent(in) :: KA, KS, KE
235  integer, intent(in) :: IA, IS, IE
236  integer, intent(in) :: JA, JS, JE
237 
238  real(RP), intent(in) :: DENS (KA,IA,JA)
239  real(RP), intent(in) :: U (KA,IA,JA)
240  real(RP), intent(in) :: V (KA,IA,JA)
241  real(RP), intent(in) :: W (KA,IA,JA)
242  real(RP), intent(in) :: TEMP (KA,IA,JA)
243  real(RP), intent(in) :: POTT (KA,IA,JA)
244  real(RP), intent(in) :: PRES (KA,IA,JA)
245  real(RP), intent(in) :: EXNER(KA,IA,JA)
246  real(RP), intent(in) :: QDRY (KA,IA,JA)
247  real(RP), intent(in) :: QV (KA,IA,JA)
248  real(RP), intent(in) :: QC (KA,IA,JA)
249  real(RP), intent(in) :: QI (KA,IA,JA)
250  real(RP), intent(in) :: us (IA,JA)
251  real(RP), intent(in) :: PBLH (IA,JA)
252  real(RP), intent(in) :: SFLX_BUOY(IA,JA)
253  real(RP), intent(in) :: CZ (KA,IA,JA)
254  real(RP), intent(in) :: FZ (0:KA,IA,JA)
255 
256  real(RP), intent(out) :: DENS_t (KA,IA,JA)
257  real(RP), intent(inout) :: RHOT_t (KA,IA,JA)
258  real(RP), intent(inout) :: RHOQV_t(KA,IA,JA)
259  real(RP), intent(inout) :: RHOQ_t (KA,IA,JA,N_HYD)
260  real(RP), intent(inout) :: w0avg (KA,IA,JA)
261  real(RP), intent(inout) :: nca (IA,JA)
262  real(RP), intent(inout) :: SFLX_rain(IA,JA)
263  real(RP), intent(inout) :: SFLX_snow(IA,JA)
264  real(RP), intent(inout) :: SFLX_prec(IA,JA)
265  real(RP), intent(inout) :: cloudtop (IA,JA)
266  real(RP), intent(inout) :: cloudbase(IA,JA)
267 
268 #ifdef JMAPPLIB
269 
270  integer, parameter :: iconv_run = 1
271 
272  real(RP) :: tend_pt(KS:KE)
273  real(RP) :: tend_qv(KS:KE)
274  real(RP) :: tend_qc(KS:KE)
275  real(RP) :: tend_qr(KS:KE)
276  real(RP) :: tend_qi(KS:KE)
277  real(RP) :: tend_qs(KS:KE)
278  real(RP) :: tend_qg(KS:KE)
279  real(RP) :: dens_lc(KS:KE)
280  real(RP) :: qv_lc (KS:KE)
281  real(RP) :: qc_lc
282  real(RP) :: qi_lc
283  real(RP) :: z_f (KS:KE)
284  real(RP) :: dz_f (KS:KE)
285 
286  real(RP) :: lcl_temp(IA,JA)
287  real(RP) :: abe_out
288  real(RP) :: tv_ex
289  real(RP) :: ptv_ex
290 
291  integer :: k, i, j
292  !---------------------------------------------------------------------------
293 
294  log_progress(*) "atmosphere / physics / cumulus / KF-JMAPPLIB"
295 
296  !$omp parallel do &
297  !$omp private(tend_pt,tend_qv,tend_qc,tend_qr,tend_qi,tend_qs,tend_qg, &
298  !$omp dens_lc,qv_lc,qc_lc,qi_lc,z_f,dz_f, &
299  !$omp ptv_ex,tv_ex,abe_out)
300  do j = js, je
301  do i = is, ie
302 
303  do k = ks, ke
304  dens_lc(k) = dens(k,i,j) * qdry(k,i,j)
305  qv_lc(k) = qv(k,i,j) / qdry(k,i,j)
306  tend_pt(k) = rhot_t(k,i,j) / dens(k,i,j)
307  tend_qv(k) = rhoqv_t(k,i,j) / dens_lc(k)
308  tend_qc(k) = rhoq_t(k,i,j,i_hc) / dens_lc(k)
309  tend_qr(k) = rhoq_t(k,i,j,i_hr) / dens_lc(k)
310  tend_qi(k) = rhoq_t(k,i,j,i_hi) / dens_lc(k)
311  tend_qs(k) = rhoq_t(k,i,j,i_hs) / dens_lc(k)
312  tend_qg(k) = rhoq_t(k,i,j,i_hg) / dens_lc(k)
313 
314  z_f(k) = cz(k,i,j) - fz(ks-1,i,j)
315  dz_f(k) = fz(k,i,j) - fz(k-1,i,j)
316  end do
317 
318  select case ( atmos_phy_cp_kf_jmapplib_type )
319  case ( "KF" )
320 
321  call conv_kf_run( &
322  iconv_run, & ! (in)
323  u(ks:ke,i,j), v(ks:ke,i,j), w(ks:ke,i,j), temp(ks:ke,i,j), qv_lc(:), & ! (in)
324  pres(ks:ke,i,j), exner(ks:ke,i,j), dens_lc(:), z_f(:), dz_f(:), & ! (in)
325  nca(i,j), & ! (inout)
326  sflx_prec(i,j), sflx_rain(i,j), sflx_snow(i,j), & ! (inout)
327  cloudtop(i,j), cloudbase(i,j), abe_out, & ! (inout)
328  ishall_counter(i,j), ideep_counter(i,j), & ! (inout)
329  lcl_temp(i,j), & ! (inout)
330  w0avg(ks:ke,i,j), & ! (inout)
331  tend_pt(:), tend_qv(:), tend_qc(:), & ! (inout)
332  tend_qi(:), tend_qr(:), tend_qs(:) ) ! (inout)
333 
334  tend_qg(ks:ke) = 0.0_rp
335 
336  case ( "KF1701" )
337 
338  qc_lc = qc(ks,i,j) / qdry(ks,i,j)
339  qi_lc = qi(ks,i,j) / qdry(ks,i,j)
340  call pbl_diag_buoyancy_excess( &
341  pott(ks,i,j), qv_lc(ks),qc_lc, qi_lc, & ! (in)
342  pblh(i,j), us(i,j), sflx_buoy(i,j), exner(ks,i,j), & ! (in)
343  ptv_ex, tv_ex ) ! (out)
344 
345  call conv_kf_run1701( &
346  iconv_run, & ! (in)
347  u(ks:ke,i,j), v(ks:ke,i,j), w(ks:ke,i,j), temp(ks:ke,i,j), qv_lc(:), & ! (in)
348  pres(ks:ke,i,j), exner(ks:ke,i,j), dens_lc(:), z_f(:), dz_f(:), & ! (in)
349  pblh(i,j), tv_ex, & ! (in)
350  nca(i,j), & ! (inout)
351  sflx_prec(i,j), sflx_rain(i,j), sflx_snow(i,j), & ! (inout)
352  cloudtop(i,j), cloudbase(i,j), abe_out, & ! (inout)
353  ishall_counter(i,j), ideep_counter(i,j), & ! (inout)
354  lcl_temp(i,j), & ! (inout)
355  w0avg(ks:ke,i,j), & ! (inout)
356  tend_pt(:), tend_qv(:), tend_qc(:), & ! (inout)
357  tend_qi(:), tend_qr(:), & ! (inout)
358  tend_qs(:), tend_qg(:) ) ! (inout)
359 
360  end select
361 
362  if ( nca(i,j) < -1 ) nca(i,j) = -1
363 
364  do k = ks, ke
365  rhot_t(k,i,j) = tend_pt(k) * dens(k,i,j)
366  rhoqv_t(k,i,j) = tend_qv(k) * dens_lc(k)
367  rhoq_t(k,i,j,i_hc) = tend_qc(k) * dens_lc(k)
368  rhoq_t(k,i,j,i_hr) = tend_qr(k) * dens_lc(k)
369  rhoq_t(k,i,j,i_hi) = tend_qi(k) * dens_lc(k)
370  rhoq_t(k,i,j,i_hs) = tend_qs(k) * dens_lc(k)
371  rhoq_t(k,i,j,i_hg) = tend_qg(k) * dens_lc(k)
372  dens_t(k,i,j) = rhoqv_t(k,i,j) &
373  + rhoq_t(k,i,j,i_hc) + rhoq_t(k,i,j,i_hr) &
374  + rhoq_t(k,i,j,i_hi) + rhoq_t(k,i,j,i_hs) + rhoq_t(k,i,j,i_hg)
375  end do
376 
377  end do
378  end do
379 
380 #endif
381 
382  return

References scale_atmos_hydrometeor::i_hc, scale_atmos_hydrometeor::i_hg, scale_atmos_hydrometeor::i_hi, scale_atmos_hydrometeor::i_hr, scale_atmos_hydrometeor::i_hs, and scale_atmos_hydrometeor::n_hyd.

Referenced by mod_atmos_phy_cp_driver::atmos_phy_cp_driver_calc_tendency().

Here is the caller graph for this function:
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::i_hs
integer, parameter, public i_hs
snow
Definition: scale_atmos_hydrometeor.F90:100
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
scale_prc
module PROCESS
Definition: scale_prc.F90:11
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_const::const_pre00
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:97
scale_atmos_hydrometeor::i_hg
integer, parameter, public i_hg
graupel
Definition: scale_atmos_hydrometeor.F90:101