SCALE-RM
Functions/Subroutines
mod_atmos_phy_bl_driver Module Reference

module atmosphere / physics / PBL More...

Functions/Subroutines

subroutine, public atmos_phy_bl_driver_tracer_setup
 Config. More...
 
subroutine, public atmos_phy_bl_driver_setup
 Setup. More...
 
subroutine, public atmos_phy_bl_driver_calc_tendency (update_flag)
 calculate tendency More...
 

Detailed Description

module atmosphere / physics / PBL

Description
Planetary boundary layer turbulence
Author
Team SCALE
NAMELIST
  • No namelist group
History Output
namedescriptionunitvariable
Kh_BL eddy diffusion m2/s Kh
Nu_BL eddy viscosity m2/s Nu
QL_BL liquid water content in partial condensation kg/kg QL
RHOT_t_BL RHOT tendency (BL) K.kg/m3/s RHOT_t_BL
RHOU_t_BL MOMX tendency (BL) kg/m2/s2 RHOU_t_BL
RHOV_t_BL MOMY tendency (BL) kg/m2/s2 RHOV_t_BL
Zi_BL depth of the boundary layer m Zi
cldfrac_BL cloud fraction in partial condensation 1 cldfrac
{TRACER_NAME}_t_BL RHO*{TRACER_NAME} tendency (BL);
{TRACER_NAME} depends on the physics schemes, e.g., QV, QC, QR.
kg/m3/s RHOQ_t_BL

Function/Subroutine Documentation

◆ atmos_phy_bl_driver_tracer_setup()

subroutine, public mod_atmos_phy_bl_driver::atmos_phy_bl_driver_tracer_setup

Config.

Definition at line 50 of file mod_atmos_phy_bl_driver.F90.

50  use scale_prc, only: &
51  prc_abort
52  use scale_tracer, only: &
54  use scale_atmos_phy_bl_mynn, only: &
60  use mod_atmos_admin, only: &
63  use mod_atmos_phy_bl_vars, only: &
64  qs, qe
65  implicit none
66  !---------------------------------------------------------------------------
67 
68  log_newline
69  log_info("ATMOS_PHY_BL_driver_tracer_setup",*) 'Setup'
70 
71  if ( atmos_sw_phy_bl ) then
72  select case ( atmos_phy_bl_type )
73  case ( 'MYNN' )
75  call tracer_regist( &
76  qs, &
82  case default
83  log_error("ATMOS_PHY_BL_driver_tracer_setup",*) 'ATMOS_PHY_BL_TYPE is invalid: ', trim(atmos_phy_bl_type)
84  call prc_abort
85  end select
86  end if
87 
88  return

References scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_desc, scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_name, scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_ntracer, scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_tracer_setup(), scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_units, mod_atmos_admin::atmos_phy_bl_type, mod_atmos_admin::atmos_sw_phy_bl, scale_prc::prc_abort(), mod_atmos_phy_bl_vars::qe, mod_atmos_phy_bl_vars::qs, and scale_tracer::tracer_regist().

Referenced by mod_atmos_driver::atmos_driver_tracer_setup().

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

◆ atmos_phy_bl_driver_setup()

subroutine, public mod_atmos_phy_bl_driver::atmos_phy_bl_driver_setup

Setup.

Definition at line 94 of file mod_atmos_phy_bl_driver.F90.

94  use scale_atmos_phy_bl_mynn, only: &
96  use mod_atmos_admin, only: &
99  use mod_atmos_phy_bl_vars, only: &
101  use scale_bulkflux, only: &
103  implicit none
104  !---------------------------------------------------------------------------
105 
106  log_newline
107  log_info("ATMOS_PHY_BL_driver_setup",*) 'Setup'
108 
109  if ( atmos_sw_phy_bl ) then
110  select case ( atmos_phy_bl_type )
111  case ( 'MYNN' )
113  ka, ks, ke, ia, is, ie, ja, js, je, &
114  bulkflux_type ) ! (in)
115  end select
116  else
117  log_info("ATMOS_PHY_BL_driver_setup",*) 'this component is never called.'
118  atmos_phy_bl_zi(:,:) = 0.0_rp
119  endif
120 
121  return

References scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_setup(), mod_atmos_admin::atmos_phy_bl_type, mod_atmos_phy_bl_vars::atmos_phy_bl_zi, mod_atmos_admin::atmos_sw_phy_bl, scale_bulkflux::bulkflux_type, scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ie, 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_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, and scale_atmos_grid_cartesc_index::ks.

Referenced by mod_atmos_driver::atmos_driver_setup().

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

◆ atmos_phy_bl_driver_calc_tendency()

subroutine, public mod_atmos_phy_bl_driver::atmos_phy_bl_driver_calc_tendency ( logical, intent(in)  update_flag)

calculate tendency

Definition at line 127 of file mod_atmos_phy_bl_driver.F90.

127  use scale_statistics, only: &
129  statistics_total
130  use scale_file_history, only: &
131  file_history_in
132  use scale_time, only: &
133  dt_bl => time_dtsec_atmos_phy_bl
134  use scale_atmos_phy_bl_mynn, only: &
136  atmos_phy_bl_mynn_tendency_tracer
137  use scale_atmos_grid_cartesc_real, only: &
142  use scale_atmos_hydrometeor, only: &
143  i_qv
144  use scale_bulkflux, only: &
146  use mod_atmos_admin, only: &
149  use mod_atmos_vars, only: &
150  dens => dens_av, &
151  qtrc => qtrc_av, &
152  u, &
153  v, &
154  pott, &
155  pres, &
156  exner, &
157  qdry, &
158  qv, &
159  qc, &
160  qi, &
161  rhou_t => rhou_tp, &
162  rhov_t => rhov_tp, &
163  rhot_t => rhot_tp, &
164  rhoq_t => rhoq_tp, &
165  atmos_vars_get_diagnostic
166  use mod_atmos_phy_bl_vars, only: &
167  qs, qe, &
168  rhou_t_bl => atmos_phy_bl_rhou_t, &
169  rhov_t_bl => atmos_phy_bl_rhov_t, &
170  rhot_t_bl => atmos_phy_bl_rhot_t, &
171  rhoq_t_bl => atmos_phy_bl_rhoq_t, &
172  zi => atmos_phy_bl_zi, &
173  ql => atmos_phy_bl_ql, &
174  cldfrac => atmos_phy_bl_cldfrac
175  use mod_atmos_phy_sf_vars, only: &
176  sfc_dens => atmos_phy_sf_sfc_dens, &
177  sflx_mu => atmos_phy_sf_sflx_mu, &
178  sflx_mv => atmos_phy_sf_sflx_mv, &
179  sflx_sh => atmos_phy_sf_sflx_sh, &
180  sflx_q => atmos_phy_sf_sflx_qtrc, &
181  sflx_qv => atmos_phy_sf_sflx_qv, &
182  ustar => atmos_phy_sf_ustar, &
183  tstar => atmos_phy_sf_tstar, &
184  qstar => atmos_phy_sf_qstar, &
185  rlmo => atmos_phy_sf_rlmo
186  implicit none
187 
188  logical, intent(in) :: update_flag
189 
190  real(RP) :: Nu(KA,IA,JA)
191  real(RP) :: Kh(KA,IA,JA)
192 
193  real(RP) :: QW(KA,IA,JA)
194 
195  real(RP) :: N2 (KA,IA,JA)
196  real(RP) :: POTL(KA,IA,JA)
197  real(RP) :: POTV(KA,IA,JA)
198 
199  real(RP), pointer :: RHOQV_t(:,:,:)
200 
201  integer :: k, i, j, iq
202  !---------------------------------------------------------------------------
203 
204  if ( update_flag ) then
205 
206  rhoq_t_bl(:,:,:,:) = 0.0_rp
207 
208  select case ( atmos_phy_bl_type )
209  case ( 'MYNN' )
210  call atmos_vars_get_diagnostic( "N2", n2 )
211  call atmos_vars_get_diagnostic( "POTL", potl )
212  call atmos_vars_get_diagnostic( "POTV", potv )
213  do j = jsb, jeb
214  do i = isb, ieb
215  do k = ks, ke
216  qw(k,i,j) = qv(k,i,j) + qc(k,i,j) + qi(k,i,j)
217  end do
218  end do
219  end do
220  if ( i_qv > 0 ) then
221  rhoqv_t => rhoq_t_bl(:,:,:,i_qv)
222  else
223  allocate( rhoqv_t(ka,ia,ja) )
224  end if
226  ka, ks, ke, ia, is, ie, ja, js, je, &
227  dens(:,:,:), u(:,:,:), v(:,:,:), & ! (in)
228  pott(:,:,:), qtrc(:,:,:,qs:qe), & ! (in)
229  pres(:,:,:), exner(:,:,:), n2(:,:,:), & ! (in)
230  qdry(:,:,:), qv(:,:,:), qw(:,:,:), & ! (in)
231  potl(:,:,:), potv(:,:,:), & ! (in)
232  sfc_dens(:,:), & ! (in)
233  sflx_mu(:,:), sflx_mv(:,:), sflx_sh(:,:), sflx_qv(:,:), & ! (in)
234  ustar(:,:), tstar(:,:), qstar(:,:), rlmo(:,:), & ! (in)
235  cz(:,:,:), fz(:,:,:), dt_bl, & ! (in)
236  bulkflux_type, & ! (in)
237  rhou_t_bl(:,:,:), rhov_t_bl(:,:,:), rhot_t_bl(:,:,:), & ! (out)
238  rhoqv_t(:,:,:), rhoq_t_bl(:,:,:,qs:qe), & ! (out)
239  nu(:,:,:), kh(:,:,:), & ! (out)
240  ql(:,:,:), cldfrac(:,:,:), zi(:,:) ) ! (out)
241  if ( i_qv <= 0 ) deallocate( rhoqv_t )
242  do iq = 1, qa
243  if ( ( .not. tracer_advc(iq) ) .or. iq==i_qv .or. (iq>=qs .and. iq<=qe) ) cycle
244  call atmos_phy_bl_mynn_tendency_tracer( &
245  ka, ks, ke, ia, is, ie, ja, js, je, &
246  dens(:,:,:), qtrc(:,:,:,iq), & ! (in)
247  sflx_q(:,:,iq), & ! (in)
248  kh(:,:,:), & ! (in)
249  tracer_mass(iq), & ! (in)
250  cz(:,:,:), fz(:,:,:), & ! (in)
251  dt_bl, tracer_name(iq), & ! (in)
252  rhoq_t_bl(:,:,:,iq) ) ! (out)
253  end do
254  end select
255 
256  call file_history_in( nu(:,:,:), 'Nu_BL', 'eddy viscosity', 'm2/s', fill_halo=.true., dim_type="ZHXY" )
257  call file_history_in( kh(:,:,:), 'Kh_BL', 'eddy diffusion', 'm2/s', fill_halo=.true., dim_type="ZHXY" )
258 
259  call file_history_in( ql(:,:,:), 'QL_BL', 'liquid water content in partial condensation', 'kg/kg', fill_halo=.true. )
260  call file_history_in( cldfrac(:,:,:), 'cldfrac_BL', 'cloud fraction in partial condensation', '1', fill_halo=.true. )
261 
262  call file_history_in( zi(:,:), 'Zi_BL', 'depth of the boundary layer', 'm', fill_halo=.true. )
263 
264  call file_history_in( rhou_t_bl(:,:,:), 'RHOU_t_BL', 'MOMX tendency (BL)', 'kg/m2/s2', fill_halo=.true. )
265  call file_history_in( rhov_t_bl(:,:,:), 'RHOV_t_BL', 'MOMY tendency (BL)', 'kg/m2/s2', fill_halo=.true. )
266  call file_history_in( rhot_t_bl(:,:,:), 'RHOT_t_BL', 'RHOT tendency (BL)', 'K.kg/m3/s', fill_halo=.true. )
267 
268  do iq = 1, qa
269  if ( .not. tracer_advc(iq) ) cycle
270  call file_history_in( rhoq_t_bl(:,:,:,iq), trim(tracer_name(iq))//'_t_BL', &
271  'RHO*'//trim(tracer_name(iq))//' tendency (BL)', 'kg/m3/s', fill_halo=.true. )
272  enddo
273 
274  if ( statistics_checktotal ) then
275  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
276  rhou_t_bl(:,:,:), 'RHOU_t_BL', &
279  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
280  rhov_t_bl(:,:,:), 'RHOV_t_BL', &
283  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
284  rhot_t_bl(:,:,:), 'RHOT_t_BL', &
287  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
288  nu(:,:,:), 'Nu_BL', &
291  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
292  kh(:,:,:), 'Kh_BL', &
295 
296  do iq = 1, qa
297  if ( .not. tracer_advc(iq) ) cycle
298  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
299  rhoq_t_bl(:,:,:,iq), trim(tracer_name(iq))//'_t_BL', &
302  enddo
303  endif
304 
305  endif
306 
307  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
308  !$omp shared(JS,JE,IS,IE,KS,KE,RHOU_t,RHOU_t_BL,RHOV_t,RHOV_t_BL,RHOT_t,RHOT_t_BL)
309  do j = js, je
310  do i = is, ie
311  do k = ks, ke
312  rhou_t(k,i,j) = rhou_t(k,i,j) + rhou_t_bl(k,i,j)
313  rhov_t(k,i,j) = rhov_t(k,i,j) + rhov_t_bl(k,i,j)
314  rhot_t(k,i,j) = rhot_t(k,i,j) + rhot_t_bl(k,i,j)
315  enddo
316  enddo
317  enddo
318 
319  do iq = 1, qa
320  if ( .not. tracer_advc(iq) ) cycle
321  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(3)
322  do j = js, je
323  do i = is, ie
324  do k = ks, ke
325  rhoq_t(k,i,j,iq) = rhoq_t(k,i,j,iq) + rhoq_t_bl(k,i,j,iq)
326  enddo
327  enddo
328  enddo
329  enddo
330 
331  return

References scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_fz, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvol, scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_vol, mod_atmos_phy_bl_vars::atmos_phy_bl_cldfrac, scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_tendency(), mod_atmos_phy_bl_vars::atmos_phy_bl_ql, mod_atmos_phy_bl_vars::atmos_phy_bl_rhoq_t, mod_atmos_phy_bl_vars::atmos_phy_bl_rhot_t, mod_atmos_phy_bl_vars::atmos_phy_bl_rhou_t, mod_atmos_phy_bl_vars::atmos_phy_bl_rhov_t, mod_atmos_admin::atmos_phy_bl_type, mod_atmos_phy_bl_vars::atmos_phy_bl_zi, mod_atmos_phy_sf_vars::atmos_phy_sf_qstar, mod_atmos_phy_sf_vars::atmos_phy_sf_rlmo, mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_dens, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mu, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mv, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_qtrc, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_qv, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_sh, mod_atmos_phy_sf_vars::atmos_phy_sf_tstar, mod_atmos_phy_sf_vars::atmos_phy_sf_ustar, mod_atmos_admin::atmos_sw_phy_bl, scale_bulkflux::bulkflux_type, mod_atmos_vars::dens, mod_atmos_vars::dens_av, mod_atmos_vars::exner, scale_atmos_hydrometeor::i_qv, scale_atmos_grid_cartesc_index::ia, scale_atmos_grid_cartesc_index::ie, scale_atmos_grid_cartesc_index::ieb, scale_atmos_grid_cartesc_index::is, scale_atmos_grid_cartesc_index::isb, scale_atmos_grid_cartesc_index::ja, scale_atmos_grid_cartesc_index::je, scale_atmos_grid_cartesc_index::jeb, scale_atmos_grid_cartesc_index::js, scale_atmos_grid_cartesc_index::jsb, scale_tracer::k, scale_atmos_grid_cartesc_index::ka, scale_atmos_grid_cartesc_index::ke, scale_atmos_grid_cartesc_index::ks, mod_atmos_vars::pott, mod_atmos_vars::pres, scale_tracer::qa, mod_atmos_vars::qc, mod_atmos_vars::qdry, mod_atmos_phy_bl_vars::qe, mod_atmos_vars::qi, mod_atmos_phy_bl_vars::qs, mod_atmos_vars::qtrc, mod_atmos_vars::qtrc_av, mod_atmos_vars::qv, mod_atmos_vars::rhoq_tp, mod_atmos_vars::rhot_tp, mod_atmos_vars::rhou_tp, mod_atmos_vars::rhov_tp, scale_statistics::statistics_checktotal, scale_time::time_dtsec_atmos_phy_bl, scale_tracer::tracer_advc, scale_tracer::tracer_mass, scale_tracer::tracer_name, mod_atmos_vars::u, and mod_atmos_vars::v.

Referenced by mod_atmos_driver::atmos_driver_calc_tendency(), and mod_atmos_driver::atmos_driver_calc_tendency_from_sflux().

Here is the call graph for this function:
Here is the caller graph for this function:
mod_atmos_admin::atmos_phy_bl_type
character(len=h_short), public atmos_phy_bl_type
Definition: mod_atmos_admin.F90:42
scale_statistics
module Statistics
Definition: scale_statistics.F90:11
mod_atmos_phy_sf_vars::atmos_phy_sf_sfc_dens
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sfc_dens
Definition: mod_atmos_phy_sf_vars.F90:71
scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_desc
character(len=h_long), dimension(:), allocatable, public atmos_phy_bl_mynn_desc
Definition: scale_atmos_phy_bl_mynn.F90:52
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mu
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mu
Definition: mod_atmos_phy_sf_vars.F90:77
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
mod_atmos_vars::rhoq_tp
real(rp), dimension(:,:,:,:), allocatable, public rhoq_tp
Definition: mod_atmos_vars.F90:120
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_cz
geopotential height [m] (zxy)
Definition: scale_atmos_grid_cartesC_real.F90:38
mod_atmos_vars::pott
real(rp), dimension(:,:,:), allocatable, target, public pott
Definition: mod_atmos_vars.F90:132
scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_tendency
subroutine, public atmos_phy_bl_mynn_tendency(KA, KS, KE, IA, IS, IE, JA, JS, JE, DENS, U, V, POTT, PROG, PRES, EXNER, N2, QDRY, QV, Qw, POTL, POTV, SFC_DENS, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_QV, us, ts, qs, RLmo, CZ, FZ, dt_DP, BULKFLUX_type, RHOU_t, RHOV_t, RHOT_t, RHOQV_t, RPROG_t, Nu, Kh, Qlp, cldfrac, Zi)
ATMOS_PHY_BL_MYNN_tendency calculate tendency by the virtical eddy viscosity.
Definition: scale_atmos_phy_bl_mynn.F90:276
mod_atmos_vars::qtrc_av
real(rp), dimension(:,:,:,:), pointer, public qtrc_av
Definition: mod_atmos_vars.F90:94
scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_units
character(len=h_short), dimension(:), allocatable, public atmos_phy_bl_mynn_units
Definition: scale_atmos_phy_bl_mynn.F90:53
mod_atmos_vars::rhov_tp
real(rp), dimension(:,:,:), allocatable, public rhov_tp
Definition: mod_atmos_vars.F90:117
scale_atmos_phy_bl_mynn
module atmosphere / physics / pbl / mynn
Definition: scale_atmos_phy_bl_mynn.F90:18
mod_atmos_admin
module ATMOS admin
Definition: mod_atmos_admin.F90:11
mod_atmos_phy_sf_vars::atmos_phy_sf_qstar
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_qstar
Definition: mod_atmos_phy_sf_vars.F90:90
mod_atmos_vars::qdry
real(rp), dimension(:,:,:), allocatable, target, public qdry
Definition: mod_atmos_vars.F90:139
mod_atmos_phy_sf_vars::atmos_phy_sf_ustar
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_ustar
Definition: mod_atmos_phy_sf_vars.F90:88
scale_time::time_dtsec_atmos_phy_bl
real(dp), public time_dtsec_atmos_phy_bl
time interval of physics(pbl ) [sec]
Definition: scale_time.F90:42
mod_atmos_phy_bl_vars::atmos_phy_bl_rhou_t
real(rp), dimension(:,:,:), allocatable, public atmos_phy_bl_rhou_t
Definition: mod_atmos_phy_bl_vars.F90:58
scale_bulkflux
module Surface bulk flux
Definition: scale_bulkflux.F90:12
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mv
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mv
Definition: mod_atmos_phy_sf_vars.F90:78
mod_atmos_phy_bl_vars::atmos_phy_bl_rhoq_t
real(rp), dimension(:,:,:,:), allocatable, target, public atmos_phy_bl_rhoq_t
Definition: mod_atmos_phy_bl_vars.F90:62
mod_atmos_phy_sf_vars
module ATMOSPHERIC Surface Variables
Definition: mod_atmos_phy_sf_vars.F90:12
scale_atmos_grid_cartesc_real
module Atmosphere GRID CartesC Real(real space)
Definition: scale_atmos_grid_cartesC_real.F90:11
mod_atmos_vars::qtrc
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc
Definition: mod_atmos_vars.F90:80
scale_file_history
module file_history
Definition: scale_file_history.F90:15
scale_prc
module PROCESS
Definition: scale_prc.F90:11
mod_atmos_phy_bl_vars::atmos_phy_bl_rhov_t
real(rp), dimension(:,:,:), allocatable, public atmos_phy_bl_rhov_t
Definition: mod_atmos_phy_bl_vars.F90:59
mod_atmos_vars::rhou_tp
real(rp), dimension(:,:,:), allocatable, public rhou_tp
Definition: mod_atmos_vars.F90:116
mod_atmos_phy_bl_vars
module atmosphere / physics / PBL
Definition: mod_atmos_phy_bl_vars.F90:12
mod_atmos_vars::dens
real(rp), dimension(:,:,:), allocatable, target, public dens
Definition: mod_atmos_vars.F90:75
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_vol
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_vol
control volume (zxy) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:83
mod_atmos_phy_sf_vars::atmos_phy_sf_rlmo
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_rlmo
Definition: mod_atmos_phy_sf_vars.F90:98
mod_atmos_vars::v
real(rp), dimension(:,:,:), allocatable, target, public v
Definition: mod_atmos_vars.F90:130
mod_atmos_phy_bl_vars::atmos_phy_bl_zi
real(rp), dimension(:,:), allocatable, public atmos_phy_bl_zi
Definition: mod_atmos_phy_bl_vars.F90:64
mod_atmos_vars::exner
real(rp), dimension(:,:,:), allocatable, target, public exner
Definition: mod_atmos_vars.F90:135
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_qv
real(rp), dimension(:,:), pointer, public atmos_phy_sf_sflx_qv
Definition: mod_atmos_phy_sf_vars.F90:86
scale_bulkflux::bulkflux_type
character(len=h_short), public bulkflux_type
Definition: scale_bulkflux.F90:85
scale_time
module TIME
Definition: scale_time.F90:11
mod_atmos_vars::qv
real(rp), dimension(:,:,:), allocatable, pointer, target, public qv
Definition: mod_atmos_vars.F90:96
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_qtrc
real(rp), dimension(:,:,:), allocatable, target, public atmos_phy_sf_sflx_qtrc
Definition: mod_atmos_phy_sf_vars.F90:84
scale_tracer::tracer_regist
subroutine, public tracer_regist(QS, NQ, NAME, DESC, UNIT, CV, CP, R, ENGI0, ADVC, MASS)
Regist tracer.
Definition: scale_tracer.F90:65
mod_atmos_phy_bl_vars::qe
integer, public qe
Definition: mod_atmos_phy_bl_vars.F90:45
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
scale_atmos_hydrometeor::i_qv
integer, public i_qv
Definition: scale_atmos_hydrometeor.F90:77
mod_atmos_vars::pres
real(rp), dimension(:,:,:), allocatable, target, public pres
Definition: mod_atmos_vars.F90:134
mod_atmos_vars::dens_av
real(rp), dimension(:,:,:), pointer, public dens_av
Definition: mod_atmos_vars.F90:89
mod_atmos_vars::u
real(rp), dimension(:,:,:), allocatable, target, public u
Definition: mod_atmos_vars.F90:129
scale_statistics::statistics_checktotal
logical, public statistics_checktotal
calc&report variable totals to logfile?
Definition: scale_statistics.F90:64
mod_atmos_vars
module ATMOSPHERIC Variables
Definition: mod_atmos_vars.F90:12
mod_atmos_vars::qi
real(rp), dimension(:,:,:), pointer, public qi
Definition: mod_atmos_vars.F90:99
mod_atmos_vars::rhot_tp
real(rp), dimension(:,:,:), allocatable, public rhot_tp
Definition: mod_atmos_vars.F90:118
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_sh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_sh
Definition: mod_atmos_phy_sf_vars.F90:79
mod_atmos_admin::atmos_sw_phy_bl
logical, public atmos_sw_phy_bl
Definition: mod_atmos_admin.F90:58
scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_ntracer
integer, public atmos_phy_bl_mynn_ntracer
Definition: scale_atmos_phy_bl_mynn.F90:50
mod_atmos_vars::qc
real(rp), dimension(:,:,:), pointer, public qc
Definition: mod_atmos_vars.F90:97
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_fz
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_fz
geopotential height [m] (wxy)
Definition: scale_atmos_grid_cartesC_real.F90:42
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvol
real(rp), public atmos_grid_cartesc_real_totvol
total volume (zxy, local) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:87
mod_atmos_phy_bl_vars::atmos_phy_bl_cldfrac
real(rp), dimension(:,:,:), allocatable, public atmos_phy_bl_cldfrac
Definition: mod_atmos_phy_bl_vars.F90:67
mod_atmos_phy_bl_vars::atmos_phy_bl_rhot_t
real(rp), dimension(:,:,:), allocatable, public atmos_phy_bl_rhot_t
Definition: mod_atmos_phy_bl_vars.F90:60
scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_setup
subroutine, public atmos_phy_bl_mynn_setup(KA, KS, KE, IA, IS, IE, JA, JS, JE, BULKFLUX_type, TKE_MIN, PBL_MAX)
ATMOS_PHY_BL_MYNN_setup Setup.
Definition: scale_atmos_phy_bl_mynn.F90:190
mod_atmos_phy_bl_vars::atmos_phy_bl_ql
real(rp), dimension(:,:,:), allocatable, public atmos_phy_bl_ql
Definition: mod_atmos_phy_bl_vars.F90:66
scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_name
character(len=h_short), dimension(:), allocatable, public atmos_phy_bl_mynn_name
Definition: scale_atmos_phy_bl_mynn.F90:51
mod_atmos_phy_sf_vars::atmos_phy_sf_tstar
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_tstar
Definition: mod_atmos_phy_sf_vars.F90:89
mod_atmos_phy_bl_vars::qs
integer, public qs
Definition: mod_atmos_phy_bl_vars.F90:45
scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_tracer_setup
subroutine, public atmos_phy_bl_mynn_tracer_setup
ATMOS_PHY_BL_MYNN_tracer_setup Tracer Setup.
Definition: scale_atmos_phy_bl_mynn.F90:134