SCALE-RM
Functions/Subroutines
mod_atmos_phy_tb_driver Module Reference

module ATMOSPHERE / Physics Turbulence More...

Functions/Subroutines

subroutine, public atmos_phy_tb_driver_setup
 Setup. More...
 
subroutine, public atmos_phy_tb_driver_resume
 Resume. More...
 
subroutine, public atmos_phy_tb_driver (update_flag)
 Driver. More...
 

Detailed Description

module ATMOSPHERE / Physics Turbulence

Description
Sub-grid scale turbulence process
Author
Team SCALE
History
  • 2013-12-05 (S.Nishizawa) [new]
NAMELIST
  • No namelist group
History Output
namedescriptionunitvariable
MOMX_t_TB MOMX tendency (TB) kg/m2/s2 MOMX_t_TB
MOMY_t_TB MOMY tendency (TB) kg/m2/s2 MOMY_t_TB
MOMZ_t_TB MOMZ tendency (TB) kg/m2/s2 MOMZ_t_TB
N2 squared Brunt-Vaisala frequency 1/s2 N2
NU eddy viscosity m2/s NU
Pr Prantle number NIL Pr
RHOT_t_TB RHOT tendency (TB) K.kg/m3/s RHOT_t_TB
Ri Richardson number NIL Ri
SGS_XFLX_MOMX SGS X FLUX of MOMX kg/m/s2 QFLX_MOMX
SGS_XFLX_MOMY SGS X FLUX of MOMY kg/m/s2 QFLX_MOMY
SGS_XFLX_MOMZ SGS X FLUX of MOMZ kg/m/s2 QFLX_MOMZ
SGS_XFLX_QC SGS X FLUX of QC kg/m2/s QFLX_RHOQ
SGS_XFLX_QG SGS X FLUX of QG kg/m2/s QFLX_RHOQ
SGS_XFLX_QI SGS X FLUX of QI kg/m2/s QFLX_RHOQ
SGS_XFLX_QR SGS X FLUX of QR kg/m2/s QFLX_RHOQ
SGS_XFLX_QS SGS X FLUX of QS kg/m2/s QFLX_RHOQ
SGS_XFLX_QV SGS X FLUX of QV kg/m2/s QFLX_RHOQ
SGS_XFLX_RHOT SGS X FLUX of RHOT K*kg/m2/s QFLX_RHOT
SGS_YFLX_MOMX SGS Y FLUX of MOMX kg/m/s2 QFLX_MOMX
SGS_YFLX_MOMY SGS Y FLUX of MOMY kg/m/s2 QFLX_MOMY
SGS_YFLX_MOMZ SGS Y FLUX of MOMZ kg/m/s2 QFLX_MOMZ
SGS_YFLX_QC SGS Y FLUX of QC kg/m2/s QFLX_RHOQ
SGS_YFLX_QG SGS Y FLUX of QG kg/m2/s QFLX_RHOQ
SGS_YFLX_QI SGS Y FLUX of QI kg/m2/s QFLX_RHOQ
SGS_YFLX_QR SGS Y FLUX of QR kg/m2/s QFLX_RHOQ
SGS_YFLX_QS SGS Y FLUX of QS kg/m2/s QFLX_RHOQ
SGS_YFLX_QV SGS Y FLUX of QV kg/m2/s QFLX_RHOQ
SGS_YFLX_RHOT SGS Y FLUX of RHOT K*kg/m2/s QFLX_RHOT
SGS_ZFLX_MOMX SGS Z FLUX of MOMX kg/m/s2 QFLX_MOMX
SGS_ZFLX_MOMY SGS Z FLUX of MOMY kg/m/s2 QFLX_MOMY
SGS_ZFLX_MOMZ SGS Z FLUX of MOMZ kg/m/s2 QFLX_MOMZ
SGS_ZFLX_QC SGS Z FLUX of QC kg/m2/s QFLX_RHOQ
SGS_ZFLX_QG SGS Z FLUX of QG kg/m2/s QFLX_RHOQ
SGS_ZFLX_QI SGS Z FLUX of QI kg/m2/s QFLX_RHOQ
SGS_ZFLX_QR SGS Z FLUX of QR kg/m2/s QFLX_RHOQ
SGS_ZFLX_QS SGS Z FLUX of QS kg/m2/s QFLX_RHOQ
SGS_ZFLX_QV SGS Z FLUX of QV kg/m2/s QFLX_RHOQ
SGS_ZFLX_RHOT SGS Z FLUX of RHOT K*kg/m2/s QFLX_RHOT
TKE turburent kinetic energy m2/s2 TKE
trim(AQ_NAME(iq))//'_t_TB' RHO*'//trim(AQ_NAME(iq))//' tendency (TB) kg/m3/s RHOQ_t_TB

