SCALE-RM
Functions/Subroutines | Variables
scale_atmos_phy_bl_mynn_jmapplib Module Reference

module atmosphere / physics / pbl / mynn-jmapplib More...

Functions/Subroutines

subroutine, public atmos_phy_bl_mynn_jmapplib_setup (KA, KS, KE, CZ, dt, PBL_MAX, SHCU_MAX)
 ATMOS_PHY_BL_MYNN_JMAPPLIB_setup Setup. More...
 
subroutine, public atmos_phy_bl_mynn_jmapplib_finalize
 
subroutine, public atmos_phy_bl_mynn_jmapplib_tendency (KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, U, V, POTT, PROG, PRES, EXNER, QDRY, QV, QC, QI, SFC_DENS, SFC_PRES, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_QV, us, RLmo, CZ, FZ, F2H, dt, RHOU_t, RHOV_t, RHOT_t, RHOQV_t, RPROG_t, Nu, Kh, Zi, SFLX_BUOY)
 ATMOS_PHY_BL_MYNN_JMAPPLIB_tendency calculate tendency by the virtical eddy viscosity. More...
 

Variables

integer, parameter, public atmos_phy_bl_mynn_jmapplib_ntracer = 0
 
character(len=h_short), dimension(1), public atmos_phy_bl_mynn_jmapplib_name = (/ '' /)
 
character(len=h_long), dimension(1), public atmos_phy_bl_mynn_jmapplib_desc = (/ '' /)
 
character(len=h_short), dimension(1), public atmos_phy_bl_mynn_jmapplib_units = (/ '' /)
 

Detailed Description

module atmosphere / physics / pbl / mynn-jmapplib

Description
Boundary layer turbulence model Mellor-Yamada Nakanishi-Niino model implemented in the JMA Physics Process Library
Author
Team SCALE
NAMELIST
  • PARAM_ATMOS_PHY_BL_MYNN_JMAPPLIB
    nametypedefault valuecomment
    ATMOS_PHY_BL_MYNN_JMAPPLIB_LEVEL character(len=3) "3"
    ATMOS_PHY_BL_MYNN_JMAPPLIB_SHCU_BUOY_FLAG logical .true.
    ATMOS_PHY_BL_MYNN_JMAPPLIB_PBL_MAX real(RP) 1.E+99_RP > maximum height of the PBL
    ATMOS_PHY_BL_MYNN_JMAPPLIB_SHCU_MAX real(RP) 1.E+99_RP > maximum height of the shallow cumulus
    ATMOS_PHY_BL_MYNN_JMAPPLIB_SGM_MIN_FCT real(RP) 0.09_RP

History Output
namedescriptionunitvariable
L_mix_MYNN minxing length m l
TKE_diss_MYNN TKE dissipation m2/s3 diss

Function/Subroutine Documentation

◆ atmos_phy_bl_mynn_jmapplib_setup()

subroutine, public scale_atmos_phy_bl_mynn_jmapplib::atmos_phy_bl_mynn_jmapplib_setup ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
real(rp), dimension(ka), intent(in)  CZ,
real(dp), intent(in), optional  dt,
real(rp), intent(in), optional  PBL_MAX,
real(rp), intent(in), optional  SHCU_MAX 
)

ATMOS_PHY_BL_MYNN_JMAPPLIB_setup Setup.

Definition at line 98 of file scale_atmos_phy_bl_mynn_jmapplib.F90.

