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

module ATMOSPHERE / Physics Turbulence More...

Functions/Subroutines

subroutine, public atmos_phy_tb_hybrid_setup (TB_TYPE, CDZ, CDX, CDY, CZ)
 
subroutine, public atmos_phy_tb_hybrid (qflx_sgs_momz, qflx_sgs_momx, qflx_sgs_momy, qflx_sgs_rhot, qflx_sgs_rhoq, tke, tke_t, Nu, Ri, Pr, N2, MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_QV, GSQRT, J13G, J23G, J33G, MAPF, dt)
 

Variables

procedure(tb), pointer sgs_tb => NULL()
 

Detailed Description

module ATMOSPHERE / Physics Turbulence

Description
Boundary layer turbulence model
Author
Team SCALE
History
  • 2014-09-18 (S.Nishizawa) [new]
NAMELIST
  • PARAM_ATMOS_PHY_TB_HYBRID
    nametypedefault valuecomment
    ATMOS_PHY_TB_HYBRID_SGS_DX real(RP) 100.0_RP horizontal resolution for SGS model
    ATMOS_PHY_TB_HYBRID_PBL_DX real(RP) 500.0_RP horizontal resolution for PBL model
    ATMOS_PHY_TB_HYBRID_SGS_TYPE character(len=H_SHORT) 'SMAGORINSKY' scheme type for SGS
    ATMOS_PHY_TB_HYBRID_PBL_TYPE character(len=H_SHORT) 'MYNN' scheme type for turbulent parametarization
    ATMOS_PHY_TB_HYBRID_TKE_TYPE character(len=H_SHORT) 'PBL' SGS, PBL, or MIXED

History Output
No history output

Function/Subroutine Documentation

◆ atmos_phy_tb_hybrid_setup()

subroutine, public scale_atmos_phy_tb_hybrid::atmos_phy_tb_hybrid_setup ( character(len=*), intent(in)  TB_TYPE,
real(rp), dimension(ka), intent(in)  CDZ,
real(rp), dimension(ia), intent(in)  CDX,
real(rp), dimension(ja), intent(in)  CDY,
real(rp), dimension (ka,ia,ja), intent(in)  CZ 
)

Definition at line 115 of file scale_atmos_phy_tb_hybrid.F90.

References scale_atmos_phy_tb_mynn::atmos_phy_tb_mynn(), scale_atmos_phy_tb_mynn::atmos_phy_tb_mynn_setup(), scale_atmos_phy_tb_smg::atmos_phy_tb_smg(), scale_atmos_phy_tb_smg::atmos_phy_tb_smg_setup(), scale_grid_index::ia, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, scale_grid_index::ja, scale_process::prc_mpistop(), and sgs_tb.

Referenced by scale_atmos_phy_tb::atmos_phy_tb_setup().