Function/Subroutine Documentation

◆ atmos_phy_tb_driver_setup()

subroutine, public mod_atmos_phy_tb_driver::atmos_phy_tb_driver_setup ( )

Setup.

Definition at line 52 of file mod_atmos_phy_tb_driver.f90.

References mod_atmos_phy_tb_vars::atmos_phy_tb_momz_t, scale_atmos_phy_tb::atmos_phy_tb_setup(), mod_atmos_phy_tb_vars::atmos_phy_tb_tke_t, mod_atmos_admin::atmos_phy_tb_type, mod_atmos_admin::atmos_sw_phy_tb, scale_grid::grid_cdx, scale_grid::grid_cdy, scale_grid::grid_cdz, 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, scale_grid_index::ks, scale_process::prc_mpistop(), and scale_grid_real::real_cz.

Referenced by mod_atmos_driver::atmos_driver_setup().

52  use scale_grid, only: &
53  cdz => grid_cdz, &
54  cdx => grid_cdx, &
55  cdy => grid_cdy
56  use scale_grid_real, only: &
57  cz => real_cz
58  use scale_process, only: &
60  use mod_atmos_admin, only: &
63  use scale_atmos_phy_tb, only: &
65  use mod_atmos_phy_tb_vars, only: &
66  tke_t_tb => atmos_phy_tb_tke_t, &
67  momz_t_tb => atmos_phy_tb_momz_t
68  implicit none
69 
70 ! NAMELIST / PARAM_ATMOS_PHY_TB / &
71 
72  integer :: k, i, j
73 ! integer :: ierr
74  !---------------------------------------------------------------------------
75 
76  if( io_l ) write(io_fid_log,*)
77  if( io_l ) write(io_fid_log,*) '++++++ Module[DRIVER] / Categ[ATMOS PHY_TB] / Origin[SCALE-RM]'
78 
79 ! !--- read namelist
80 ! rewind(IO_FID_CONF)
81 ! read(IO_FID_CONF,nml=PARAM_ATMOS_PHY_TB,iostat=ierr)
82 ! if( ierr < 0 ) then !--- missing
83 ! if( IO_L ) write(IO_FID_LOG,*) '*** Not found namelist. Default used.'
84 ! elseif( ierr > 0 ) then !--- fatal error
85 ! write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_TB. Check!'
86 ! call PRC_MPIstop
87 ! endif
88 ! if( IO_LNML ) write(IO_FID_LOG,nml=PARAM_ATMOS_PHY_TB)
89 
90 
91  ! initialize
92  do j = js, je
93  do i = is, ie
94  momz_t_tb(ks-1,i,j) = 0.0_rp
95  momz_t_tb(ke ,i,j) = 0.0_rp
96  enddo
97  enddo
98 
99  do j = js, je
100  do i = is, ie
101  do k = ks, ke
102  tke_t_tb(k,i,j) = 0.0_rp
103  enddo
104  enddo
105  enddo
106 
107  if ( atmos_sw_phy_tb ) then
108 
109  ! setup library component
110  call atmos_phy_tb_setup( atmos_phy_tb_type, & ! [IN]
111  cdz, cdx, cdy, cz ) ! [IN]
112  else
113  if( io_l ) write(io_fid_log,*) '*** this component is never called.'
114  endif
115 
116  return
module ATMOS admin
subroutine, public prc_mpistop
Abort MPI.
real(rp), dimension(:,:,:), allocatable, public real_cz
geopotential height [m] (cell center)
logical, public atmos_sw_phy_tb
module GRID (real space)
module PROCESS
module Atmosphere / Physics Turbulence
subroutine, public atmos_phy_tb_setup(TB_TYPE, CDZ, CDX, CDY, CZ)
character(len=h_short), public atmos_phy_tb_type
module ATMOSPHERE / Physics Turbulence
module GRID (cartesian)
real(rp), dimension(:), allocatable, public grid_cdz
z-length of control volume [m]
real(rp), dimension(:), allocatable, public grid_cdy
y-length of control volume [m]
real(rp), dimension(:,:,:), allocatable, public atmos_phy_tb_tke_t
real(rp), dimension(:), allocatable, public grid_cdx
x-length of control volume [m]
real(rp), dimension(:,:,:), allocatable, public atmos_phy_tb_momz_t
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_tb_driver_resume()