98  use scale_prc, only: &
99  prc_abort
100  use scale_const, only: &
101  pre00 => const_pre00
102 #ifdef JMAPPLIB
103  use pbl_const, only: &
104  pbl_const_ini
105  use pbl_parm, only: &
106  pbl_parm_ini
107  use pbl_grid, only: &
108  pbl_grid_ini
109  use pbl_mym_option, only: &
110  pbl_mym_option_ini
111  use pbl_mym_option_symbol, only: &
112  mymodel3, &
113  mymodel25
114  use pbl_mym_parm, only: &
115  pbl_mym_parm_ini
116  use pbl_mym_const, only: &
117  pbl_mym_const_ini
118 #endif
119  implicit none
120 
121  integer, intent(in) :: KA, KS, KE
122 
123  real(RP), intent(in) :: CZ(KA)
124 
125  real(DP), intent(in), optional :: dt
126  real(RP), intent(in), optional :: PBL_MAX
127  real(RP), intent(in), optional :: SHCU_MAX
128 
129  real(RP) :: e_emit
130 
131  integer :: KE_shcu
132  integer :: nz_pbl, nz_shcu
133 
134  integer :: ierr
135  integer :: k
136  !---------------------------------------------------------------------------
137 
138  log_newline
139  log_info("ATMOS_PHY_BL_MYNN_JMAPPLIB_setup",*) 'Setup'
140  log_info("ATMOS_PHY_BL_MYNN_JMAPPLIB_setup",*) 'Mellor-Yamada Nakanishi-Niino scheme implemented in the JMA Physics Process Library'
141 
142 #ifdef JMAPPLIB
143  if ( present(pbl_max) ) atmos_phy_bl_mynn_jmapplib_pbl_max = pbl_max
144  if ( present(shcu_max) ) atmos_phy_bl_mynn_jmapplib_shcu_max = shcu_max
145 
146  !--- read namelist
147  rewind(io_fid_conf)
148  read(io_fid_conf,nml=param_atmos_phy_bl_mynn_jmapplib,iostat=ierr)
149  if( ierr > 0 ) then !--- fatal error
150  log_error("ATMOS_PHY_BL_MYNN_JMAPPLIB_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_BL_MYNN_JMAPPLIB. Check!'
151  call prc_abort
152  endif
153  log_nml(param_atmos_phy_bl_mynn_jmapplib)
154 
155 
156  ke_pbl = ks+1
157  ke_shcu = ks+1
158  do k = ks+2, ke-1
159  if ( atmos_phy_bl_mynn_jmapplib_pbl_max >= cz(k) ) then
160  ke_pbl = k
161  end if
162  if ( atmos_phy_bl_mynn_jmapplib_shcu_max >= cz(k) ) then
163  ke_shcu = k
164  end if
165  end do
166  nz_pbl = ke_pbl - ks + 1
167  nz_shcu = min(ke_shcu - ks + 1, nz_pbl)
168 
169  e_emit = 1.0_rp
170  call pbl_const_ini(pref_in = pre00, timestep_in = real(dt,rp), e_emit_in = e_emit)
171  call pbl_parm_ini
172  call pbl_grid_ini(nz_pbl, shcu_levels_in = nz_shcu)
173 
174  select case ( atmos_phy_bl_mynn_jmapplib_level )
175  case ( "3" )
176  call pbl_mym_option_ini(levflag_in = mymodel3, l_shcu_buoy_in = atmos_phy_bl_mynn_jmapplib_shcu_buoy_flag)
177 
178  case ( "2.5" )
179  call pbl_mym_option_ini(levflag_in = mymodel25, l_shcu_buoy_in = atmos_phy_bl_mynn_jmapplib_shcu_buoy_flag)
180  case default
181  log_error("ATMOS_PHY_BL_MYNN_JMAPPLIB_setup",*) 'only level 2.5 and 3 are supported at this moment'
182  call prc_abort
183  end select
184 
185  call pbl_mym_parm_ini(my_sgm_min_fct_in = atmos_phy_bl_mynn_jmapplib_sgm_min_fct)
186  call pbl_mym_const_ini
187 
188 #else
189 
190  log_error("ATMOS_PHY_BL_MYNN_JMAPPLIB_setup",*) 'To use "MYNN-JMAPPLIB", compile SCALE with "SCALE_ENABLE_JMAPPLIB=T" option.'
191  call prc_abort
192 
193 #endif
194 
195  return

References scale_const::const_pre00, scale_io::io_fid_conf, scale_tracer::k, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, scale_prc::prc_abort(), and scale_precision::rp.

Referenced by mod_atmos_phy_bl_driver::atmos_phy_bl_driver_setup().

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

◆ atmos_phy_bl_mynn_jmapplib_finalize()

subroutine, public scale_atmos_phy_bl_mynn_jmapplib::atmos_phy_bl_mynn_jmapplib_finalize

Definition at line 202 of file scale_atmos_phy_bl_mynn_jmapplib.F90.

202  implicit none
203 
204  return

Referenced by mod_atmos_phy_bl_driver::atmos_phy_bl_driver_finalize().

Here is the caller graph for this function:

◆ atmos_phy_bl_mynn_jmapplib_tendency()

subroutine, public scale_atmos_phy_bl_mynn_jmapplib::atmos_phy_bl_mynn_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(out)  DENS,
real(rp), dimension (ka,ia,ja), intent(out)  U,
real(rp), dimension (ka,ia,ja), intent(out)  V,
real(rp), dimension (ka,ia,ja), intent(in)  POTT,
real(rp), dimension (ka,ia,ja,atmos_phy_bl_mynn_jmapplib_ntracer), intent(in)  PROG,
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(out)  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)  SFC_DENS,
real(rp), dimension( ia,ja), intent(in)  SFC_PRES,
real(rp), dimension ( ia,ja), intent(in)  SFLX_MU,
real(rp), dimension ( ia,ja), intent(in)  SFLX_MV,
real(rp), dimension ( ia,ja), intent(in)  SFLX_SH,
real(rp), dimension ( ia,ja), intent(in)  SFLX_QV,
real(rp), dimension ( ia,ja), intent(in)  us,
real(rp), dimension ( ia,ja), intent(in)  RLmo,
real(rp), dimension( ka,ia,ja), intent(in)  CZ,
real(rp), dimension(0:ka,ia,ja), intent(in)  FZ,
real(rp), dimension(ka,2,ia,ja), intent(in)  F2H,
real(dp), intent(in)  dt,
real(rp), dimension (ka,ia,ja), intent(out)  RHOU_t,
real(rp), dimension (ka,ia,ja), intent(out)  RHOV_t,
real(rp), dimension (ka,ia,ja), intent(out)  RHOT_t,
real(rp), dimension(ka,ia,ja), intent(out)  RHOQV_t,
real(rp), dimension(ka,ia,ja,atmos_phy_bl_mynn_jmapplib_ntracer), intent(out)  RPROG_t,
real(rp), dimension (ka,ia,ja), intent(out)  Nu,
real(rp), dimension (ka,ia,ja), intent(out)  Kh,
real(rp), dimension (ia,ja), intent(out)  Zi,
real(rp), dimension (ia,ja), intent(out)  SFLX_BUOY 
)

ATMOS_PHY_BL_MYNN_JMAPPLIB_tendency calculate tendency by the virtical eddy viscosity.

Definition at line 224 of file scale_atmos_phy_bl_mynn_jmapplib.F90.

224  use scale_const, only: &
225  cpdry => const_cpdry
226  use scale_atmos_hydrometeor, only: &
227  cp_vapor
228  use scale_file_history, only: &
229  file_history_in
230 #ifdef JMAPPLIB
231  use pbl_mym_main, only: &
232  pbl_mym_main_level3, &
233  pbl_mym_main_level25
234  use pbl_coupler, only: &
235  pbl_coupler_flx_force_tend_run
236  use pbl_diag, only: &
237  pbl_diag_pbl_height
238 #endif
239  implicit none
240 
241  integer, intent(in) :: KA, KS, KE
242  integer, intent(in) :: IA, IS, IE
243  integer, intent(in) :: JA, JS, JE
244 
245  real(RP), intent(in) :: DENS (KA,IA,JA)
246  real(RP), intent(in) :: U (KA,IA,JA)
247  real(RP), intent(in) :: V (KA,IA,JA)
248  real(RP), intent(in) :: POTT (KA,IA,JA)
249  real(RP), intent(in) :: PROG (KA,IA,JA,ATMOS_PHY_BL_MYNN_JMAPPLIB_ntracer)
250  real(RP), intent(in) :: PRES (KA,IA,JA)
251  real(RP), intent(in) :: EXNER (KA,IA,JA)
252  real(RP), intent(in) :: QDRY (KA,IA,JA)
253  real(RP), intent(in) :: QV (KA,IA,JA)
254  real(RP), intent(in) :: QC (KA,IA,JA)
255  real(RP), intent(in) :: QI (KA,IA,JA)
256  real(RP), intent(in) :: SFC_DENS( IA,JA)
257  real(RP), intent(in) :: SFC_PRES( IA,JA)
258  real(RP), intent(in) :: SFLX_MU ( IA,JA)
259  real(RP), intent(in) :: SFLX_MV ( IA,JA)
260  real(RP), intent(in) :: SFLX_SH ( IA,JA)
261  real(RP), intent(in) :: SFLX_QV ( IA,JA)
262  real(RP), intent(in) :: us ( IA,JA)
263  real(RP), intent(in) :: RLmo ( IA,JA)
264 
265  real(RP), intent(in) :: CZ( KA,IA,JA)
266  real(RP), intent(in) :: FZ(0:KA,IA,JA)
267  real(RP), intent(in) :: F2H(KA,2,IA,JA)
268  real(DP), intent(in) :: dt
269 
270  real(RP), intent(out) :: RHOU_t (KA,IA,JA)
271  real(RP), intent(out) :: RHOV_t (KA,IA,JA)
272  real(RP), intent(out) :: RHOT_t (KA,IA,JA)
273  real(RP), intent(out) :: RHOQV_t(KA,IA,JA)
274  real(RP), intent(out) :: RPROG_t(KA,IA,JA,ATMOS_PHY_BL_MYNN_JMAPPLIB_ntracer)
275  real(RP), intent(out) :: Nu (KA,IA,JA)
276  real(RP), intent(out) :: Kh (KA,IA,JA)
277  real(RP), intent(out) :: Zi (IA,JA)
278  real(RP), intent(out) :: SFLX_BUOY (IA,JA)
279 
280 #ifdef JMAPPLIB
281  real(RP) :: l(KA,IA,JA)
282 
283  real(RP) :: qke(KS:KE_PBL)
284  real(RP) :: qv_lc(KS:KE_PBL)
285  real(RP) :: qc_lc(KS:KE_PBL)
286  real(RP) :: qi_lc(KS:KE_PBL)
287  real(RP) :: dens_lc(KS:KE_PBL)
288  real(RP) :: dfm(KS:KE_PBL)
289  real(RP) :: dfh(KS:KE_PBL)
290  real(RP) :: taux_ex(KS:KE_PBL)
291  real(RP) :: tauy_ex(KS:KE_PBL)
292  real(RP) :: ftl_ex(KS:KE_PBL)
293  real(RP) :: fqw_ex(KS:KE_PBL)
294  real(RP) :: tend_qke(KS:KE_PBL)
295  real(RP) :: tend_tsq(KS:KE_PBL)
296  real(RP) :: tend_qsq(KS:KE_PBL)
297  real(RP) :: tend_cov(KS:KE_PBL)
298  real(RP) :: tend_u(KS:KE_PBL)
299  real(RP) :: tend_v(KS:KE_PBL)
300  real(RP) :: tend_pt(KS:KE_PBL)
301  real(RP) :: tend_qv(KS:KE_PBL)
302  real(RP) :: z_f(KS:KE_PBL)
303  real(RP) :: dz_f(KS:KE_PBL)
304  real(RP) :: rdz_f(KS:KE_PBL)
305  real(RP) :: rdz_h(KS:KE_PBL)
306  real(RP) :: h2f_m(KS:KE_PBL)
307  real(RP) :: h2f_p(KS:KE_PBL)
308  real(RP) :: fb_surf
309  real(RP) :: rho_ov_rhoa
310  real(RP) :: SFLX_U
311  real(RP) :: SFLX_V
312  real(RP) :: SFLX_PT
313  real(RP) :: SFLX_Q
314  real(RP) :: CPtot
315 
316  real(RP) :: diss(KA,IA,JA)
317  real(RP) :: dummy(KA,IA,JA)
318 
319  integer :: k, i, j
320  !---------------------------------------------------------------------------
321 
322  log_progress(*) "atmosphere / physics / pbl / MYNN-JMAPPLIB"
323 
324  h2f_m(:) = 0.5_rp
325  h2f_p(:) = 0.5_rp
326 
327  !$omp parallel do &
328  !$omp private(qke,qv_lc,qc_lc,qi_lc,dens_lc,taux_ex,tauy_ex,ftl_ex,fqw_ex,&
329  !$omp tend_qke,tend_tsq,tend_qsq,tend_cov,tend_u,tend_v,tend_pt,tend_qv, &
330  !$omp z_f,dz_f,rdz_f,rdz_h, &
331  !$omp rho_ov_rhoa,SFLX_U,SFLX_V,SFLX_PT,SFLX_Q,CPtot)
332  do j = js, je
333  do i = is, ie
334 
335  do k = ks, ke_pbl
336  qke(k) = prog(k,i,j,i_tke) * 2.0_rp
337  rho_ov_rhoa = 1.0_rp / ( qdry(k,i,j) + qv(k,i,j) )
338  qv_lc(k) = qv(k,i,j) * rho_ov_rhoa
339  qc_lc(k) = qc(k,i,j) * rho_ov_rhoa
340  qi_lc(k) = qi(k,i,j) * rho_ov_rhoa
341  z_f(k) = cz(k,i,j) - fz(ks-1,i,j)
342  dz_f(k) = fz(k,i,j) - fz(k-1,i,j)
343  rdz_f(k) = 1.0_rp / dz_f(k)
344  rdz_h(k) = 1.0_rp / ( cz(k+1,i,j) - cz(k,i,j) )
345  end do
346 
347  cptot = cpdry + sflx_qv(i,j) * ( cp_vapor - cpdry )
348 
349  sflx_u = sflx_mu(i,j) / sfc_dens(i,j)
350  sflx_v = sflx_mv(i,j) / sfc_dens(i,j)
351  sflx_pt = sflx_sh(i,j) / ( cptot * exner(ks,i,j) * sfc_dens(i,j) )
352  sflx_q = sflx_qv(i,j) / sfc_dens(i,j)
353 
354  select case ( atmos_phy_bl_mynn_jmapplib_level )
355  case ( "3" )
356  call pbl_mym_main_level3( &
357  sflx_u, sflx_v, sflx_pt, sflx_q, & ! (in)
358  us(i,j), rlmo(i,j), sfc_pres(i,j), & ! (in)
359  u(ks:ke_pbl,i,j), v(ks:ke_pbl,i,j), pott(ks:ke_pbl,i,j), & ! (in)
360  qv_lc, qc_lc, qi_lc, pres(ks:ke_pbl,i,j), exner(ks:ke_pbl,i,j), & ! (in)
361  qke(:), prog(ks:ke_pbl,i,j,i_tsq), prog(ks:ke_pbl,i,j,i_qsq), prog(ks:ke_pbl,i,j,i_cov), & ! (in)
362  z_f, dz_f, rdz_f, rdz_h, f2h(ks:ke_pbl,2,i,j), f2h(ks:ke_pbl,1,i,j), h2f_m, h2f_p, & ! (in)
363  sflx_buoy(i,j), & ! (out)
364  dfm(:), dfh(:), tend_qke(:), tend_tsq(:), tend_qsq(:), tend_cov(:), & ! (out)
365  taux_ex, tauy_ex, ftl_ex, fqw_ex, & ! (out)
366  l_e = l(ks:ke_pbl,i,j), tke_dis = diss(ks:ke_pbl,i,j), tsq_dis=dummy(ks:ke_pbl,i,j) ) ! (out)
367  case ( "2.5" )
368  call pbl_mym_main_level25( &
369  sflx_u, sflx_v, sflx_pt, sflx_q, & ! (in)
370  us(i,j), rlmo(i,j), sfc_pres(i,j), & ! (in)
371  u(ks:ke_pbl,i,j), v(ks:ke_pbl,i,j), pott(ks:ke_pbl,i,j), & ! (in)
372  qv_lc, qc_lc, qi_lc, pres(ks:ke_pbl,i,j), exner(ks:ke_pbl,i,j), & ! (in)
373  qke(:), prog(ks:ke_pbl,i,j,i_tsq), prog(ks:ke_pbl,i,j,i_qsq), prog(ks:ke_pbl,i,j,i_cov), & ! (in)
374  z_f, dz_f, rdz_f, rdz_h, f2h(ks:ke_pbl,2,i,j), f2h(ks:ke_pbl,1,i,j), h2f_m, h2f_p, & ! (in)
375  sflx_buoy(i,j), & ! (out)
376  dfm(:), dfh(:), tend_qke(:), tend_tsq(:), tend_qsq(:), tend_cov(:), & ! (out)
377  taux_ex, tauy_ex, ftl_ex, fqw_ex, & ! (out)
378  l_e = l(ks:ke_pbl,i,j), tke_dis = diss(ks:ke_pbl,i,j), tsq_dis=dummy(ks:ke_pbl,i,j) ) ! (out)
379  do k = ks, ke_pbl
380  tend_tsq(k) = ( tend_tsq(k) - prog(k,i,j,i_tsq) ) / dt
381  tend_qsq(k) = ( tend_qsq(k) - prog(k,i,j,i_qsq) ) / dt
382  tend_cov(k) = ( tend_cov(k) - prog(k,i,j,i_cov) ) / dt
383  end do
384  end select
385 
386  do k = ks, ke_pbl
387  dens_lc(k) = dens(k,i,j) * ( qdry(k,i,j) + qv(k,i,j) )
388  end do
389 
390  call pbl_coupler_flx_force_tend_run( &
391  sfc_dens(i,j), sflx_u, sflx_v, sflx_pt, sflx_q, & ! (in)
392  dfm(:), dfh(:), dens_lc(:), & ! (in)
393  taux_ex, tauy_ex, ftl_ex, fqw_ex, & ! (in)
394  f2h(ks:ke_pbl,2,i,j), f2h(ks:ke_pbl,1,i,j), rdz_f, rdz_h, & ! (in)
395  tend_u, tend_v, tend_pt, tend_qv ) ! (out)
396 
397  do k = ks, ke_pbl
398  rhou_t(k,i,j) = tend_u(k) * dens(k,i,j)
399  rhov_t(k,i,j) = tend_v(k) * dens(k,i,j)
400  rhot_t(k,i,j) = tend_pt(k) * dens(k,i,j)
401  rhoqv_t(k,i,j) = tend_qv(k) * dens(k,i,j) * ( qdry(k,i,j) + qv(k,i,j) )
402  rprog_t(k,i,j,i_tke) = tend_qke(k) * dens(k,i,j) * 0.5_rp
403  rprog_t(k,i,j,i_tsq) = tend_tsq(k) * dens(k,i,j)
404  rprog_t(k,i,j,i_qsq) = tend_qsq(k) * dens(k,i,j)
405  rprog_t(k,i,j,i_cov) = tend_cov(k) * dens(k,i,j)
406  end do
407  rhou_t(ks,i,j) = rhou_t(ks,i,j) - sflx_mu(i,j) * rdz_f(ks)
408  rhov_t(ks,i,j) = rhov_t(ks,i,j) - sflx_mv(i,j) * rdz_f(ks)
409  rhot_t(ks,i,j) = rhot_t(ks,i,j) - sflx_pt * sfc_dens(i,j) * rdz_f(ks)
410  rhoqv_t(ks,i,j) = rhoqv_t(ks,i,j) - sflx_qv(i,j) * rdz_f(ks)
411  do k = ke_pbl+1, ke
412  rhou_t(k,i,j) = 0.0_rp
413  rhov_t(k,i,j) = 0.0_rp
414  rhot_t(k,i,j) = 0.0_rp
415  rhoqv_t(k,i,j) = 0.0_rp
416  rprog_t(k,i,j,i_tke) = 0.0_rp
417  rprog_t(k,i,j,i_tsq) = 0.0_rp
418  rprog_t(k,i,j,i_qsq) = 0.0_rp
419  rprog_t(k,i,j,i_cov) = 0.0_rp
420  end do
421 
422  nu(ks-1,i,j) = 0.0_rp
423  kh(ks-1,i,j) = 0.0_rp
424  do k = ks, ke_pbl-1
425  nu(k,i,j) = dfm(k+1) * f2h(k,1,i,j) + dfm(k) * f2h(k,2,i,j)
426  kh(k,i,j) = dfh(k+1) * f2h(k,1,i,j) + dfh(k) * f2h(k,2,i,j)
427  end do
428  do k = ke_pbl, ke
429  nu(k,i,j) = 0.0_rp
430  kh(k,i,j) = 0.0_rp
431  end do
432 
433  call pbl_diag_pbl_height( &
434  sflx_pt, us(i,j), rlmo(i,j), & ! (in)
435  u(ks:ke_pbl,i,j), v(ks:ke_pbl,i,j), pott(ks:ke_pbl,i,j), qv_lc, z_f, & ! (in)
436  zi(i,j) ) ! (out)
437 
438  end do
439  end do
440 
441 
442  l(ke_pbl+1:ke,:,:) = undef
443  call file_history_in(l(:,:,:), 'L_mix_MYNN', 'minxing length', 'm', fill_halo=.true.)
444 
445  diss(ks:ke_pbl,:,:) = - diss(ks:ke_pbl,:,:)
446  diss(ke_pbl+1:ke,:,:) = undef
447  call file_history_in(diss(:,:,:), 'TKE_diss_MYNN', 'TKE dissipation', 'm2/s3', fill_halo=.true.)
448 #endif
449 
450  return

