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