subroutine, public mod_atmos_phy_tb_driver::atmos_phy_tb_driver_resume ( )

Resume.

Definition at line 122 of file mod_atmos_phy_tb_driver.f90.

References atmos_phy_tb_driver(), mod_atmos_admin::atmos_sw_phy_tb, scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

Referenced by mod_atmos_driver::atmos_driver_resume2().

122  use mod_atmos_admin, only: &
124  implicit none
125 
126  if ( atmos_sw_phy_tb ) then
127 
128  ! run once (only for the diagnostic value)
129  call prof_rapstart('ATM_Turbulence', 1)
130  call atmos_phy_tb_driver( update_flag = .true. )
131  call prof_rapend ('ATM_Turbulence', 1)
132 
133  end if
134 
135  return
module ATMOS admin
logical, public atmos_sw_phy_tb
subroutine, public atmos_phy_tb_driver(update_flag)
Driver.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_tb_driver()

subroutine, public mod_atmos_phy_tb_driver::atmos_phy_tb_driver ( logical, intent(in)  update_flag)

Driver.

Definition at line 141 of file mod_atmos_phy_tb_driver.f90.

References scale_tracer::aq_name, 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_mw, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_qtrc, mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_sh, scale_atmos_phy_tb::atmos_phy_tb, scale_atmos_phy_tb_common::atmos_phy_tb_calc_tend_momx(), scale_atmos_phy_tb_common::atmos_phy_tb_calc_tend_momy(), scale_atmos_phy_tb_common::atmos_phy_tb_calc_tend_momz(), scale_atmos_phy_tb_common::atmos_phy_tb_calc_tend_phi(), mod_atmos_phy_tb_vars::atmos_phy_tb_momx_t, mod_atmos_phy_tb_vars::atmos_phy_tb_momy_t, mod_atmos_phy_tb_vars::atmos_phy_tb_momz_t, mod_atmos_phy_tb_vars::atmos_phy_tb_nu, mod_atmos_phy_tb_vars::atmos_phy_tb_rhoq_t, mod_atmos_phy_tb_vars::atmos_phy_tb_rhot_t, mod_atmos_phy_tb_vars::atmos_phy_tb_tke, mod_atmos_phy_tb_vars::atmos_phy_tb_tke_t, mod_atmos_vars::dens, mod_atmos_vars::dens_av, scale_gridtrans::gtrans_gsqrt, scale_gridtrans::gtrans_j13g, scale_gridtrans::gtrans_j23g, scale_gridtrans::gtrans_j33g, scale_gridtrans::gtrans_mapf, scale_tracer::i_qc, scale_tracer::i_qg, scale_tracer::i_qi, scale_tracer::i_qr, scale_tracer::i_qs, scale_tracer::i_qv, scale_grid_index::iblock, scale_grid_index::ie, scale_grid_index::is, scale_grid_index::jblock, scale_grid_index::je, scale_grid_index::js, scale_grid_index::ka, scale_grid_index::ke, scale_grid_index::ks, mod_atmos_vars::momx, mod_atmos_vars::momx_av, mod_atmos_vars::momx_tp, mod_atmos_vars::momy, mod_atmos_vars::momy_av, mod_atmos_vars::momy_tp, mod_atmos_vars::momz, mod_atmos_vars::momz_av, mod_atmos_vars::momz_tp, scale_tracer::qa, mod_atmos_vars::qtrc, mod_atmos_vars::qtrc_av, mod_atmos_vars::rhoq_tp, mod_atmos_vars::rhot, mod_atmos_vars::rhot_av, mod_atmos_vars::rhot_tp, scale_rm_statistics::statistics_checktotal, scale_time::time_dtsec_atmos_phy_tb, scale_grid_index::xdir, scale_grid_index::ydir, and scale_grid_index::zdir.

Referenced by mod_atmos_driver::atmos_driver(), and atmos_phy_tb_driver_resume().