115  use scale_process, only: &
117  use scale_atmos_phy_tb_smg, only: &
120  use scale_atmos_phy_tb_mynn, only: &
123  implicit none
124 
125  character(len=*), intent(in) :: tb_type
126 
127  real(RP), intent(in) :: cdz(ka)
128  real(RP), intent(in) :: cdx(ia)
129  real(RP), intent(in) :: cdy(ja)
130  real(RP), intent(in) :: cz (ka,ia,ja)
131 
132  real(RP) :: atmos_phy_tb_hybrid_sgs_dx = 100.0_rp
133  real(RP) :: atmos_phy_tb_hybrid_pbl_dx = 500.0_rp
134  character(len=H_SHORT) :: atmos_phy_tb_hybrid_sgs_type = 'SMAGORINSKY'
135  character(len=H_SHORT) :: atmos_phy_tb_hybrid_pbl_type = 'MYNN'
136 
137  character(len=H_SHORT) :: atmos_phy_tb_hybrid_tke_type = 'PBL'
138 
139  namelist / param_atmos_phy_tb_hybrid / &
140  atmos_phy_tb_hybrid_sgs_dx, &
141  atmos_phy_tb_hybrid_pbl_dx, &
142  atmos_phy_tb_hybrid_sgs_type, &
143  atmos_phy_tb_hybrid_pbl_type, &
144  atmos_phy_tb_hybrid_tke_type
145 
146 
147  real(RP) :: dxy
148 
149  integer :: i, j
150  integer :: ierr
151  !---------------------------------------------------------------------------
152 
153  if( io_l ) write(io_fid_log,*)
154  if( io_l ) write(io_fid_log,*) '++++++ Module[TURBULENCE] / Categ[ATMOS PHYSICS] / Origin[SCALElib]'
155  if( io_l ) write(io_fid_log,*) '+++ SGS-parameterization hybrid Model'
156 
157  if ( tb_type /= 'HYBRID' ) then
158  write(*,*) 'xxx ATMOS_PHY_TB_TYPE is not HYBRID. Check!'
159  call prc_mpistop
160  endif
161 
162 
163  !--- read namelist
164  rewind(io_fid_conf)
165  read(io_fid_conf,nml=param_atmos_phy_tb_hybrid,iostat=ierr)
166  if( ierr < 0 ) then !--- missing
167  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
168  elseif( ierr > 0 ) then !--- fatal error
169  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_TB_HYBRID. Check!'
170  call prc_mpistop
171  endif
172  if( io_lnml ) write(io_fid_log,nml=param_atmos_phy_tb_hybrid)
173 
174  select case ( atmos_phy_tb_hybrid_sgs_type )
175  case ('SMAGORINSKY')
176  call atmos_phy_tb_smg_setup( &
177  atmos_phy_tb_hybrid_sgs_type, &
178  cdz, cdx, cdy, &
179  cz )
180  sgs_tb => atmos_phy_tb_smg
181  case default
182  write(*,*) 'xxx ATMOS_PHY_TB_HYBRID_SGS_TYPE is invalid'
183  call prc_mpistop
184  end select
185 
186  select case ( atmos_phy_tb_hybrid_pbl_type )
187  case ('MYNN')
189  atmos_phy_tb_hybrid_pbl_type, &
190  cdz, cdx, cdy, &
191  cz )
192  pbl_tb => atmos_phy_tb_mynn
193  case default
194  write(*,*) 'xxx ATMOS_PHY_TB_HYBRID_PBL_TYPE is invalid'
195  call prc_mpistop
196  end select
197 
198  allocate( frac_sgs(ia,ja) )
199  allocate( frac_pbl(ia,ja) )
200  allocate( frac_sgs_tke(ia,ja) )
201  allocate( frac_pbl_tke(ia,ja) )
202 
203  do j = 1, ja
204  do i = 1, ia
205  dxy = sqrt( ( cdx(i)**2 + cdy(j)**2 )*0.5_rp )
206  frac_pbl(i,j) = &
207  min( 1.0_rp, &
208  max( 0.0_rp, &
209  ( dxy - atmos_phy_tb_hybrid_sgs_dx ) &
210  / ( atmos_phy_tb_hybrid_pbl_dx - atmos_phy_tb_hybrid_sgs_dx ) ) )
211  frac_sgs(i,j) = 1.0_rp - frac_pbl(i,j)
212  end do
213  end do
214 
215  select case ( atmos_phy_tb_hybrid_tke_type )
216  case ('SGS')
217  frac_pbl_tke(:,:) = 0.0_rp
218  frac_sgs_tke(:,:) = 1.0_rp
219  case ('PBL')
220  frac_pbl_tke(:,:) = 1.0_rp
221  frac_sgs_tke(:,:) = 0.0_rp
222  case ('MIXED')
223  frac_pbl_tke(:,:) = frac_pbl(:,:)
224  frac_sgs_tke(:,:) = frac_sgs(:,:)
225  case default
226  write(*,*) 'xxx ATMOS_PHY_TB_HYBRID_TKE_TYPE is invalid'
227  call prc_mpistop
228  end select
229 
230  return
subroutine, public prc_mpistop
Abort MPI.
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
subroutine, public atmos_phy_tb_mynn(qflx_sgs_momz, qflx_sgs_momx, qflx_sgs_momy, qflx_sgs_rhot, qflx_sgs_rhoq, tke, tke_t, Nu, Ri, Pr, N2, MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_QV, GSQRT, J13G, J23G, J33G, MAPF, dt)
subroutine, public atmos_phy_tb_mynn_setup(TYPE_TB, CDZ, CDX, CDY, CZ)
module ATMOSPHERE / Physics Turbulence
integer, public ia
of x whole cells (local, with HALO)
integer, public ka
of z whole cells (local, with HALO)
subroutine, public atmos_phy_tb_smg(qflx_sgs_momz, qflx_sgs_momx, qflx_sgs_momy, qflx_sgs_rhot, qflx_sgs_rhoq, tke, tke_t, nu, Ri, Pr, N2, MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_QV, GSQRT, J13G, J23G, J33G, MAPF, dt)
subroutine, public atmos_phy_tb_smg_setup(TYPE_TB, CDZ, CDX, CDY, CZ)
module PROCESS
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
module ATMOSPHERE / Physics Turbulence
integer, public ja
of y whole cells (local, with HALO)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ atmos_phy_tb_hybrid()

subroutine, public scale_atmos_phy_tb_hybrid::atmos_phy_tb_hybrid ( real(rp), dimension(ka,ia,ja,3), intent(out)  qflx_sgs_momz,
real(rp), dimension(ka,ia,ja,3), intent(out)  qflx_sgs_momx,
real(rp), dimension(ka,ia,ja,3), intent(out)  qflx_sgs_momy,
real(rp), dimension(ka,ia,ja,3), intent(out)  qflx_sgs_rhot,
real(rp), dimension(ka,ia,ja,3,qa), intent(out)  qflx_sgs_rhoq,
real(rp), dimension (ka,ia,ja), intent(inout)  tke,
real(rp), dimension(ka,ia,ja), intent(out)  tke_t,
real(rp), dimension(ka,ia,ja), intent(out)  Nu,
real(rp), dimension(ka,ia,ja), intent(out)  Ri,
real(rp), dimension(ka,ia,ja), intent(out)  Pr,
real(rp), dimension(ka,ia,ja), intent(out)  N2,
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), intent(in)  DENS,
real(rp), dimension(ka,ia,ja,qa), intent(in)  QTRC,
real(rp), dimension(ia,ja), intent(in)  SFLX_MW,
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 (ka,ia,ja,7), intent(in)  GSQRT,
real(rp), dimension (ka,ia,ja,7), intent(in)  J13G,
real(rp), dimension (ka,ia,ja,7), intent(in)  J23G,
real(rp), intent(in)  J33G,
real(rp), dimension (ia,ja,2,4), intent(in)  MAPF,
real(dp), intent(in)  dt 
)
Parameters
[in]gsqrtvertical metrics {G}^1/2
[in]j13g(1,3) element of Jacobian matrix
[in]j23g(1,3) element of Jacobian matrix
[in]j33g(3,3) element of Jacobian matrix
[in]mapfmap factor

Definition at line 242 of file scale_atmos_phy_tb_hybrid.F90.

References scale_const::const_grav, scale_grid_index::ia, scale_grid_index::ja, scale_grid_index::ke, scale_grid_index::ks, scale_tracer::qa, sgs_tb, scale_grid_index::xdir, scale_grid_index::ydir, and scale_grid_index::zdir.

Referenced by scale_atmos_phy_tb::atmos_phy_tb_setup().

