SCALE-RM
scale_atmos_phy_tb_dns.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
13 !-------------------------------------------------------------------------------
15  !-----------------------------------------------------------------------------
16  !
17  !++ used modules
18  !
19  use scale_precision
20  use scale_stdio
21  use scale_prof
23  use scale_tracer
24 
25 #ifdef DEBUG
26  use scale_debug, only: &
27  check
28  use scale_const, only: &
29  undef => const_undef, &
30  iundef => const_undef2
31 #endif
32  !-----------------------------------------------------------------------------
33  implicit none
34  private
35  !-----------------------------------------------------------------------------
36  !
37  !++ Public procedure
38  !
39  public :: atmos_phy_tb_dns_config
40  public :: atmos_phy_tb_dns_setup
41  public :: atmos_phy_tb_dns
42 
43  !-----------------------------------------------------------------------------
44  !
45  !++ Public parameters & variables
46  !
47  !-----------------------------------------------------------------------------
48  !
49  !++ Private procedure
50  !
51  !-----------------------------------------------------------------------------
52  !
53  !++ Private parameters & variables
54  !
55  real(RP), private :: ATMOS_PHY_TB_DNS_NU = 1.512e-5_rp ! [m2/s] kinematic viscosity coefficient for air at 20degC
56 ! real(RP), private :: mu = 1.8E-5_RP ! [m2/s] molecular diffusive coefficient for air at 20degC
57  real(RP), private :: ATMOS_PHY_TB_DNS_MU = 1.512e-5_rp ! same as NU (needed based on hyposes. see Mellado 2010)
58 
59  !-----------------------------------------------------------------------------
60 contains
61  !-----------------------------------------------------------------------------
63  subroutine atmos_phy_tb_dns_config( &
64  TYPE_TB, &
65  I_TKE_out )
66  use scale_process, only: &
68  implicit none
69 
70  character(len=*), intent(in) :: type_tb
71  integer, intent(out) :: i_tke_out
72  !---------------------------------------------------------------------------
73 
74  if( io_l ) write(io_fid_log,*)
75  if( io_l ) write(io_fid_log,*) '++++++ Module[Turbulence Tracer] / Categ[ATMOS PHYSICS] / Origin[SCALElib]'
76  if( io_l ) write(io_fid_log,*) '*** Tracers for Deardorff (1980) 1.5th TKE Model'
77 
78  if ( type_tb /= 'DNS' ) then
79  write(*,*) 'xxx ATMOS_PHY_TB_TYPE is not DNS. Check!'
80  call prc_mpistop
81  endif
82 
83  i_tke_out = -1
84 
85  return
86  end subroutine atmos_phy_tb_dns_config
87  !-----------------------------------------------------------------------------
88  subroutine atmos_phy_tb_dns_setup( &
89  CDZ, CDX, CDY, CZ )
90  use scale_process, only: &
92  implicit none
93 
94  real(RP), intent(in) :: cdz(ka)
95  real(RP), intent(in) :: cdx(ia)
96  real(RP), intent(in) :: cdy(ja)
97  real(RP), intent(in) :: cz (ka,ia,ja)
98 
99  namelist / param_atmos_phy_tb_dns / &
100  atmos_phy_tb_dns_nu, &
101  atmos_phy_tb_dns_mu
102 
103  integer :: ierr
104  !---------------------------------------------------------------------------
105 
106  if( io_l ) write(io_fid_log,*)
107  if( io_l ) write(io_fid_log,*) '++++++ Module[Turbulence] / Categ[ATMOS PHYSICS] / Origin[SCALElib]'
108  if( io_l ) write(io_fid_log,*) '*** Eddy Viscocity Model for DNS'
109 
110  !--- read namelist
111  rewind(io_fid_conf)
112  read(io_fid_conf,nml=param_atmos_phy_tb_dns,iostat=ierr)
113  if( ierr < 0 ) then !--- missing
114  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
115  elseif( ierr > 0 ) then !--- fatal error
116  write(*,*) 'xxx Not appropriate names in namelist PARAM_ATMOS_PHY_TB_DNS. Check!'
117  call prc_mpistop
118  endif
119  if( io_nml ) write(io_fid_nml,nml=param_atmos_phy_tb_dns)
120 
121  return
122  end subroutine atmos_phy_tb_dns_setup
123 
124  !-----------------------------------------------------------------------------
125  subroutine atmos_phy_tb_dns( &
126  qflx_sgs_MOMZ, qflx_sgs_MOMX, qflx_sgs_MOMY, &
127  qflx_sgs_rhot, qflx_sgs_rhoq, &
128  RHOQ_t, nu, Ri, Pr, &
129  MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, N2, &
130  SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_Q, &
131  GSQRT, J13G, J23G, J33G, MAPF, dt )
133  use scale_tracer
134  use scale_const, only: &
135  grav => const_grav
136  use scale_grid, only: &
137  rcdz => grid_rcdz, &
138  rcdx => grid_rcdx, &
139  rcdy => grid_rcdy, &
140  rfdz => grid_rfdz, &
141  rfdx => grid_rfdx, &
142  rfdy => grid_rfdy
143  use scale_gridtrans, only: &
144  i_xy, &
145  i_uy, &
146  i_xv
147  implicit none
148 
149  ! SGS flux
150  real(RP), intent(out) :: qflx_sgs_momz(ka,ia,ja,3)
151  real(RP), intent(out) :: qflx_sgs_momx(ka,ia,ja,3)
152  real(RP), intent(out) :: qflx_sgs_momy(ka,ia,ja,3)
153  real(RP), intent(out) :: qflx_sgs_rhot(ka,ia,ja,3)
154  real(RP), intent(out) :: qflx_sgs_rhoq(ka,ia,ja,3,qa)
155 
156  real(RP), intent(inout) :: rhoq_t(ka,ia,ja,qa) ! tendency of rho * QTRC
157 
158  real(RP), intent(out) :: nu(ka,ia,ja) ! eddy viscosity (center)
159  real(RP), intent(out) :: ri(ka,ia,ja) ! Richardson number
160  real(RP), intent(out) :: pr(ka,ia,ja) ! Prantle number
161 
162  real(RP), intent(in) :: momz(ka,ia,ja)
163  real(RP), intent(in) :: momx(ka,ia,ja)
164  real(RP), intent(in) :: momy(ka,ia,ja)
165  real(RP), intent(in) :: rhot(ka,ia,ja)
166  real(RP), intent(in) :: dens(ka,ia,ja)
167  real(RP), intent(in) :: qtrc(ka,ia,ja,qa)
168  real(RP), intent(in) :: n2(ka,ia,ja)
169 
170  real(RP), intent(in) :: sflx_mw(ia,ja)
171  real(RP), intent(in) :: sflx_mu(ia,ja)
172  real(RP), intent(in) :: sflx_mv(ia,ja)
173  real(RP), intent(in) :: sflx_sh(ia,ja)
174  real(RP), intent(in) :: sflx_q (ia,ja,qa)
175 
176  real(RP), intent(in) :: gsqrt (ka,ia,ja,7)
177  real(RP), intent(in) :: j13g (ka,ia,ja,7)
178  real(RP), intent(in) :: j23g (ka,ia,ja,7)
179  real(RP), intent(in) :: j33g
180  real(RP), intent(in) :: mapf (ia,ja,2,4)
181  real(DP), intent(in) :: dt
182 
183  real(RP) :: pott(ka,ia,ja)
184 
185  integer :: iis, iie
186  integer :: jjs, jje
187 
188  integer :: k, i, j, iq
189  !---------------------------------------------------------------------------
190 
191 #ifdef DEBUG
192  qflx_sgs_momz(:,:,:,:) = undef
193  qflx_sgs_momx(:,:,:,:) = undef
194  qflx_sgs_momy(:,:,:,:) = undef
195  qflx_sgs_rhot(:,:,:,:) = undef
196  qflx_sgs_rhoq(:,:,:,:,:) = undef
197 
198  pott(:,:,:) = undef
199 #endif
200 
201  nu(:,:,:) = 0.0_rp
202  ri(:,:,:) = 0.0_rp
203  pr(:,:,:) = 1.0_rp
204 
205  ! potential temperature
206  do j = js-1, je+1
207  do i = is-1, ie+1
208  do k = ks, ke
209  pott(k,i,j) = rhot(k,i,j) / dens(k,i,j)
210  enddo
211  enddo
212  enddo
213 
214  !##### Start Upadate #####
215  do jjs = js, je, jblock
216  jje = jjs+jblock-1
217  do iis = is, ie, iblock
218  iie = iis+iblock-1
219 
220  !##### momentum equation (z) #####
221  ! (cell center)
222  do j = jjs, jje
223  do i = iis, iie
224  do k = ks+1, ke-1
225  qflx_sgs_momz(k,i,j,zdir) = -atmos_phy_tb_dns_nu * ( momz(k,i,j)-momz(k-1,i,j) ) * rcdz(k)
226  enddo
227  enddo
228  enddo
229 
230  do j = jjs, jje
231  do i = iis, iie
232  qflx_sgs_momz(ks,i,j,zdir) = 0.0_rp ! bottom boundary
233  qflx_sgs_momz(ke,i,j,zdir) = 0.0_rp ! top boundary
234  enddo
235  enddo
236 
237  ! (y edge)
238  do j = jjs, jje
239  do i = iis-1, iie
240  do k = ks, ke-1
241  qflx_sgs_momz(k,i,j,xdir) = -atmos_phy_tb_dns_nu * ( momz(k,i+1,j)-momz(k,i,j) ) * rfdx(i) * mapf(i,j,1,i_xy)
242  enddo
243  enddo
244  enddo
245 
246  ! (x edge)
247  do j = jjs-1, jje
248  do i = iis, iie
249  do k = ks, ke-1
250  qflx_sgs_momz(k,i,j,ydir) = -atmos_phy_tb_dns_nu * ( momz(k,i,j+1)-momz(k,i,j) ) * rfdy(j) * mapf(i,j,2,i_xy)
251  enddo
252  enddo
253  enddo
254 
255  !##### momentum equation (x) #####
256  ! (y edge)
257  do j = jjs, jje
258  do i = iis, iie
259  do k = ks, ke-1
260  qflx_sgs_momx(k,i,j,zdir) = -atmos_phy_tb_dns_nu * ( momx(k+1,i,j)-momx(k,i,j) ) * rfdz(k)
261  enddo
262  enddo
263  enddo
264 
265  do j = jjs, jje
266  do i = iis, iie
267  qflx_sgs_momx(ks-1,i,j,zdir) = 0.0_rp ! bottom boundary
268  qflx_sgs_momx(ke ,i,j,zdir) = 0.0_rp ! top boundary
269  enddo
270  enddo
271 
272  ! (cell center)
273  do j = jjs, jje
274  do i = iis, iie+1
275  do k = ks, ke
276  qflx_sgs_momx(k,i,j,xdir) = -atmos_phy_tb_dns_nu * ( momx(k,i,j)-momx(k,i-1,j) ) * rcdx(i) * mapf(i,j,1,i_uy)
277  enddo
278  enddo
279  enddo
280 
281  ! (z edge)
282  do j = jjs-1, jje
283  do i = iis, iie
284  do k = ks, ke
285  qflx_sgs_momx(k,i,j,ydir) = -atmos_phy_tb_dns_nu * ( momx(k,i,j+1)-momx(k,i,j) ) * rfdy(j) * mapf(i,j,2,i_uy)
286  enddo
287  enddo
288  enddo
289 
290  !##### momentum equation (y) #####
291 
292  ! (x edge)
293  do j = jjs, jje
294  do i = iis, iie
295  do k = ks, ke-1
296  qflx_sgs_momy(k,i,j,zdir) = -atmos_phy_tb_dns_nu * ( momy(k+1,i,j)-momy(k,i,j) ) * rfdz(k)
297  enddo
298  enddo
299  enddo
300 
301  do j = jjs, jje
302  do i = iis, iie
303  qflx_sgs_momy(ks-1,i,j,zdir) = 0.0_rp ! bottom boundary
304  qflx_sgs_momy(ke ,i,j,zdir) = 0.0_rp ! top boundary
305  enddo
306  enddo
307 
308  ! (z edge)
309  do j = jjs, jje
310  do i = iis-1, iie
311  do k = ks, ke
312  qflx_sgs_momy(k,i,j,xdir) = -atmos_phy_tb_dns_nu * ( momy(k,i+1,j)-momy(k,i,j) ) * rfdx(i) * mapf(i,j,1,i_xv)
313  enddo
314  enddo
315  enddo
316 
317  ! (z-x plane)
318  do j = jjs, jje+1
319  do i = iis, iie
320  do k = ks, ke
321  qflx_sgs_momy(k,i,j,ydir) = -atmos_phy_tb_dns_nu * ( momy(k,i,j)-momy(k,i,j-1) ) * rcdy(j) * mapf(i,j,2,i_xv)
322  enddo
323  enddo
324  enddo
325 
326  !##### Thermodynamic Equation #####
327 
328  ! at x, y ,w
329  do j = jjs, jje
330  do i = iis, iie
331  do k = ks, ke-1
332  qflx_sgs_rhot(k,i,j,zdir) = -0.5_rp * ( dens(k+1,i,j)+dens(k,i,j) ) &
333  * atmos_phy_tb_dns_mu * ( pott(k+1,i,j)-pott(k,i,j) ) * rfdz(k)
334  enddo
335  enddo
336  enddo
337 
338  do j = jjs, jje
339  do i = iis, iie
340  qflx_sgs_rhot(ks-1,i,j,zdir) = 0.0_rp
341  qflx_sgs_rhot(ke ,i,j,zdir) = 0.0_rp
342  enddo
343  enddo
344 
345  ! at u, y, z
346  do j = jjs, jje
347  do i = iis-1, iie
348  do k = ks, ke
349  qflx_sgs_rhot(k,i,j,xdir) = -0.5_rp * ( dens(k,i+1,j)+dens(k,i,j) ) &
350  * atmos_phy_tb_dns_mu * ( pott(k,i+1,j)-pott(k,i,j) ) * rfdx(i) * mapf(i,j,1,i_xy)
351  enddo
352  enddo
353  enddo
354 
355  ! at x, v, z
356  do j = jjs-1, jje
357  do i = iis, iie
358  do k = ks, ke
359  qflx_sgs_rhot(k,i,j,ydir) = -0.5_rp * ( dens(k,i,j+1)+dens(k,i,j) ) &
360  * atmos_phy_tb_dns_mu * ( pott(k,i,j+1)-pott(k,i,j) ) * rfdy(j) * mapf(i,j,2,i_xy)
361  enddo
362  enddo
363  enddo
364 
365  enddo
366  enddo
367 
368  !##### Tracers #####
369  do iq = 1, qa
370 
371  if ( .not. tracer_advc(iq) ) cycle
372 
373  do jjs = js, je, jblock
374  jje = jjs+jblock-1
375  do iis = is, ie, iblock
376  iie = iis+iblock-1
377 
378  ! at x, y ,w
379  do j = jjs, jje
380  do i = iis, iie
381  do k = ks, ke-1
382  qflx_sgs_rhoq(k,i,j,zdir,iq) = -0.5_rp * ( dens(k+1,i,j)+dens(k,i,j) ) &
383  * atmos_phy_tb_dns_mu * ( qtrc(k+1,i,j,iq)-qtrc(k,i,j,iq) ) * rfdz(k)
384  enddo
385  enddo
386  enddo
387  do j = jjs, jje
388  do i = iis, iie
389  qflx_sgs_rhoq(ks-1,i,j,zdir,iq) = 0.0_rp
390  qflx_sgs_rhoq(ke ,i,j,zdir,iq) = 0.0_rp
391  enddo
392  enddo
393 
394  ! at u, y, z
395  do j = jjs, jje
396  do i = iis-1, iie
397  do k = ks, ke
398  qflx_sgs_rhoq(k,i,j,xdir,iq) = -0.5_rp * ( dens(k,i+1,j)+dens(k,i,j) ) &
399  * atmos_phy_tb_dns_mu * ( qtrc(k,i+1,j,iq)-qtrc(k,i,j,iq) ) * rfdx(i) * mapf(i,j,1,i_xy)
400  enddo
401  enddo
402  enddo
403 
404  ! at x, v, z
405  do j = jjs-1, jje
406  do i = iis, iie
407  do k = ks, ke
408  qflx_sgs_rhoq(k,i,j,ydir,iq) = -0.5_rp * ( dens(k,i,j+1)+dens(k,i,j) ) &
409  * atmos_phy_tb_dns_mu * ( qtrc(k,i,j+1,iq)-qtrc(k,i,j,iq) ) * rfdy(j) * mapf(i,j,2,i_xy)
410  enddo
411  enddo
412  enddo
413 
414  enddo
415  enddo
416 
417  enddo ! scalar quantities loop
418 
419  return
420  end subroutine atmos_phy_tb_dns
421 
422 end module scale_atmos_phy_tb_dns
integer, public is
start point of inner domain: x, local
module DEBUG
Definition: scale_debug.F90:13
integer, public je
end point of inner domain: y, local
real(rp), dimension(:), allocatable, public grid_rcdy
reciprocal of center-dy
subroutine, public prc_mpistop
Abort MPI.
subroutine, public atmos_phy_tb_dns(qflx_sgs_MOMZ, qflx_sgs_MOMX, qflx_sgs_MOMY, qflx_sgs_rhot, qflx_sgs_rhoq, RHOQ_t, nu, Ri, Pr, MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, N2, SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_Q, GSQRT, J13G, J23G, J33G, MAPF, dt)
integer, public iblock
block size for cache blocking: x
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:61
real(rp), dimension(:), allocatable, public grid_rcdx
reciprocal of center-dx
integer, parameter, public zdir
module STDIO
Definition: scale_stdio.F90:12
integer, parameter, public ydir
integer, public ke
end point of inner domain: z, local
integer, parameter, public xdir
integer, public qa
logical, dimension(qa_max), public tracer_advc
integer, public i_xy
real(rp), dimension(:), allocatable, public grid_rfdy
reciprocal of face-dy
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:58
real(rp), dimension(:), allocatable, public grid_rcdz
reciprocal of center-dz
real(rp), public const_undef
Definition: scale_const.F90:43
module grid index
logical, public io_nml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:62
module TRACER
subroutine, public atmos_phy_tb_dns_setup(CDZ, CDX, CDY, CZ)
integer, public ia
of whole cells: x, local, with HALO
module GRIDTRANS
integer, public ka
of whole cells: z, local, with HALO
integer, public i_uy
integer, public jblock
block size for cache blocking: y
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:48
subroutine, public atmos_phy_tb_dns_config(TYPE_TB, I_TKE_out)
Config.
integer, public js
start point of inner domain: y, local
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:40
module ATMOSPHERE / Physics Turbulence
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
module GRID (cartesian)
integer, public i_xv
module profiler
Definition: scale_prof.F90:10
integer, public ie
end point of inner domain: x, local
module PRECISION
real(rp), dimension(:), allocatable, public grid_rfdx
reciprocal of face-dx
real(rp), dimension(:), allocatable, public grid_rfdz
reciprocal of face-dz
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
integer, public io_fid_nml
Log file ID (only for output namelist)
Definition: scale_stdio.F90:57
integer, public ja
of whole cells: y, local, with HALO