141  use scale_gridtrans, only: &
142  gsqrt => gtrans_gsqrt, &
143  j13g => gtrans_j13g, &
144  j23g => gtrans_j23g, &
145  j33g => gtrans_j33g, &
146  mapf => gtrans_mapf
147  use scale_rm_statistics, only: &
149  stat_total
150  use scale_comm, only: &
151  comm_vars8, &
152  comm_wait
153  use scale_history, only: &
154  hist_in
155  use scale_time, only: &
156  dt_tb => time_dtsec_atmos_phy_tb
157  use scale_atmos_phy_tb, only: &
159  use scale_atmos_phy_tb_common, only: &
160  calc_tend_momz => atmos_phy_tb_calc_tend_momz, &
161  calc_tend_momx => atmos_phy_tb_calc_tend_momx, &
162  calc_tend_momy => atmos_phy_tb_calc_tend_momy, &
163  calc_tend_phi => atmos_phy_tb_calc_tend_phi
164  use mod_atmos_vars, only: &
165  dens => dens_av, &
166  momz => momz_av, &
167  momx => momx_av, &
168  momy => momy_av, &
169  rhot => rhot_av, &
170  qtrc => qtrc_av, &
171  momz_t => momz_tp, &
172  momx_t => momx_tp, &
173  momy_t => momy_tp, &
174  rhot_t => rhot_tp, &
175  rhoq_t => rhoq_tp
176  use mod_atmos_phy_tb_vars, only: &
177  momz_t_tb => atmos_phy_tb_momz_t, &
178  momx_t_tb => atmos_phy_tb_momx_t, &
179  momy_t_tb => atmos_phy_tb_momy_t, &
180  rhot_t_tb => atmos_phy_tb_rhot_t, &
181  rhoq_t_tb => atmos_phy_tb_rhoq_t, &
182  tke_t_tb => atmos_phy_tb_tke_t, &
183  tke => atmos_phy_tb_tke, &
184  nu => atmos_phy_tb_nu
185  use mod_atmos_phy_sf_vars, only: &
186  sflx_mw => atmos_phy_sf_sflx_mw, &
187  sflx_mu => atmos_phy_sf_sflx_mu, &
188  sflx_mv => atmos_phy_sf_sflx_mv, &
189  sflx_sh => atmos_phy_sf_sflx_sh, &
190  sflx_q => atmos_phy_sf_sflx_qtrc
191  implicit none
192 
193  logical, intent(in) :: update_flag
194 
195  ! eddy viscosity/diffusion flux
196  real(RP) :: qflx_momz(ka,ia,ja,3)
197  real(RP) :: qflx_momx(ka,ia,ja,3)
198  real(RP) :: qflx_momy(ka,ia,ja,3)
199  real(RP) :: qflx_rhot(ka,ia,ja,3)
200  real(RP) :: qflx_rhoq(ka,ia,ja,3,qa)
201 
202  real(RP) :: ri(ka,ia,ja) ! Richardson number
203  real(RP) :: pr(ka,ia,ja) ! Prandtl number
204  real(RP) :: n2(ka,ia,ja) ! squared Brunt-Vaisala frequency
205 
206  integer :: jjs, jje
207  integer :: iis, iie
208 
209  real(RP) :: total ! dummy
210 
211  integer :: k, i, j, iq
212  !---------------------------------------------------------------------------
213 
214  if ( update_flag ) then
215 
216  do j = js, je
217  do i = is, ie
218  do k = ks, ke
219  tke(k,i,j) = max( tke(k,i,j) + tke_t_tb(k,i,j) * dt_tb, 0.0_rp )
220  end do
221  end do
222  end do
223 
224  do j = js, je
225  do i = is, ie
226  tke(1:ks-1,i,j) = tke(ks,i,j)
227  tke(ke+1:ka,i,j) = tke(ke,i,j)
228  end do
229  end do
230 
231  call comm_vars8( tke(:,:,:), 1 )
232  call comm_wait ( tke(:,:,:), 1 )
233 
234  call atmos_phy_tb( qflx_momz, & ! [OUT]
235  qflx_momx, & ! [OUT]
236  qflx_momy, & ! [OUT]
237  qflx_rhot, & ! [OUT]
238  qflx_rhoq, & ! [OUT]
239  tke, & ! [INOUT]
240  tke_t_tb, & ! [OUT]
241  nu, & ! [OUT]
242  ri, & ! [OUT]
243  pr, & ! [OUT]
244  n2, & ! [OUT]
245  momz, & ! [IN]
246  momx, & ! [IN]
247  momy, & ! [IN]
248  rhot, & ! [IN]
249  dens, & ! [IN]
250  qtrc, & ! [IN]
251  sflx_mw, & ! [IN]
252  sflx_mu, & ! [IN]
253  sflx_mv, & ! [IN]
254  sflx_sh, & ! [IN]
255  sflx_q(:,:,i_qv), & ! [IN]
256  gsqrt, & ! [IN]
257  j13g, & ! [IN]
258  j23g, & ! [IN]
259  j33g, & ! [IN]
260  mapf, & ! [IN]
261  dt_tb ) ! [IN]
262 
263  do jjs = js, je, jblock
264  jje = jjs+jblock-1
265  do iis = is, ie, iblock
266  iie = iis+iblock-1
267  call calc_tend_momz( momz_t_tb, & ! (out)
268  qflx_momz, & ! (in)
269  gsqrt, j13g, j23g, j33g, mapf, & ! (in)
270  iis, iie, jjs, jje ) ! (in)
271  end do
272  end do
273 
274  do jjs = js, je, jblock
275  jje = jjs+jblock-1
276  do iis = is, ie, iblock
277  iie = iis+iblock-1
278  call calc_tend_momx( momx_t_tb, & ! (out)
279  qflx_momx, & ! (in)
280  gsqrt, j13g, j23g, j33g, mapf, & ! (in)
281  iis, iie, jjs, jje ) ! (in)
282  end do
283  end do
284 
285  do jjs = js, je, jblock
286  jje = jjs+jblock-1
287  do iis = is, ie, iblock
288  iie = iis+iblock-1
289  call calc_tend_momy( momy_t_tb, & ! (out)
290  qflx_momy, & ! (in)
291  gsqrt, j13g, j23g, j33g, mapf, & ! (in)
292  iis, iie, jjs, jje ) ! (in)
293  end do
294  end do
295 
296  do jjs = js, je, jblock
297  jje = jjs+jblock-1
298  do iis = is, ie, iblock
299  iie = iis+iblock-1
300  call calc_tend_phi ( rhot_t_tb, & ! (out)
301  qflx_rhot, & ! (in)
302  gsqrt, j13g, j23g, j33g, mapf, & ! (in)
303  iis, iie, jjs, jje ) ! (in)
304  end do
305  end do
306 
307  do iq = 1, qa
308  do jjs = js, je, jblock
309  jje = jjs+jblock-1
310  do iis = is, ie, iblock
311  iie = iis+iblock-1
312  call calc_tend_phi ( rhoq_t_tb(:,:,:,iq), & ! (out)
313  qflx_rhoq(:,:,:,:,iq), & ! (in)
314  gsqrt, j13g, j23g, j33g, mapf, & ! (in)
315  iis, iie, jjs, jje ) ! (in)
316  end do
317  end do
318  end do
319 
320  call hist_in( tke(:,:,:), 'TKE', 'turburent kinetic energy', 'm2/s2', nohalo=.true. )
321  call hist_in( nu(:,:,:), 'NU', 'eddy viscosity', 'm2/s' , nohalo=.true. )
322  call hist_in( ri(:,:,:), 'Ri', 'Richardson number', 'NIL' , nohalo=.true. )
323  call hist_in( pr(:,:,:), 'Pr', 'Prantle number', 'NIL' , nohalo=.true. )
324  call hist_in( n2(:,:,:), 'N2', 'squared Brunt-Vaisala frequency', '1/s2', nohalo=.true. )
325 
326  call hist_in( momz_t_tb(:,:,:), 'MOMZ_t_TB', 'MOMZ tendency (TB)', 'kg/m2/s2', nohalo=.true. )
327  call hist_in( momx_t_tb(:,:,:), 'MOMX_t_TB', 'MOMX tendency (TB)', 'kg/m2/s2', nohalo=.true. )
328  call hist_in( momy_t_tb(:,:,:), 'MOMY_t_TB', 'MOMY tendency (TB)', 'kg/m2/s2', nohalo=.true. )
329  call hist_in( rhot_t_tb(:,:,:), 'RHOT_t_TB', 'RHOT tendency (TB)', 'K.kg/m3/s', nohalo=.true. )
330 
331  do iq = 1, qa
332  call hist_in( rhoq_t_tb(:,:,:,iq), trim(aq_name(iq))//'_t_TB', &
333  'RHO*'//trim(aq_name(iq))//' tendency (TB)', 'kg/m3/s', nohalo=.true. )
334  enddo
335 
336  call hist_in( qflx_momz(:,:,:,zdir), 'SGS_ZFLX_MOMZ', 'SGS Z FLUX of MOMZ', 'kg/m/s2', &
337  nohalo=.true.)
338  call hist_in( qflx_momz(:,:,:,xdir), 'SGS_XFLX_MOMZ', 'SGS X FLUX of MOMZ', 'kg/m/s2', &
339  xdim='half', zdim='half', nohalo=.true.)
340  call hist_in( qflx_momz(:,:,:,ydir), 'SGS_YFLX_MOMZ', 'SGS Y FLUX of MOMZ', 'kg/m/s2', &
341  ydim='half', zdim='half', nohalo=.true.)
342 
343  call hist_in( qflx_momx(:,:,:,zdir), 'SGS_ZFLX_MOMX', 'SGS Z FLUX of MOMX', 'kg/m/s2', &
344  xdim='half', zdim='half', nohalo=.true.)
345  call hist_in( qflx_momx(:,:,:,xdir), 'SGS_XFLX_MOMX', 'SGS X FLUX of MOMX', 'kg/m/s2', &
346  nohalo=.true.)
347  call hist_in( qflx_momx(:,:,:,ydir), 'SGS_YFLX_MOMX', 'SGS Y FLUX of MOMX', 'kg/m/s2', &
348  xdim='half', ydim='half', nohalo=.true.)
349 
350  call hist_in( qflx_momy(:,:,:,zdir), 'SGS_ZFLX_MOMY', 'SGS Z FLUX of MOMY', 'kg/m/s2', &
351  ydim='half', zdim='half', nohalo=.true.)
352  call hist_in( qflx_momy(:,:,:,xdir), 'SGS_XFLX_MOMY', 'SGS X FLUX of MOMY', 'kg/m/s2', &
353  xdim='half', ydim='half', nohalo=.true.)
354  call hist_in( qflx_momy(:,:,:,ydir), 'SGS_YFLX_MOMY', 'SGS Y FLUX of MOMY', 'kg/m/s2', &
355  nohalo=.true.)
356 
357  call hist_in( qflx_rhot(:,:,:,zdir), 'SGS_ZFLX_RHOT', 'SGS Z FLUX of RHOT', 'K*kg/m2/s', &
358  zdim='half', nohalo=.true.)
359  call hist_in( qflx_rhot(:,:,:,xdir), 'SGS_XFLX_RHOT', 'SGS X FLUX of RHOT', 'K*kg/m2/s', &
360  xdim='half', nohalo=.true.)
361  call hist_in( qflx_rhot(:,:,:,ydir), 'SGS_YFLX_RHOT', 'SGS Y FLUX of RHOT', 'K*kg/m2/s', &
362  ydim='half', nohalo=.true.)
363 
364  if ( i_qv > 0 ) then
365  call hist_in( qflx_rhoq(:,:,:,zdir,i_qv), 'SGS_ZFLX_QV', 'SGS Z FLUX of QV', 'kg/m2/s', &
366  zdim='half', nohalo=.true.)
367  call hist_in( qflx_rhoq(:,:,:,xdir,i_qv), 'SGS_XFLX_QV', 'SGS X FLUX of QV', 'kg/m2/s', &
368  xdim='half', nohalo=.true.)
369  call hist_in( qflx_rhoq(:,:,:,ydir,i_qv), 'SGS_YFLX_QV', 'SGS Y FLUX of QV', 'kg/m2/s', &
370  ydim='half', nohalo=.true.)
371  endif
372 
373 #ifndef DRY
374  if ( i_qc > 0 ) then
375  call hist_in( qflx_rhoq(:,:,:,zdir,i_qc), 'SGS_ZFLX_QC', 'SGS Z FLUX of QC', 'kg/m2/s', &
376  zdim='half', nohalo=.true.)
377  call hist_in( qflx_rhoq(:,:,:,xdir,i_qc), 'SGS_XFLX_QC', 'SGS X FLUX of QC', 'kg/m2/s', &
378  xdim='half', nohalo=.true.)
379  call hist_in( qflx_rhoq(:,:,:,ydir,i_qc), 'SGS_YFLX_QC', 'SGS Y FLUX of QC', 'kg/m2/s', &
380  ydim='half', nohalo=.true.)
381  endif
382 
383  if ( i_qr > 0 ) then
384  call hist_in( qflx_rhoq(:,:,:,zdir,i_qr), 'SGS_ZFLX_QR', 'SGS Z FLUX of QR', 'kg/m2/s', &
385  zdim='half', nohalo=.true.)
386  call hist_in( qflx_rhoq(:,:,:,xdir,i_qr), 'SGS_XFLX_QR', 'SGS X FLUX of QR', 'kg/m2/s', &
387  xdim='half', nohalo=.true.)
388  call hist_in( qflx_rhoq(:,:,:,ydir,i_qr), 'SGS_YFLX_QR', 'SGS Y FLUX of QR', 'kg/m2/s', &
389  ydim='half', nohalo=.true.)
390  endif
391 
392  if ( i_qi > 0 ) then
393  call hist_in( qflx_rhoq(:,:,:,zdir,i_qi), 'SGS_ZFLX_QI', 'SGS Z FLUX of QI', 'kg/m2/s', &
394  zdim='half', nohalo=.true.)
395  call hist_in( qflx_rhoq(:,:,:,xdir,i_qi), 'SGS_XFLX_QI', 'SGS X FLUX of QI', 'kg/m2/s', &
396  xdim='half', nohalo=.true.)
397  call hist_in( qflx_rhoq(:,:,:,ydir,i_qi), 'SGS_YFLX_QI', 'SGS Y FLUX of QI', 'kg/m2/s', &
398  ydim='half', nohalo=.true.)
399  endif
400 
401  if ( i_qs > 0 ) then
402  call hist_in( qflx_rhoq(:,:,:,zdir,i_qs), 'SGS_ZFLX_QS', 'SGS Z FLUX of QS', 'kg/m2/s', &
403  zdim='half', nohalo=.true.)
404  call hist_in( qflx_rhoq(:,:,:,xdir,i_qs), 'SGS_XFLX_QS', 'SGS X FLUX of QS', 'kg/m2/s', &
405  xdim='half', nohalo=.true.)
406  call hist_in( qflx_rhoq(:,:,:,ydir,i_qs), 'SGS_YFLX_QS', 'SGS Y FLUX of QS', 'kg/m2/s', &
407  ydim='half', nohalo=.true.)
408  endif
409 
410  if ( i_qg > 0 ) then
411  call hist_in( qflx_rhoq(:,:,:,zdir,i_qg), 'SGS_ZFLX_QG', 'SGS Z FLUX of QG', 'kg/m2/s', &
412  zdim='half', nohalo=.true.)
413  call hist_in( qflx_rhoq(:,:,:,xdir,i_qg), 'SGS_XFLX_QG', 'SGS X FLUX of QG', 'kg/m2/s', &
414  xdim='half', nohalo=.true.)
415  call hist_in( qflx_rhoq(:,:,:,ydir,i_qg), 'SGS_YFLX_QG', 'SGS Y FLUX of QG', 'kg/m2/s', &
416  ydim='half', nohalo=.true.)
417  endif
418 #endif
419 
420  if ( statistics_checktotal ) then
421  call stat_total( total, momz_t_tb(:,:,:), 'MOMZ_t_TB' )
422  call stat_total( total, momx_t_tb(:,:,:), 'MOMX_t_TB' )
423  call stat_total( total, momy_t_tb(:,:,:), 'MOMY_t_TB' )
424  call stat_total( total, rhot_t_tb(:,:,:), 'RHOT_t_TB' )
425  call stat_total( total, tke(:,:,:), 'TKE' )
426  call stat_total( total, nu(:,:,:), 'Nu' )
427  call stat_total( total, ri(:,:,:), 'Ri' )
428  call stat_total( total, pr(:,:,:), 'Pr' )
429 
430  do iq = 1, qa
431  call stat_total( total, rhoq_t_tb(:,:,:,iq), trim(aq_name(iq))//'_t_TB' )
432  enddo
433  endif
434 
435  endif
436 
437  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
438  do j = js, je
439  do i = is, ie
440  do k = ks, ke
441  momz_t(k,i,j) = momz_t(k,i,j) + momz_t_tb(k,i,j)
442  momx_t(k,i,j) = momx_t(k,i,j) + momx_t_tb(k,i,j)
443  momy_t(k,i,j) = momy_t(k,i,j) + momy_t_tb(k,i,j)
444  rhot_t(k,i,j) = rhot_t(k,i,j) + rhot_t_tb(k,i,j)
445  enddo
446  enddo
447  enddo
448 
449  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(3)
450  do iq = 1, qa
451  do j = js, je
452  do i = is, ie
453  do k = ks, ke
454  rhoq_t(k,i,j,iq) = rhoq_t(k,i,j,iq) + rhoq_t_tb(k,i,j,iq)
455  enddo
456  enddo
457  enddo
458  enddo
459 
460  return
logical, public statistics_checktotal
calc&report variable totals to logfile?
real(rp), dimension(:,:,:), allocatable, target, public momz
subroutine, public atmos_phy_tb_calc_tend_momz(MOMZ_t_TB, QFLX_MOMZ, GSQRT, J13G, J23G, J33G, MAPF, IIS, IIE, JJS, JJE)
real(rp), dimension(:,:,:), allocatable, public atmos_phy_tb_momx_t
real(rp), dimension(:,:,:), allocatable, public atmos_phy_sf_sflx_qtrc
real(rp), dimension(:,:,:), allocatable, public atmos_phy_tb_rhot_t
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mw
real(rp), dimension(:,:,:), allocatable, target, public rhot
real(rp), dimension(:,:,:,:), allocatable, public gtrans_j23g
(2,3) element of Jacobian matrix * {G}^1/2
real(rp), dimension(:,:,:), allocatable, public momy_tp
subroutine, public atmos_phy_tb_calc_tend_momy(MOMY_t_TB, QFLX_MOMY, GSQRT, J13G, J23G, J33G, MAPF, IIS, IIE, JJS, JJE)
procedure(tb), pointer, public atmos_phy_tb
module ATMOSPHERIC Variables
real(rp), dimension(:,:,:,:), pointer, public qtrc_av
real(rp), dimension(:,:,:), allocatable, target, public momx
real(rp), dimension(:,:,:), allocatable, public atmos_phy_tb_tke
real(rp), public gtrans_j33g
(3,3) element of Jacobian matrix * {G}^1/2
real(rp), dimension(:,:,:), allocatable, public rhot_tp
real(rp), dimension(:,:,:), allocatable, target, public dens
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mv
module Statistics
subroutine, public atmos_phy_tb_calc_tend_phi(phi_t_TB, QFLX_phi, GSQRT, J13G, J23G, J33G, MAPF, IIS, IIE, JJS, JJE)
module ATMOSPHERIC Surface Variables
real(rp), dimension(:,:,:), pointer, public momx_av
real(rp), dimension(:,:,:,:), allocatable, public gtrans_mapf
map factor
real(rp), dimension(:,:,:), allocatable, public atmos_phy_tb_nu
subroutine, public atmos_phy_tb_calc_tend_momx(MOMX_t_TB, QFLX_MOMX, GSQRT, J13G, J23G, J33G, MAPF, IIS, IIE, JJS, JJE)
module GRIDTRANS
real(rp), dimension(:,:,:,:), allocatable, public rhoq_tp
module COMMUNICATION
Definition: scale_comm.F90:23
module TIME
Definition: scale_time.F90:15
real(dp), public time_dtsec_atmos_phy_tb
time interval of physics(turbulence ) [sec]
Definition: scale_time.F90:44
module Atmosphere / Physics Turbulence
real(rp), dimension(:,:,:), pointer, public dens_av
real(rp), dimension(:,:,:,:), allocatable, public gtrans_j13g
(1,3) element of Jacobian matrix * {G}^1/2
module ATMOSPHERE / Physics Turbulence
module ATMOSPHERE / Physics Turbulence
real(rp), dimension(:,:,:), allocatable, public atmos_phy_tb_momy_t
real(rp), dimension(:,:,:), allocatable, public momx_tp
real(rp), dimension(:,:,:), allocatable, target, public momy
real(rp), dimension(:,:,:), allocatable, public momz_tp
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_sh
real(rp), dimension(:,:,:,:), allocatable, public gtrans_gsqrt
transformation metrics from Z to Xi, {G}^1/2
real(rp), dimension(:,:,:,:), allocatable, public atmos_phy_tb_rhoq_t
module HISTORY
real(rp), dimension(:,:,:), pointer, public momz_av
real(rp), dimension(:,:,:), pointer, public rhot_av
real(rp), dimension(:,:,:), pointer, public momy_av
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mu
real(rp), dimension(:,:,:), allocatable, public atmos_phy_tb_tke_t
real(rp), dimension(:,:,:), allocatable, public atmos_phy_tb_momz_t
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc
Here is the call graph for this function:
Here is the caller graph for this function: