SCALE-RM
scale_atmos_phy_tb_d1980.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
16 !-------------------------------------------------------------------------------
17 #include "scalelib.h"
19  !-----------------------------------------------------------------------------
20  !
21  !++ used modules
22  !
23  use scale_precision
24  use scale_io
25  use scale_prof
27  use scale_tracer
28 
29 #ifdef DEBUG
30  use scale_debug, only: &
31  check
32  use scale_const, only: &
33  undef => const_undef, &
34  iundef => const_undef2
35 #endif
36  !-----------------------------------------------------------------------------
37  implicit none
38  private
39  !-----------------------------------------------------------------------------
40  !
41  !++ Public procedure
42  !
44  public :: atmos_phy_tb_d1980_setup
46  public :: atmos_phy_tb_d1980
47 
48  !-----------------------------------------------------------------------------
49  !
50  !++ Public parameters & variables
51  !
52  !-----------------------------------------------------------------------------
53  !
54  !++ Private procedure
55  !
56  !-----------------------------------------------------------------------------
57  !
58  !++ Private parameters & variables
59  !
60  real(RP), private, parameter :: C1 = 0.19_rp
61  real(RP), private, parameter :: C2 = 0.51_rp
62 
63  real(RP), private, parameter :: OneOverThree = 1.0_rp / 3.0_rp
64  real(RP), private, parameter :: TwoOverThree = 2.0_rp / 3.0_rp
65  real(RP), private, parameter :: FourOverThree = 4.0_rp / 3.0_rp
66 
67  real(RP), private, allocatable :: delta(:,:,:) ! (dx*dy*dz)^(1/3)
68 
69  integer, private :: I_TKE
70 
71  real(RP), private :: ATMOS_PHY_TB_D1980_Km_MAX = 1000.0_rp
72  logical, private :: ATMOS_PHY_TB_D1980_implicit = .false.
73 
74  !-----------------------------------------------------------------------------
75 contains
76  !-----------------------------------------------------------------------------
78  subroutine atmos_phy_tb_d1980_config( &
79  TYPE_TB, &
80  I_TKE_out )
81  use scale_prc, only: &
82  prc_abort
83  use scale_tracer, only: &
85  implicit none
86 
87  character(len=*), intent(in) :: type_tb
88  integer, intent(out) :: i_tke_out
89  !---------------------------------------------------------------------------
90 
91  log_newline
92  log_info("ATMOS_PHY_TB_d1980_config",*) 'Setup'
93  log_info("ATMOS_PHY_TB_d1980_config",*) 'Tracers for Deardorff (1980) 1.5th TKE Model'
94 
95  if ( type_tb /= 'D1980' ) then
96  log_error("ATMOS_PHY_TB_d1980_config",*) 'ATMOS_PHY_TB_TYPE is not D1980. Check!'
97  call prc_abort
98  endif
99 
100  call tracer_regist( i_tke, & ! [OUT]
101  1, & ! [IN]
102  (/ 'TKE_D1980' /), & ! [IN]
103  (/ 'turbulent kinetic energy (D1980)' /), & ! [IN]
104  (/ 'm2/s2' /) ) ! [IN]
105 
106  i_tke_out = i_tke
107 
108  return
109  end subroutine atmos_phy_tb_d1980_config
110 
111  !-----------------------------------------------------------------------------
113  subroutine atmos_phy_tb_d1980_setup( &
114  CDZ, CDX, CDY, CZ )
115  use scale_prc, only: &
116  prc_abort
117  implicit none
118 
119  real(rp), intent(in) :: cdz(ka)
120  real(rp), intent(in) :: cdx(ia)
121  real(rp), intent(in) :: cdy(ja)
122  real(rp), intent(in) :: cz (ka,ia,ja)
123 
124  namelist / param_atmos_phy_tb_d1980 / &
125  atmos_phy_tb_d1980_km_max, &
126  atmos_phy_tb_d1980_implicit
127 
128  integer :: ierr
129  integer :: k, i, j
130  !---------------------------------------------------------------------------
131 
132  log_newline
133  log_info("ATMOS_PHY_TB_d1980_setup",*) 'Setup'
134  log_info("ATMOS_PHY_TB_d1980_setup",*) 'Deardorff (1980) 1.5th TKE Model'
135 
136  !--- read namelist
137  rewind(io_fid_conf)
138  read(io_fid_conf,nml=param_atmos_phy_tb_d1980,iostat=ierr)
139  if( ierr < 0 ) then !--- missing
140  log_info("ATMOS_PHY_TB_d1980_setup",*) 'Not found namelist. Default used.'
141  elseif( ierr > 0 ) then !--- fatal error
142  log_error("ATMOS_PHY_TB_d1980_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_TB_D1980. Check!'
143  call prc_abort
144  endif
145  log_nml(param_atmos_phy_tb_d1980)
146 
147  allocate( delta(ka,ia,ja) )
148 
149 #ifdef DEBUG
150  delta(:,:,:) = undef
151 #endif
152  !$omp parallel do
153  do j = js-1, je+1
154  do i = is-1, ie+1
155  do k = ks, ke
156  delta(k,i,j) = ( cdz(k) * cdx(i) * cdy(j) )**(1.0_rp/3.0_rp)
157  enddo
158  enddo
159  enddo
160 
161  return
162  end subroutine atmos_phy_tb_d1980_setup
163 
164  !-----------------------------------------------------------------------------
166  subroutine atmos_phy_tb_d1980_finalize
167 
168  deallocate( delta )
169 
170  return
171  end subroutine atmos_phy_tb_d1980_finalize
172 
173  !-----------------------------------------------------------------------------
174  subroutine atmos_phy_tb_d1980( &
175  qflx_sgs_momz, qflx_sgs_momx, qflx_sgs_momy, &
176  qflx_sgs_rhot, qflx_sgs_rhoq, &
177  RHOQ_t, &
178  Km, Ri, Pr, &
179  MOMZ, MOMX, MOMY, RHOT, DENS, QTRC, N2, &
180  SFLX_MW, SFLX_MU, SFLX_MV, SFLX_SH, SFLX_Q, &
181  GSQRT, J13G, J23G, J33G, MAPF, dt )
184  use scale_tracer
185  use scale_const, only: &
186  grav => const_grav
187  use scale_atmos_grid_cartesc, only: &
188  fdz => atmos_grid_cartesc_fdz, &
189  rcdz => atmos_grid_cartesc_rcdz, &
190  rcdx => atmos_grid_cartesc_rcdx, &
191  rcdy => atmos_grid_cartesc_rcdy, &
192  rfdz => atmos_grid_cartesc_rfdz, &
193  rfdx => atmos_grid_cartesc_rfdx, &
195  use scale_matrix, only: &
196  matrix_solver_tridiagonal
197  use scale_atmos_phy_tb_common, only: &
198  calc_strain_tensor => atmos_phy_tb_calc_strain_tensor, &
199  calc_tend_momz => atmos_phy_tb_calc_tend_momz, &
200  calc_tend_momx => atmos_phy_tb_calc_tend_momx, &
201  calc_tend_momy => atmos_phy_tb_calc_tend_momy, &
202  calc_tend_phi => atmos_phy_tb_calc_tend_phi, &
203  calc_flux_phi => atmos_phy_tb_calc_flux_phi
204  implicit none
205 
206  ! SGS flux
207  real(rp), intent(out) :: qflx_sgs_momz(ka,ia,ja,3)
208  real(rp), intent(out) :: qflx_sgs_momx(ka,ia,ja,3)
209  real(rp), intent(out) :: qflx_sgs_momy(ka,ia,ja,3)
210  real(rp), intent(out) :: qflx_sgs_rhot(ka,ia,ja,3)
211  real(rp), intent(out) :: qflx_sgs_rhoq(ka,ia,ja,3,qa)
212 
213  real(rp), intent(inout) :: rhoq_t (ka,ia,ja,qa) ! tendency of rho * QTRC
214 
215  real(rp), intent(out) :: km (ka,ia,ja) ! eddy viscosity (center)
216  real(rp), intent(out) :: ri (ka,ia,ja) ! Richardson number
217  real(rp), intent(out) :: pr (ka,ia,ja) ! Prantle number
218 
219  real(rp), intent(in) :: momz (ka,ia,ja)
220  real(rp), intent(in) :: momx (ka,ia,ja)
221  real(rp), intent(in) :: momy (ka,ia,ja)
222  real(rp), intent(in) :: rhot (ka,ia,ja)
223  real(rp), intent(in) :: dens (ka,ia,ja)
224  real(rp), intent(in) :: qtrc (ka,ia,ja,qa)
225  real(rp), intent(in) :: n2 (ka,ia,ja)
226 
227  real(rp), intent(in) :: sflx_mw (ia,ja)
228  real(rp), intent(in) :: sflx_mu (ia,ja)
229  real(rp), intent(in) :: sflx_mv (ia,ja)
230  real(rp), intent(in) :: sflx_sh (ia,ja)
231  real(rp), intent(in) :: sflx_q (ia,ja,qa)
232 
233  real(rp), intent(in) :: gsqrt (ka,ia,ja,7)
234  real(rp), intent(in) :: j13g (ka,ia,ja,7)
235  real(rp), intent(in) :: j23g (ka,ia,ja,7)
236  real(rp), intent(in) :: j33g
237  real(rp), intent(in) :: mapf (ia,ja,2,4)
238  real(dp), intent(in) :: dt
239 
240  ! diagnostic variables
241  real(rp) :: pott (ka,ia,ja)
242 
243  ! deformation rate tensor
244  real(rp) :: s33_c (ka,ia,ja) ! (cell center)
245  real(rp) :: s11_c (ka,ia,ja)
246  real(rp) :: s22_c (ka,ia,ja)
247  real(rp) :: s31_c (ka,ia,ja)
248  real(rp) :: s12_c (ka,ia,ja)
249  real(rp) :: s23_c (ka,ia,ja)
250  real(rp) :: s12_z (ka,ia,ja) ! (z edge or x-y plane)
251  real(rp) :: s23_x (ka,ia,ja) ! (x edge or y-z plane)
252  real(rp) :: s31_y (ka,ia,ja) ! (y edge or z-x plane)
253  real(rp) :: s2 (ka,ia,ja) ! |S|^2
254 
255  real(rp) :: kh (ka,ia,ja) ! eddy diffusion
256  real(rp) :: l (ka,ia,ja) ! mixing length
257 
258  real(rp) :: qflx_tke(ka,ia,ja,3)
259 
260  real(rp) :: tend(ka,ia,ja)
261  real(rp) :: tend2(ka,ia,ja)
262  real(rp) :: a(ka,ia,ja)
263  real(rp) :: b(ka,ia,ja)
264  real(rp) :: c(ka,ia,ja)
265  real(rp) :: ap
266 
267  integer :: iis, iie
268  integer :: jjs, jje
269 
270  integer :: k, i, j, iq
271  !---------------------------------------------------------------------------
272 
273  log_progress(*) 'atmosphere / physics / turbulence / D1980'
274 
275 #ifdef DEBUG
276  qflx_sgs_momz(:,:,:,:) = undef
277  qflx_sgs_momx(:,:,:,:) = undef
278  qflx_sgs_momy(:,:,:,:) = undef
279  qflx_sgs_rhot(:,:,:,:) = undef
280  qflx_sgs_rhoq(:,:,:,:,:) = undef
281 
282  pr(:,:,:) = undef
283  ri(:,:,:) = undef
284  kh(:,:,:) = undef
285  km(:,:,:) = undef
286 
287  pott(:,:,:) = undef
288 #endif
289 
290  ! potential temperature
291  !$omp parallel do
292  do j = js-1, je+1
293  do i = is-1, ie+1
294  do k = ks, ke
295 #ifdef DEBUG
296  call check( __line__, rhot(k,i,j) )
297  call check( __line__, dens(k,i,j) )
298 #endif
299  pott(k,i,j) = rhot(k,i,j) / dens(k,i,j)
300  enddo
301  enddo
302  enddo
303 #ifdef DEBUG
304  i = iundef; j = iundef; k = iundef
305 #endif
306 
307  !##### Start Upadate #####
308 
309  call calc_strain_tensor( &
310  s33_c, s11_c, s22_c, & ! (out)
311  s31_c, s12_c, s23_c, & ! (out)
312  s12_z, s23_x, s31_y, & ! (out)
313  s2 , & ! (out)
314  dens, momz, momx, momy, & ! (in)
315  gsqrt, j13g, j23g, j33g, mapf ) ! (in)
316 
317  do jjs = js, je, jblock
318  jje = jjs+jblock-1
319  do iis = is, ie, iblock
320  iie = iis+iblock-1
321 
322  ! Ri = N^2 / |S|^2, N^2 = g / theta * dtheta/dz
323  !$omp parallel do
324  do j = jjs-1, jje+1
325  do i = iis-1, iie+1
326  do k = ks, ke
327 #ifdef DEBUG
328  call check( __line__, pott(k+1,i,j) )
329  call check( __line__, pott(k,i,j) )
330  call check( __line__, pott(k-1,i,j) )
331  call check( __line__, fdz(k) )
332  call check( __line__, fdz(k-1) )
333  call check( __line__, s2(k,i,j) )
334 #endif
335  ri(k,i,j) = n2(k,i,j) / s2(k,i,j)
336  enddo
337  enddo
338  enddo
339 #ifdef DEBUG
340  i = iundef; j = iundef; k = iundef
341 #endif
342 
343  ! mixing length
344  !$omp parallel do
345  do j = jjs-1, jje+1
346  do i = iis-1, iie+1
347  do k = ks, ke
348  if ( n2(k,i,j) > 1e-10_rp ) then
349  l(k,i,j) = max( min( 0.76_rp * sqrt(qtrc(k,i,j,i_tke)/n2(k,i,j)), delta(k,i,j) ), 1e-10_rp )
350  else
351  l(k,i,j) = delta(k,i,j)
352  end if
353  enddo
354  enddo
355  enddo
356 #ifdef DEBUG
357  i = iundef; j = iundef; k = iundef
358 #endif
359 
360  !$omp parallel do
361  do j = jjs-1, jje+1
362  do i = iis-1, iie+1
363  do k = ks, ke
364  km(k,i,j) = min( 0.10_rp * l(k,i,j) * sqrt(qtrc(k,i,j,i_tke)), atmos_phy_tb_d1980_km_max )
365  pr(k,i,j) = 1.0_rp / ( 1.0_rp + 2.0_rp * l(k,i,j) / delta(k,i,j) )
366  kh(k,i,j) = km(k,i,j) / pr(k,i,j)
367  enddo
368  enddo
369  enddo
370 #ifdef DEBUG
371  i = iundef; j = iundef; k = iundef
372 #endif
373 
374  !##### momentum equation (z) #####
375  ! (cell center)
376  !$omp parallel do
377  do j = jjs, jje
378  do i = iis, iie
379  do k = ks+1, ke-1
380 #ifdef DEBUG
381  call check( __line__, dens(k,i,j) )
382  call check( __line__, km(k,i,j) )
383  call check( __line__, s33_c(k,i,j) )
384 #endif
385  qflx_sgs_momz(k,i,j,zdir) = dens(k,i,j) * ( - 2.0_rp * km(k,i,j) * s33_c(k,i,j) )
386  enddo
387  enddo
388  enddo
389 #ifdef DEBUG
390  i = iundef; j = iundef; k = iundef
391 #endif
392  !$omp parallel do
393  do j = jjs, jje
394  do i = iis, iie
395  qflx_sgs_momz(ks,i,j,zdir) = 0.0_rp ! just above bottom boundary
396  qflx_sgs_momz(ke,i,j,zdir) = 0.0_rp ! just below top boundary
397  enddo
398  enddo
399 #ifdef DEBUG
400  i = iundef; j = iundef; k = iundef
401 #endif
402  ! (y edge)
403  !$omp parallel do
404  do j = jjs, jje
405  do i = iis-1, iie
406  do k = ks, ke-1
407 #ifdef DEBUG
408  call check( __line__, dens(k,i,j) )
409  call check( __line__, dens(k,i+1,j) )
410  call check( __line__, dens(k+1,i,j) )
411  call check( __line__, dens(k+1,i+1,j) )
412  call check( __line__, km(k,i,j) )
413  call check( __line__, km(k,i+1,j) )
414  call check( __line__, km(k+1,i,j) )
415  call check( __line__, km(k+1,i+1,j) )
416  call check( __line__, s31_y(k,i,j) )
417 #endif
418  qflx_sgs_momz(k,i,j,xdir) = - 0.125_rp & ! 2.0 / 4 / 4
419  * ( dens(k,i,j)+dens(k+1,i,j)+dens(k,i+1,j)+dens(k+1,i+1,j) ) &
420  * ( km(k,i,j)+km(k+1,i,j)+km(k,i+1,j)+km(k+1,i+1,j)) &
421  * s31_y(k,i,j)
422  enddo
423  enddo
424  enddo
425 #ifdef DEBUG
426  i = iundef; j = iundef; k = iundef
427 #endif
428  ! (x edge)
429  !$omp parallel do
430  do j = jjs-1, jje
431  do i = iis, iie
432  do k = ks, ke-1
433 #ifdef DEBUG
434  call check( __line__, dens(k,i,j) )
435  call check( __line__, dens(k,i,j+1) )
436  call check( __line__, dens(k+1,i,j) )
437  call check( __line__, dens(k+1,i,j+1) )
438  call check( __line__, km(k,i,j) )
439  call check( __line__, km(k,i,j+1) )
440  call check( __line__, km(k+1,i,j) )
441  call check( __line__, km(k+1,i,j+1) )
442  call check( __line__, s23_x(k,i,j) )
443 #endif
444  qflx_sgs_momz(k,i,j,ydir) = - 0.125_rp & ! 2/4/4
445  * ( dens(k,i,j)+dens(k+1,i,j)+dens(k,i,j+1)+dens(k+1,i,j+1) ) &
446  * ( km(k,i,j)+km(k+1,i,j)+km(k,i,j+1)+km(k+1,i,j+1) ) &
447  * s23_x(k,i,j)
448  enddo
449  enddo
450  enddo
451 #ifdef DEBUG
452  i = iundef; j = iundef; k = iundef
453 #endif
454 
455  if ( atmos_phy_tb_d1980_implicit ) then
456  call calc_tend_momz( tend, & ! (out)
457  qflx_sgs_momz, & ! (in)
458  gsqrt, j13g, j23g, j33g, mapf, & ! (in)
459  iis, iie, jjs, jje ) ! (in)
460 
461  !$omp parallel do &
462  !$omp private(ap)
463  do j = jjs, jje
464  do i = iis, iie
465 
466  ap = - fouroverthree * dt &
467  * dens(ks+1,i,j)*km(ks+1,i,j) &
468  * rcdz(ks+1) / gsqrt(ks+1,i,j,i_xyz)
469  a(ks,i,j) = ap * rfdz(ks) / gsqrt(ks,i,j,i_xyw)
470  c(ks,i,j) = 0.0_rp
471  b(ks,i,j) = - a(ks,i,j) + 0.5_rp * ( dens(ks,i,j)+dens(ks+1,i,j) )
472  do k = ks+1, ke-2
473  c(k,i,j) = ap * rfdz(k+1) / gsqrt(k+1,i,j,i_xyw)
474  ap = - fouroverthree * dt &
475  * dens(k+1,i,j)*km(k+1,i,j) &
476  * rcdz(k+1) / gsqrt(k+1,i,j,i_xyz)
477  a(k,i,j) = ap * rfdz(k) / gsqrt(k,i,j,i_xyw)
478  b(k,i,j) = - a(k,i,j) - c(k,i,j) + 0.5_rp * ( dens(k,i,j)+dens(k+1,i,j) )
479  end do
480  a(ke-1,i,j) = 0.0_rp
481  c(ke-1,i,j) = ap * rfdz(ke) / gsqrt(ke,i,j,i_xyw)
482  b(ke-1,i,j) = - c(ke-1,i,j) + 0.5_rp * ( dens(ke-1,i,j)+dens(ke,i,j) )
483  end do
484  end do
485 
486  call matrix_solver_tridiagonal( ka, ks, ke-1, ia, iis, iie, ja, jjs, jje, &
487  a(:,:,:), b(:,:,:), c(:,:,:), tend(:,:,:), & ! (in)
488  tend2(:,:,:) ) ! (out)
489 
490  !$omp parallel do
491  do j = jjs, jje
492  do i = iis, iie
493  do k = ks+1, ke-1
494  qflx_sgs_momz(k,i,j,zdir) = qflx_sgs_momz(k,i,j,zdir) &
495  - fouroverthree * dens(k,i,j) * km(k,i,j) * dt &
496  * ( tend2(k,i,j) - tend2(k-1,i,j) ) * rcdz(k) / gsqrt(k,i,j,i_xyz)
497  end do
498 
499  end do
500  end do
501 
502  end if
503 
504  !##### momentum equation (x) #####
505  ! (y edge)
506  !$omp parallel do
507  do j = jjs, jje
508  do i = iis, iie
509  do k = ks, ke-1
510 #ifdef DEBUG
511  call check( __line__, dens(k,i,j) )
512  call check( __line__, dens(k,i+1,j) )
513  call check( __line__, dens(k+1,i,j) )
514  call check( __line__, dens(k+1,i+1,j) )
515  call check( __line__, km(k,i,j) )
516  call check( __line__, km(k,i+1,j) )
517  call check( __line__, km(k+1,i,j) )
518  call check( __line__, km(k+1,i+1,j) )
519  call check( __line__, s31_y(k,i,j) )
520 #endif
521  qflx_sgs_momx(k,i,j,zdir) = - 0.125_rp & ! 2/4/4
522  * ( dens(k,i,j)+dens(k+1,i,j)+dens(k,i+1,j)+dens(k+1,i+1,j) ) &
523  * ( km(k,i,j)+km(k+1,i,j)+km(k,i+1,j)+km(k+1,i+1,j) ) &
524  * s31_y(k,i,j)
525  enddo
526  enddo
527  enddo
528 #ifdef DEBUG
529  i = iundef; j = iundef; k = iundef
530 #endif
531  !$omp parallel do
532  do j = jjs, jje
533  do i = iis, iie
534  qflx_sgs_momx(ks-1,i,j,zdir) = 0.0_rp ! bottom boundary
535  qflx_sgs_momx(ke ,i,j,zdir) = 0.0_rp ! top boundary
536  enddo
537  enddo
538 #ifdef DEBUG
539  i = iundef; j = iundef; k = iundef
540 #endif
541  ! (cell center)
542  !$omp parallel do
543  do j = jjs, jje
544  do i = iis, iie+1
545  do k = ks, ke
546 #ifdef DEBUG
547  call check( __line__, dens(k,i,j) )
548  call check( __line__, km(k,i,j) )
549  call check( __line__, s11_c(k,i,j) )
550 #endif
551  qflx_sgs_momx(k,i,j,xdir) = dens(k,i,j) * ( - 2.0_rp * km(k,i,j) * s11_c(k,i,j) )
552  enddo
553  enddo
554  enddo
555 #ifdef DEBUG
556  i = iundef; j = iundef; k = iundef
557 #endif
558  ! (z edge)
559  !$omp parallel do
560  do j = jjs-1, jje
561  do i = iis, iie
562  do k = ks, ke
563 #ifdef DEBUG
564  call check( __line__, dens(k,i,j) )
565  call check( __line__, dens(k,i+1,j) )
566  call check( __line__, dens(k,i,j+1) )
567  call check( __line__, dens(k,i+1,j+1) )
568  call check( __line__, km(k,i,j) )
569  call check( __line__, km(k,i+1,j) )
570  call check( __line__, km(k,i,j+1) )
571  call check( __line__, km(k,i+1,j+1) )
572  call check( __line__, s12_z(k,i,j) )
573 #endif
574  qflx_sgs_momx(k,i,j,ydir) = - 0.125_rp & ! 2/4/4
575  * ( dens(k,i,j)+dens(k,i+1,j)+dens(k,i,j+1)+dens(k,i+1,j+1) ) &
576  * ( km(k,i,j)+km(k,i+1,j)+km(k,i,j+1)+km(k,i+1,j+1) ) &
577  * s12_z(k,i,j)
578  enddo
579  enddo
580  enddo
581 #ifdef DEBUG
582  i = iundef; j = iundef; k = iundef
583 #endif
584 
585  if ( atmos_phy_tb_d1980_implicit ) then
586  call calc_tend_momx( tend, & ! (out)
587  qflx_sgs_momx, & ! (in)
588  gsqrt, j13g, j23g, j33g, mapf, & ! (in)
589  iis, iie, jjs, jje ) ! (in)
590 
591  !$omp parallel do &
592  !$omp private(ap)
593  do j = jjs, jje
594  do i = iis, iie
595 
596  ap = - dt * 0.25_rp * ( dens(ks ,i ,j)*km(ks ,i ,j) &
597  + dens(ks+1,i ,j)*km(ks+1,i ,j) &
598  + dens(ks ,i+1,j)*km(ks ,i+1,j) &
599  + dens(ks+1,i+1,j)*km(ks+1,i+1,j) ) &
600  * rfdz(ks) / gsqrt(ks,i,j,i_uyw)
601  a(ks,i,j) = ap * rcdz(ks) / gsqrt(ks,i,j,i_uyz)
602  c(ks,i,j) = 0.0_rp
603  b(ks,i,j) = - a(ks,i,j) + 0.5_rp * ( dens(ks,i,j)+dens(ks,i+1,j) )
604  do k = ks+1, ke-1
605  c(k,i,j) = ap * rcdz(k) / gsqrt(k,i,j,i_uyz)
606  ap = - dt * 0.25_rp * ( dens(k ,i ,j)*km(k ,i ,j) &
607  + dens(k+1,i ,j)*km(k+1,i ,j) &
608  + dens(k ,i+1,j)*km(k ,i+1,j) &
609  + dens(k+1,i+1,j)*km(k+1,i+1,j) ) &
610  * rfdz(k) / gsqrt(k,i,j,i_uyw)
611  a(k,i,j) = ap * rcdz(k) / gsqrt(k,i,j,i_uyz)
612  b(k,i,j) = - a(k,i,j) - c(k,i,j) + 0.5_rp * ( dens(k,i,j)+dens(k,i+1,j) )
613  end do
614  a(ke,i,j) = 0.0_rp
615  c(ke,i,j) = ap * rcdz(ke) / gsqrt(ke,i,j,i_uyz)
616  b(ke,i,j) = - c(ke,i,j) + 0.5_rp * ( dens(ke,i,j)+dens(ke,i+1,j) )
617  end do
618  end do
619 
620  call matrix_solver_tridiagonal( ka, ks, ke, ia, iis, iie, ja, jjs, jje, &
621  a(:,:,:), b(:,:,:), c(:,:,:), tend(:,:,:), & ! (in)
622  tend2(:,:,:) ) ! (out)
623 
624  !$omp parallel do
625  do j = jjs, jje
626  do i = iis, iie
627  do k = ks, ke-1
628  qflx_sgs_momx(k,i,j,zdir) = qflx_sgs_momx(k,i,j,zdir) &
629  - 0.25_rp * ( dens(k ,i ,j)*km(k ,i ,j) &
630  + dens(k+1,i ,j)*km(k+1,i ,j) &
631  + dens(k ,i+1,j)*km(k ,i+1,j) &
632  + dens(k+1,i+1,j)*km(k+1,i+1,j) ) &
633  * dt * ( tend2(k+1,i,j) - tend2(k,i,j) ) * rfdz(k) / gsqrt(k,i,j,i_uyw)
634  end do
635 
636  end do
637  end do
638 
639  end if
640 
641  !##### momentum equation (y) #####
642  ! (x edge)
643  !$omp parallel do
644  do j = jjs, jje
645  do i = iis, iie
646  do k = ks, ke-1
647 #ifdef DEBUG
648  call check( __line__, dens(k,i,j) )
649  call check( __line__, dens(k,i,j+1) )
650  call check( __line__, dens(k+1,i,j) )
651  call check( __line__, dens(k+1,i,j+1) )
652  call check( __line__, km(k,i,j) )
653  call check( __line__, km(k,i,j+1) )
654  call check( __line__, km(k+1,i,j) )
655  call check( __line__, km(k+1,i,j+1) )
656  call check( __line__, s23_x(k,i,j) )
657 #endif
658  qflx_sgs_momy(k,i,j,zdir) = - 0.125_rp & ! 2/4/4
659  * ( dens(k,i,j)+dens(k+1,i,j)+dens(k,i,j+1)+dens(k+1,i,j+1) ) &
660  * ( km(k,i,j)+km(k+1,i,j)+km(k,i,j+1)+km(k+1,i,j+1) ) &
661  * s23_x(k,i,j)
662  enddo
663  enddo
664  enddo
665 #ifdef DEBUG
666  i = iundef; j = iundef; k = iundef
667 #endif
668  !$omp parallel do
669  do j = jjs, jje
670  do i = iis, iie
671  qflx_sgs_momy(ks-1,i,j,zdir) = 0.0_rp ! bottom boundary
672  qflx_sgs_momy(ke ,i,j,zdir) = 0.0_rp ! top boundary
673  enddo
674  enddo
675 #ifdef DEBUG
676  i = iundef; j = iundef; k = iundef
677 #endif
678 
679  ! (z edge)
680  !$omp parallel do
681  do j = jjs, jje
682  do i = iis-1, iie
683  do k = ks, ke
684 #ifdef DEBUG
685  call check( __line__, dens(k,i,j) )
686  call check( __line__, dens(k,i+1,j) )
687  call check( __line__, dens(k,i,j+1) )
688  call check( __line__, dens(k,i+1,j+1) )
689  call check( __line__, km(k,i,j) )
690  call check( __line__, km(k,i+1,j) )
691  call check( __line__, km(k,i,j+1) )
692  call check( __line__, km(k,i+1,j+1) )
693  call check( __line__, s12_z(k,i,j) )
694 #endif
695  qflx_sgs_momy(k,i,j,xdir) = - 0.125_rp & !
696  * ( dens(k,i,j)+dens(k,i+1,j)+dens(k,i,j+1)+dens(k,i+1,j+1) ) &
697  * ( km(k,i,j)+km(k,i+1,j)+km(k,i,j+1)+km(k,i+1,j+1) ) &
698  * s12_z(k,i,j)
699  enddo
700  enddo
701  enddo
702 #ifdef DEBUG
703  i = iundef; j = iundef; k = iundef
704 #endif
705 
706  ! (z-x plane)
707  !$omp parallel do
708  do j = jjs, jje+1
709  do i = iis, iie
710  do k = ks, ke
711 #ifdef DEBUG
712  call check( __line__, dens(k,i,j) )
713  call check( __line__, km(k,i,j) )
714  call check( __line__, s22_c(k,i,j) )
715 #endif
716  qflx_sgs_momy(k,i,j,ydir) = dens(k,i,j) * ( - 2.0_rp * km(k,i,j) * s22_c(k,i,j) )
717  enddo
718  enddo
719  enddo
720 #ifdef DEBUG
721  i = iundef; j = iundef; k = iundef
722 #endif
723 
724  if ( atmos_phy_tb_d1980_implicit ) then
725  call calc_tend_momy( tend, & ! (out)
726  qflx_sgs_momy, & ! (in)
727  gsqrt, j13g, j23g, j33g, mapf, & ! (in)
728  iis, iie, jjs, jje ) ! (in)
729 
730  !$omp parallel do &
731  !$omp private(ap)
732  do j = jjs, jje
733  do i = iis, iie
734 
735  ap = - dt * 0.25_rp * ( dens(ks ,i,j )*km(ks ,i,j ) &
736  + dens(ks+1,i,j )*km(ks+1,i,j ) &
737  + dens(ks ,i,j+1)*km(ks ,i,j+1) &
738  + dens(ks+1,i,j+1)*km(ks+1,i,j+1) ) &
739  * rfdz(ks) / gsqrt(ks,i,j,i_xvw)
740  a(ks,i,j) = ap * rcdz(ks) / gsqrt(ks,i,j,i_xvz)
741  c(ks,i,j) = 0.0_rp
742  b(ks,i,j) = - a(ks,i,j) + 0.5_rp * ( dens(ks,i,j)+dens(ks,i,j+1) )
743  do k = ks+1, ke-1
744  c(k,i,j) = ap * rcdz(k) / gsqrt(k,i,j,i_xvz)
745  ap = - dt * 0.25_rp * ( dens(k ,i,j )*km(k ,i,j ) &
746  + dens(k+1,i,j )*km(k+1,i,j ) &
747  + dens(k ,i,j+1)*km(k ,i,j+1) &
748  + dens(k+1,i,j+1)*km(k+1,i,j+1) ) &
749  * rfdz(k) / gsqrt(k,i,j,i_xvw)
750  a(k,i,j) = ap * rcdz(k) / gsqrt(k,i,j,i_xvz)
751  b(k,i,j) = - a(k,i,j) - c(k,i,j) + 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) )
752  end do
753  a(ke,i,j) = 0.0_rp
754  c(ke,i,j) = ap * rcdz(ke) / gsqrt(ke,i,j,i_xvz)
755  b(ke,i,j) = - c(ke,i,j) + 0.5_rp * ( dens(ke,i,j)+dens(ke,i,j+1) )
756  end do
757  end do
758 
759  call matrix_solver_tridiagonal( ka, ks, ke, ia, iis, iie, ja, jjs, jje, &
760  a(:,:,:), b(:,:,:), c(:,:,:), tend(:,:,:), & ! (in)
761  tend2(:,:,:) ) ! (out)
762 
763  !$omp parallel do
764  do j = jjs, jje
765  do i = iis, iie
766  do k = ks, ke-1
767  qflx_sgs_momy(k,i,j,zdir) = qflx_sgs_momy(k,i,j,zdir) &
768  - 0.25_rp * ( dens(k ,i,j )*km(k ,i,j ) &
769  + dens(k+1,i,j )*km(k+1,i,j ) &
770  + dens(k ,i,j+1)*km(k ,i,j+1) &
771  + dens(k+1,i,j+1)*km(k+1,i,j+1) ) &
772  * dt * ( tend2(k+1,i,j) - tend2(k,i,j) ) * rfdz(k) / gsqrt(k,i,j,i_xvw)
773  end do
774 
775  end do
776  end do
777 
778  end if
779 
780  !##### Thermodynamic Equation #####
781 
782  if ( atmos_phy_tb_d1980_implicit ) then
783 
784  !$omp parallel do &
785  !$omp private(ap)
786  do j = jjs, jje
787  do i = iis, iie
788 
789  ap = - dt * 0.25_rp * ( dens(ks,i,j)+dens(ks+1,i,j) ) &
790  * ( kh(ks,i,j)+kh(ks+1,i,j) ) &
791  * rfdz(ks) / gsqrt(ks,i,j,i_xyw)
792  a(ks,i,j) = ap * rcdz(ks) / gsqrt(ks,i,j,i_xyz)
793  c(ks,i,j) = 0.0_rp
794  b(ks,i,j) = - a(ks,i,j) + dens(ks,i,j)
795  do k = ks+1, ke-1
796  c(k,i,j) = ap * rcdz(k) / gsqrt(k,i,j,i_xyz)
797  ap = - dt * 0.25_rp * ( dens(k,i,j)+dens(k+1,i,j) ) &
798  * ( kh(k,i,j)+kh(k+1,i,j) ) &
799  * rfdz(k) / gsqrt(k,i,j,i_xyw)
800  a(k,i,j) = ap * rcdz(k) / gsqrt(k,i,j,i_xyz)
801  b(k,i,j) = - a(k,i,j) - c(k,i,j) + dens(k,i,j)
802  end do
803  a(ke,i,j) = 0.0_rp
804  c(ke,i,j) = ap * rcdz(ke) / gsqrt(ke,i,j,i_xyz)
805  b(ke,i,j) = - c(ke,i,j) + dens(ke,i,j)
806 
807  end do
808  end do
809 
810  end if
811 
812  call calc_flux_phi( &
813  qflx_sgs_rhot, &
814  dens, pott, kh, 1.0_rp, &
815  gsqrt, j13g, j23g, j33g, mapf, &
816  .false., &
817  atmos_phy_tb_d1980_implicit, &
818  a, b, c, dt, &
819  iis, iie, jjs, jje )
820 
821  enddo
822  enddo
823 
824 
825  !##### Tracers #####
826  do iq = 1, qa
827 
828  if ( iq == i_tke ) then
829  qflx_sgs_rhoq(:,:,:,:,iq) = 0.0_rp
830  cycle
831  end if
832  if ( .not. tracer_advc(iq) ) cycle
833 
834  do jjs = js, je, jblock
835  jje = jjs+jblock-1
836  do iis = is, ie, iblock
837  iie = iis+iblock-1
838 
839  call calc_flux_phi( &
840  qflx_sgs_rhoq(:,:,:,:,iq), &
841  dens, qtrc(:,:,:,iq), kh, 1.0_rp, &
842  gsqrt, j13g, j23g, j33g, mapf, &
843  .false., &
844  atmos_phy_tb_d1980_implicit, &
845  a, b, c, dt, &
846  iis, iie, jjs, jje )
847 
848  enddo
849  enddo
850 #ifdef DEBUG
851  iis = iundef; iie = iundef; jjs = iundef; jje = iundef
852 #endif
853 
854  enddo ! scalar quantities loop
855 #ifdef DEBUG
856  iq = iundef
857 #endif
858 
859  ! TKE
860  do jjs = js, je, jblock
861  jje = jjs+jblock-1
862  do iis = is, ie, iblock
863  iie = iis+iblock-1
864 
865  if ( atmos_phy_tb_d1980_implicit ) then
866 
867  !$omp parallel do &
868  !$omp private(ap)
869  do j = jjs, jje
870  do i = iis, iie
871 
872  ap = - dt * 0.25_rp * ( dens(ks,i,j)+dens(ks+1,i,j) ) &
873  * 2.0_rp * ( km(ks,i,j)+km(ks+1,i,j) ) &
874  * rfdz(ks) / gsqrt(ks,i,j,i_xyw)
875  a(ks,i,j) = ap * rcdz(ks) / gsqrt(ks,i,j,i_xyz)
876  c(ks,i,j) = 0.0_rp
877  b(ks,i,j) = - a(ks,i,j) + dens(ks,i,j)
878  do k = ks+1, ke-1
879  c(k,i,j) = ap * rcdz(k) / gsqrt(k,i,j,i_xyz)
880  ap = - dt * 0.25_rp * ( dens(k,i,j)+dens(k+1,i,j) ) &
881  * 2.0_rp * ( km(k,i,j)+km(k+1,i,j) ) &
882  * rfdz(k) / gsqrt(k,i,j,i_xyw)
883  a(k,i,j) = ap * rcdz(k) / gsqrt(k,i,j,i_xyz)
884  b(k,i,j) = - a(k,i,j) - c(k,i,j) + dens(k,i,j)
885  end do
886  a(ke,i,j) = 0.0_rp
887  c(ke,i,j) = ap * rcdz(ke) / gsqrt(ke,i,j,i_xyz)
888  b(ke,i,j) = - c(ke,i,j) + dens(ke,i,j)
889 
890  end do
891  end do
892 
893  end if
894 
895  call calc_flux_phi( &
896  qflx_tke, & ! (out)
897  dens, qtrc(:,:,:,i_tke), km, 2.0_rp, & ! (in)
898  gsqrt, j13g, j23g, j33g, mapf, & ! (in)
899  .false., & ! (in)
900  atmos_phy_tb_d1980_implicit, & ! (in)
901  a, b, c, dt, & ! (in)
902  iis, iie, jjs, jje ) ! (in)
903 
904  call calc_tend_phi( rhoq_t(:,:,:,i_tke), & ! (out)
905  qflx_tke, & ! (in)
906  gsqrt, j13g, j23g, j33g, mapf, & ! (in)
907  iis, iie, jjs, jje ) ! (in)
908  !$omp parallel do
909  do j = jjs, jje
910  do i = iis, iie
911  do k = ks, ke
912  rhoq_t(k,i,j,i_tke) = max( &
913  rhoq_t(k,i,j,i_tke) &
914  + ( km(k,i,j) * s2(k,i,j) &
915  - kh(k,i,j) * n2(k,i,j) &
916  - ( c1 + c2*l(k,i,j)/delta(k,i,j) ) * qtrc(k,i,j,i_tke)**(1.5_rp) / l(k,i,j) ) * dens(k,i,j), &
917  - qtrc(k,i,j,i_tke) * dens(k,i,j) / real(dt,kind=rp) )
918  enddo
919  enddo
920  enddo
921 #ifdef DEBUG
922  i = iundef; j = iundef; k = iundef
923 #endif
924 
925  end do
926  end do
927 
928 
929  return
930  end subroutine atmos_phy_tb_d1980
931 
932 end module scale_atmos_phy_tb_d1980
scale_const::const_grav
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:49
scale_atmos_grid_cartesc_index::ke
integer, public ke
end point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:52
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_tracer::qa
integer, public qa
Definition: scale_tracer.F90:35
scale_atmos_grid_cartesc::atmos_grid_cartesc_rfdx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rfdx
reciprocal of face-dx
Definition: scale_atmos_grid_cartesC.F90:68
scale_const::const_undef2
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:40
scale_atmos_phy_tb_common::atmos_phy_tb_calc_strain_tensor
subroutine, public atmos_phy_tb_calc_strain_tensor(S33_C, S11_C, S22_C, S31_C, S12_C, S23_C, S12_Z, S23_X, S31_Y, S2, DENS, MOMZ, MOMX, MOMY, GSQRT, J13G, J23G, J33G, MAPF)
Definition: scale_atmos_phy_tb_common.F90:66
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_atmos_grid_cartesc_index::ka
integer, public ka
Definition: scale_atmos_grid_cartesC_index.F90:47
scale_atmos_grid_cartesc_index::jblock
integer, public jblock
block size for cache blocking: y
Definition: scale_atmos_grid_cartesC_index.F90:41
scale_atmos_grid_cartesc::atmos_grid_cartesc_rcdx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rcdx
reciprocal of center-dx
Definition: scale_atmos_grid_cartesC.F90:66
scale_tracer::tracer_advc
logical, dimension(qa_max), public tracer_advc
Definition: scale_tracer.F90:46
scale_atmos_phy_tb_common::atmos_phy_tb_calc_tend_momx
subroutine, public atmos_phy_tb_calc_tend_momx(MOMX_t_TB, QFLX_MOMX, GSQRT, J13G, J23G, J33G, MAPF, IIS, IIE, JJS, JJE)
Definition: scale_atmos_phy_tb_common.F90:1679
scale_atmos_grid_cartesc_index::iblock
integer, public iblock
block size for cache blocking: x
Definition: scale_atmos_grid_cartesC_index.F90:40
scale_atmos_phy_tb_common::atmos_phy_tb_calc_tend_momy
subroutine, public atmos_phy_tb_calc_tend_momy(MOMY_t_TB, QFLX_MOMY, GSQRT, J13G, J23G, J33G, MAPF, IIS, IIE, JJS, JJE)
Definition: scale_atmos_phy_tb_common.F90:1741
scale_atmos_phy_tb_d1980::atmos_phy_tb_d1980_config
subroutine, public atmos_phy_tb_d1980_config(TYPE_TB, I_TKE_out)
Config.
Definition: scale_atmos_phy_tb_d1980.F90:81
scale_atmos_grid_cartesc::atmos_grid_cartesc_rcdz
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rcdz
reciprocal of center-dz
Definition: scale_atmos_grid_cartesC.F90:45
scale_atmos_phy_tb_d1980::atmos_phy_tb_d1980_finalize
subroutine, public atmos_phy_tb_d1980_finalize
finalize
Definition: scale_atmos_phy_tb_d1980.F90:167
scale_atmos_grid_cartesc_index::i_xyz
integer, public i_xyz
Definition: scale_atmos_grid_cartesC_index.F90:91
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_atmos_phy_tb_d1980::atmos_phy_tb_d1980
subroutine, public atmos_phy_tb_d1980(qflx_sgs_momz, qflx_sgs_momx, qflx_sgs_momy, qflx_sgs_rhot, qflx_sgs_rhoq, RHOQ_t, Km, 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)
Definition: scale_atmos_phy_tb_d1980.F90:182
scale_atmos_grid_cartesc_index::ie
integer, public ie
end point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:54
scale_atmos_grid_cartesc_index::i_uyz
integer, public i_uyz
Definition: scale_atmos_grid_cartesC_index.F90:95
scale_io
module STDIO
Definition: scale_io.F90:10
scale_atmos_grid_cartesc::atmos_grid_cartesc_rcdy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rcdy
reciprocal of center-dy
Definition: scale_atmos_grid_cartesC.F90:67
scale_atmos_grid_cartesc_index::i_xvw
integer, public i_xvw
Definition: scale_atmos_grid_cartesC_index.F90:94
scale_tracer::k
real(rp), public k
Definition: scale_tracer.F90:45
scale_atmos_grid_cartesc_index
module atmosphere / grid / cartesC index
Definition: scale_atmos_grid_cartesC_index.F90:12
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
scale_debug::check
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:59
scale_atmos_grid_cartesc_index::zdir
integer, parameter, public zdir
Definition: scale_atmos_grid_cartesC_index.F90:32
scale_atmos_grid_cartesc_index::i_xvz
integer, public i_xvz
Definition: scale_atmos_grid_cartesC_index.F90:96
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_atmos_phy_tb_common::atmos_phy_tb_calc_flux_phi
subroutine, public atmos_phy_tb_calc_flux_phi(qflx_phi, DENS, PHI, Kh, FACT, GSQRT, J13G, J23G, J33G, MAPF, horizontal, implicit, a, b, c, dt, IIS, IIE, JJS, JJE)
Definition: scale_atmos_phy_tb_common.F90:1337
scale_atmos_grid_cartesc_index::i_uyw
integer, public i_uyw
Definition: scale_atmos_grid_cartesC_index.F90:93
scale_atmos_grid_cartesc_index::is
integer, public is
start point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:53
scale_precision::dp
integer, parameter, public dp
Definition: scale_precision.F90:32
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
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:68
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
scale_atmos_phy_tb_common::atmos_phy_tb_calc_tend_phi
subroutine, public atmos_phy_tb_calc_tend_phi(phi_t_TB, QFLX_phi, GSQRT, J13G, J23G, J33G, MAPF, IIS, IIE, JJS, JJE)
Definition: scale_atmos_phy_tb_common.F90:1804
scale_atmos_grid_cartesc_index::xdir
integer, parameter, public xdir
Definition: scale_atmos_grid_cartesC_index.F90:33
scale_atmos_phy_tb_d1980::atmos_phy_tb_d1980_setup
subroutine, public atmos_phy_tb_d1980_setup(CDZ, CDX, CDY, CZ)
Setup.
Definition: scale_atmos_phy_tb_d1980.F90:115
scale_atmos_grid_cartesc_index::ks
integer, public ks
start point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:51
scale_debug
module DEBUG
Definition: scale_debug.F90:11
scale_atmos_phy_tb_common
module ATMOSPHERE / Physics Turbulence
Definition: scale_atmos_phy_tb_common.F90:12
scale_atmos_grid_cartesc_index::js
integer, public js
start point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:55
scale_atmos_grid_cartesc_index::ydir
integer, parameter, public ydir
Definition: scale_atmos_grid_cartesC_index.F90:34
scale_matrix
module MATRIX
Definition: scale_matrix.F90:17
scale_atmos_grid_cartesc::atmos_grid_cartesc_rfdz
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rfdz
reciprocal of face-dz
Definition: scale_atmos_grid_cartesC.F90:46
scale_atmos_grid_cartesc
module atmosphere / grid / cartesC
Definition: scale_atmos_grid_cartesC.F90:12
scale_atmos_grid_cartesc::atmos_grid_cartesc_fdz
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fdz
z-length of grid(i+1) to grid(i) [m]
Definition: scale_atmos_grid_cartesC.F90:44
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_atmos_phy_tb_d1980
module ATMOSPHERE / Physics Turbulence
Definition: scale_atmos_phy_tb_d1980.F90:18
scale_atmos_phy_tb_common::atmos_phy_tb_calc_tend_momz
subroutine, public atmos_phy_tb_calc_tend_momz(MOMZ_t_TB, QFLX_MOMZ, GSQRT, J13G, J23G, J33G, MAPF, IIS, IIE, JJS, JJE)
Definition: scale_atmos_phy_tb_common.F90:1616
scale_atmos_grid_cartesc::atmos_grid_cartesc_rfdy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rfdy
reciprocal of face-dy
Definition: scale_atmos_grid_cartesC.F90:69
scale_io::io_fid_conf
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:57
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56
scale_atmos_grid_cartesc_index::i_xyw
integer, public i_xyw
Definition: scale_atmos_grid_cartesC_index.F90:92