References scale_const::const_cpdry, scale_atmos_hydrometeor::cp_vapor, 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::ke, and scale_atmos_grid_cartesc_index::ks.

Referenced by mod_atmos_phy_bl_driver::atmos_phy_bl_driver_calc_tendency().

Here is the caller graph for this function:

Variable Documentation

◆ atmos_phy_bl_mynn_jmapplib_ntracer

integer, parameter, public scale_atmos_phy_bl_mynn_jmapplib::atmos_phy_bl_mynn_jmapplib_ntracer = 0

Definition at line 52 of file scale_atmos_phy_bl_mynn_jmapplib.F90.

52  integer, public, parameter :: ATMOS_PHY_BL_MYNN_JMAPPLIB_NTRACER = 0

Referenced by mod_atmos_phy_bl_driver::atmos_phy_bl_driver_tracer_setup().

◆ atmos_phy_bl_mynn_jmapplib_name

character(len=h_short), dimension(1), public scale_atmos_phy_bl_mynn_jmapplib::atmos_phy_bl_mynn_jmapplib_name = (/ '' /)

Definition at line 53 of file scale_atmos_phy_bl_mynn_jmapplib.F90.

53  character(len=H_SHORT), public :: ATMOS_PHY_BL_MYNN_JMAPPLIB_NAME(1) = (/ '' /)

Referenced by mod_atmos_phy_bl_driver::atmos_phy_bl_driver_tracer_setup().

◆ atmos_phy_bl_mynn_jmapplib_desc

character(len=h_long), dimension(1), public scale_atmos_phy_bl_mynn_jmapplib::atmos_phy_bl_mynn_jmapplib_desc = (/ '' /)

Definition at line 54 of file scale_atmos_phy_bl_mynn_jmapplib.F90.

54  character(len=H_LONG), public :: ATMOS_PHY_BL_MYNN_JMAPPLIB_DESC(1) = (/ '' /)

Referenced by mod_atmos_phy_bl_driver::atmos_phy_bl_driver_tracer_setup().

◆ atmos_phy_bl_mynn_jmapplib_units

character(len=h_short), dimension(1), public scale_atmos_phy_bl_mynn_jmapplib::atmos_phy_bl_mynn_jmapplib_units = (/ '' /)

Definition at line 55 of file scale_atmos_phy_bl_mynn_jmapplib.F90.

55  character(len=H_SHORT), public :: ATMOS_PHY_BL_MYNN_JMAPPLIB_UNITS(1) = (/ '' /)

Referenced by mod_atmos_phy_bl_driver::atmos_phy_bl_driver_tracer_setup().

scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
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_const
module CONSTANT
Definition: scale_const.F90:11
scale_const::const_cpdry
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:60
mod_atmos_phy_tb_vars::i_tke
integer, public i_tke
Definition: mod_atmos_phy_tb_vars.F90:62
scale_const::const_pre00
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:97
scale_atmos_hydrometeor::cp_vapor
real(rp), public cp_vapor
CP for vapor [J/kg/K].
Definition: scale_atmos_hydrometeor.F90:150