242  use scale_precision
243  use scale_grid_index
244  use scale_tracer
245  use scale_const, only: &
246  grav => const_grav
247  implicit none
248 
249  real(RP), intent(out) :: qflx_sgs_momz(ka,ia,ja,3)
250  real(RP), intent(out) :: qflx_sgs_momx(ka,ia,ja,3)
251  real(RP), intent(out) :: qflx_sgs_momy(ka,ia,ja,3)
252  real(RP), intent(out) :: qflx_sgs_rhot(ka,ia,ja,3)
253  real(RP), intent(out) :: qflx_sgs_rhoq(ka,ia,ja,3,qa)
254 
255  real(RP), intent(inout) :: tke (ka,ia,ja) ! TKE
256  real(RP), intent(out) :: tke_t(ka,ia,ja) ! tendency TKE
257  real(RP), intent(out) :: nu(ka,ia,ja) ! eddy viscosity (center)
258  real(RP), intent(out) :: pr(ka,ia,ja) ! Plandtle number
259  real(RP), intent(out) :: ri(ka,ia,ja) ! Richardson number
260  real(RP), intent(out) :: n2(ka,ia,ja) ! squared Brunt-Vaisala frequency
261 
262  real(RP), intent(in) :: momz(ka,ia,ja)
263  real(RP), intent(in) :: momx(ka,ia,ja)
264  real(RP), intent(in) :: momy(ka,ia,ja)
265  real(RP), intent(in) :: rhot(ka,ia,ja)
266  real(RP), intent(in) :: dens(ka,ia,ja)
267  real(RP), intent(in) :: qtrc(ka,ia,ja,qa)
268 
269  real(RP), intent(in) :: sflx_mw(ia,ja)
270  real(RP), intent(in) :: sflx_mu(ia,ja)
271  real(RP), intent(in) :: sflx_mv(ia,ja)
272  real(RP), intent(in) :: sflx_sh(ia,ja)
273  real(RP), intent(in) :: sflx_qv(ia,ja)
274 
275  real(RP), intent(in) :: gsqrt (ka,ia,ja,7)
276  real(RP), intent(in) :: j13g (ka,ia,ja,7)
277  real(RP), intent(in) :: j23g (ka,ia,ja,7)
278  real(RP), intent(in) :: j33g
279  real(RP), intent(in) :: mapf (ia,ja,2,4)
280  real(DP), intent(in) :: dt
281 
282  real(RP) :: w_qflx_sgs_momz(ka,ia,ja,3,2)
283  real(RP) :: w_qflx_sgs_momx(ka,ia,ja,3,2)
284  real(RP) :: w_qflx_sgs_momy(ka,ia,ja,3,2)
285  real(RP) :: w_qflx_sgs_rhot(ka,ia,ja,3,2)
286  real(RP) :: w_qflx_sgs_rhoq(ka,ia,ja,3,qa,2)
287 
288  real(RP) :: w_tke(ka,ia,ja,2)
289  real(RP) :: w_tke_t(ka,ia,ja,2)
290  real(RP) :: w_nu(ka,ia,ja,2)
291  real(RP) :: w_ri(ka,ia,ja,2)
292  real(RP) :: w_pr(ka,ia,ja,2)
293  real(RP) :: w_n2(ka,ia,ja,2)
294 
295  integer :: k, i, j, iq
296 
297  do j = 1, ja
298  do i = 1, ia
299  do k = ks, ke
300  w_tke(k,i,j,1) = tke(k,i,j)
301  w_tke(k,i,j,2) = tke(k,i,j)
302  end do
303  end do
304  end do
305 
306  call sgs_tb( &
307  w_qflx_sgs_momz(:,:,:,:,1), w_qflx_sgs_momx(:,:,:,:,1), & ! (out)
308  w_qflx_sgs_momy(:,:,:,:,1), w_qflx_sgs_rhot(:,:,:,:,1), & ! (out
309  w_qflx_sgs_rhoq(:,:,:,:,:,1), & ! (out)
310  w_tke(:,:,:,1), & ! (inout)
311  w_tke_t(:,:,:,1), & ! (out)
312  w_nu(:,:,:,1), w_ri(:,:,:,1), w_pr(:,:,:,1), w_n2(:,:,:,1), & ! (out)
313  momz, momx, momy, rhot, dens, qtrc, & ! (in)
314  sflx_mw, sflx_mu, sflx_mv, sflx_sh, sflx_qv, & ! (in)
315  gsqrt, j13g, j23g, j33g, mapf, dt ) ! (in)
316 
317  call pbl_tb( &
318  w_qflx_sgs_momz(:,:,:,:,2), w_qflx_sgs_momx(:,:,:,:,2), & ! (out)
319  w_qflx_sgs_momy(:,:,:,:,2), w_qflx_sgs_rhot(:,:,:,:,2), & ! (out
320  w_qflx_sgs_rhoq(:,:,:,:,:,2), & ! (out)
321  w_tke(:,:,:,2), & ! (inout)
322  w_tke_t(:,:,:,2), & ! (out)
323  w_nu(:,:,:,2), w_ri(:,:,:,2), w_pr(:,:,:,2), w_n2(:,:,:,2), & ! (out)
324  momz, momx, momy, rhot, dens, qtrc, & ! (in)
325  sflx_mw, sflx_mu, sflx_mv, sflx_sh, sflx_qv, & ! (in)
326  gsqrt, j13g, j23g, j33g, mapf, dt ) ! (in)
327 
328  do j = 1, ja
329  do i = 1, ia
330  do k = ks, ke
331  qflx_sgs_momz(k,i,j,zdir) = w_qflx_sgs_momz(k,i,j,zdir,1) * frac_sgs(i,j) &
332  + w_qflx_sgs_momz(k,i,j,zdir,2) * frac_pbl(i,j)
333  qflx_sgs_momz(k,i,j,xdir) = w_qflx_sgs_momz(k,i,j,xdir,1)
334  qflx_sgs_momz(k,i,j,ydir) = w_qflx_sgs_momz(k,i,j,ydir,1)
335  end do
336  end do
337  end do
338 
339  do j = 1, ja
340  do i = 1, ia
341  do k = ks, ke
342  qflx_sgs_momx(k,i,j,zdir) = w_qflx_sgs_momx(k,i,j,zdir,1) * frac_sgs(i,j) &
343  + w_qflx_sgs_momx(k,i,j,zdir,2) * frac_pbl(i,j)
344  qflx_sgs_momx(k,i,j,xdir) = w_qflx_sgs_momx(k,i,j,xdir,1)
345  qflx_sgs_momx(k,i,j,ydir) = w_qflx_sgs_momx(k,i,j,ydir,1)
346  end do
347  end do
348  end do
349 
350  do j = 1, ja
351  do i = 1, ia
352  do k = ks, ke
353  qflx_sgs_momy(k,i,j,zdir) = w_qflx_sgs_momy(k,i,j,zdir,1) * frac_sgs(i,j) &
354  + w_qflx_sgs_momy(k,i,j,zdir,2) * frac_pbl(i,j)
355  qflx_sgs_momy(k,i,j,xdir) = w_qflx_sgs_momy(k,i,j,xdir,1)
356  qflx_sgs_momy(k,i,j,ydir) = w_qflx_sgs_momy(k,i,j,ydir,1)
357  end do
358  end do
359  end do
360 
361  do j = 1, ja
362  do i = 1, ia
363  do k = ks, ke
364  qflx_sgs_rhot(k,i,j,zdir) = w_qflx_sgs_rhot(k,i,j,zdir,1) * frac_sgs(i,j) &
365  + w_qflx_sgs_rhot(k,i,j,zdir,2) * frac_pbl(i,j)
366  qflx_sgs_rhot(k,i,j,xdir) = w_qflx_sgs_rhot(k,i,j,xdir,1)
367  qflx_sgs_rhot(k,i,j,ydir) = w_qflx_sgs_rhot(k,i,j,ydir,1)
368  end do
369  end do
370  end do
371 
372  do iq = 1, qa
373  do j = 1, ja
374  do i = 1, ia
375  do k = ks, ke
376  qflx_sgs_rhoq(k,i,j,zdir,iq) = w_qflx_sgs_rhoq(k,i,j,zdir,iq,1) * frac_sgs(i,j) &
377  + w_qflx_sgs_rhoq(k,i,j,zdir,iq,2) * frac_pbl(i,j)
378  qflx_sgs_rhoq(k,i,j,xdir,iq) = w_qflx_sgs_rhoq(k,i,j,xdir,iq,1)
379  qflx_sgs_rhoq(k,i,j,ydir,iq) = w_qflx_sgs_rhoq(k,i,j,ydir,iq,1)
380  end do
381  end do
382  end do
383  end do
384 
385  do j = 1, ja
386  do i = 1, ia
387  do k = ks, ke
388  tke(k,i,j) = w_tke(k,i,j,1) * frac_sgs_tke(i,j) &
389  + w_tke(k,i,j,2) * frac_pbl_tke(i,j)
390  end do
391  end do
392  end do
393 
394  do j = 1, ja
395  do i = 1, ia
396  do k = ks, ke
397  tke_t(k,i,j) = w_tke_t(k,i,j,1) * frac_sgs_tke(i,j) &
398  + w_tke_t(k,i,j,2) * frac_pbl_tke(i,j)
399  end do
400  end do
401  end do
402 
403  do j = 1, ja
404  do i = 1, ia
405  do k = ks, ke
406  nu(k,i,j) = w_nu(k,i,j,1) * frac_sgs(i,j) &
407  + w_nu(k,i,j,2) * frac_pbl(i,j)
408  end do
409  end do
410  end do
411 
412  do j = 1, ja
413  do i = 1, ia
414  do k = ks, ke
415  ri(k,i,j) = w_ri(k,i,j,1) * frac_sgs(i,j) &
416  + w_ri(k,i,j,2) * frac_pbl(i,j)
417  end do
418  end do
419  end do
420 
421  do j = 1, ja
422  do i = 1, ia
423  do k = ks, ke
424  pr(k,i,j) = w_pr(k,i,j,1) * frac_sgs(i,j) &
425  + w_pr(k,i,j,2) * frac_pbl(i,j)
426  end do
427  end do
428  end do
429 
430  do j = 1, ja
431  do i = 1, ia
432  do k = ks, ke
433  n2(k,i,j) = w_n2(k,i,j,1) * frac_sgs(i,j) &
434  + w_n2(k,i,j,2) * frac_pbl(i,j)
435  end do
436  end do
437  end do
438 
439  return
integer, parameter, public zdir
integer, parameter, public ydir
integer, public ke
end point of inner domain: z, local
integer, parameter, public xdir
integer, public qa
module grid index
module TRACER
integer, public ia
of x whole cells (local, with HALO)
integer, public ka
of z whole cells (local, with HALO)
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:48
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
module PRECISION
integer, public ja
of y whole cells (local, with HALO)
Here is the caller graph for this function:

Variable Documentation

◆ sgs_tb

procedure(tb), pointer scale_atmos_phy_tb_hybrid::sgs_tb => NULL()

Definition at line 102 of file scale_atmos_phy_tb_hybrid.F90.

Referenced by atmos_phy_tb_hybrid(), and atmos_phy_tb_hybrid_setup().

102  procedure(tb), pointer :: sgs_tb => null()