SCALE-RM
mod_atmos_phy_tb_driver.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
11 #include "scalelib.h"
13  !-----------------------------------------------------------------------------
14  !
15  !++ used modules
16  !
17  use scale_precision
18  use scale_io
19  use scale_prof
21  use scale_tracer
22  !-----------------------------------------------------------------------------
23  implicit none
24  private
25  !-----------------------------------------------------------------------------
26  !
27  !++ Public procedure
28  !
32 
33  !-----------------------------------------------------------------------------
34  !
35  !++ Public parameters & variables
36  !
37  !-----------------------------------------------------------------------------
38  !
39  !++ Private procedure
40  !
41  !-----------------------------------------------------------------------------
42  !
43  !++ Private parameters & variables
44  !
45  integer, private :: monit_west
46  integer, private :: monit_east
47  integer, private :: monit_south
48  integer, private :: monit_north
49 
50  !-----------------------------------------------------------------------------
51 contains
52  !-----------------------------------------------------------------------------
55  use scale_prc, only: &
56  prc_abort
57  use scale_atmos_phy_tb_d1980, only: &
59  use scale_atmos_phy_tb_dns, only: &
61  use mod_atmos_admin, only: &
64  use mod_atmos_phy_tb_vars, only: &
65  i_tke
66  implicit none
67  !---------------------------------------------------------------------------
68 
69  if ( .not. atmos_sw_phy_tb ) return
70 
71  log_newline
72  log_info("ATMOS_PHY_TB_driver_tracer_setup",*) 'Setup'
73 
74  select case( atmos_phy_tb_type )
75  case( 'SMAGORINSKY' )
76  i_tke = -1
77  case( 'D1980' )
79  i_tke ) ! [OUT]
80  case( 'DNS' )
82  i_tke ) ! [OUT]
83  case('OFF')
84  ! do nothing
85  case default
86  log_error("ATMOS_PHY_TB_driver_tracer_setup",*) 'ATMOS_PHY_TB_TYPE is invalid: ', atmos_phy_tb_type
87  call prc_abort
88  end select
89 
90  return
92 
93  !-----------------------------------------------------------------------------
95  subroutine atmos_phy_tb_driver_setup
97  cdz => atmos_grid_cartesc_cdz, &
98  cdx => atmos_grid_cartesc_cdx, &
100  use scale_atmos_grid_cartesc_real, only: &
101  real_cz => atmos_grid_cartesc_real_cz, &
102  real_fz => atmos_grid_cartesc_real_fz
105  use scale_atmos_phy_tb_smg, only: &
107  use scale_atmos_phy_tb_d1980, only: &
109  use scale_atmos_phy_tb_dns, only: &
111  use mod_atmos_admin, only: &
113  atmos_sw_phy_tb, &
115  use mod_atmos_phy_tb_vars, only: &
116  momz_t_tb => atmos_phy_tb_momz_t
117  use scale_monitor, only: &
119  implicit none
120 
121  integer :: i, j
122  logical :: horizontal = .false.
123  !---------------------------------------------------------------------------
124 
125  if ( .not. atmos_sw_phy_tb ) return
126 
127  log_newline
128  log_info("ATMOS_PHY_TB_driver_setup",*) 'Setup'
129 
130  ! initialize
131  !$omp parallel do
132  do j = js, je
133  do i = is, ie
134  momz_t_tb(ks-1,i,j) = 0.0_rp
135  momz_t_tb(ke ,i,j) = 0.0_rp
136  enddo
137  enddo
138 
139  ! setup library component
140  select case( atmos_phy_tb_type )
141  case( 'SMAGORINSKY' )
142  if ( atmos_sw_phy_bl ) horizontal = .true.
143  call atmos_phy_tb_smg_setup( real_fz, real_cz, cdx, cdy, mapf(:,:,:,i_xy), horizontal ) ! [IN]
144  case( 'D1980' )
145  call atmos_phy_tb_d1980_setup( cdz, cdx, cdy, real_cz ) ! [IN]
146  case( 'DNS' )
147  call atmos_phy_tb_dns_setup( cdz, cdx, cdy, real_cz ) ! [IN]
148  end select
149 
150 
151  ! monitor
152  call monitor_reg( "QTOTFLX_TB_WEST", "water mass flux at the western boundary", "kg", & ! [IN]
153  monit_west, & ! [OUT]
154  dim_type="ZY-W", is_tendency=.true. ) ! [IN]
155  call monitor_reg( "QTOTFLX_TB_EAST", "water mass flux at the eastern boundary", "kg", & ! [IN]
156  monit_east, & ! [OUT]
157  dim_type="ZY-E", is_tendency=.true. ) ! [IN]
158  call monitor_reg( "QTOTFLX_TB_SOUTH", "water mass flux at the southern boundary", "kg", & ! [IN]
159  monit_south, & ! [OUT]
160  dim_type="ZX-S", is_tendency=.true. ) ! [IN]
161  call monitor_reg( "QTOTFLX_TB_NORTH", "water mass flux at the northern boundary", "kg", & ! [IN]
162  monit_north, & ! [OUT]
163  dim_type="ZX-N", is_tendency=.true. ) ! [IN]
164 
165 
166  return
167  end subroutine atmos_phy_tb_driver_setup
168 
169  !-----------------------------------------------------------------------------
171  subroutine atmos_phy_tb_driver_calc_tendency( update_flag )
172  use scale_prc_cartesc, only: &
173  prc_has_w, &
174  prc_has_e, &
175  prc_has_s, &
176  prc_has_n
177  use scale_statistics, only: &
179  statistics_total
180  use scale_atmos_grid_cartesc, only: &
181  fdz => atmos_grid_cartesc_fdz, &
182  rcdz => atmos_grid_cartesc_rcdz, &
183  rfdz => atmos_grid_cartesc_rfdz, &
184  cdx => atmos_grid_cartesc_cdx, &
185  fdx => atmos_grid_cartesc_fdx, &
186  cdy => atmos_grid_cartesc_cdy, &
188  use scale_atmos_grid_cartesc_real, only: &
204  use scale_file_history, only: &
205  file_history_in
206  use scale_monitor, only: &
207  monitor_put
208  use scale_time, only: &
209  dt_tb => time_dtsec_atmos_phy_tb
210  use scale_atmos_phy_tb_smg, only: &
212  use scale_atmos_phy_tb_d1980, only: &
214  use scale_atmos_phy_tb_dns, only: &
216  use scale_atmos_phy_tb_common, only: &
217  calc_tend_momz => atmos_phy_tb_calc_tend_momz, &
218  calc_tend_momx => atmos_phy_tb_calc_tend_momx, &
219  calc_tend_momy => atmos_phy_tb_calc_tend_momy, &
220  calc_tend_phi => atmos_phy_tb_calc_tend_phi
221  use mod_atmos_admin, only: &
223  use mod_atmos_vars, only: &
224  atmos_vars_get_diagnostic, &
225  dens => dens_av, &
226  momz => momz_av, &
227  momx => momx_av, &
228  momy => momy_av, &
229  rhot => rhot_av, &
230  pott, &
231  qtrc => qtrc_av, &
232  momz_t => momz_tp, &
233  momx_t => momx_tp, &
234  momy_t => momy_tp, &
235  rhot_t => rhot_tp, &
236  rhoq_t => rhoq_tp
237  use mod_atmos_phy_tb_vars, only: &
238  i_tke, &
239  momz_t_tb => atmos_phy_tb_momz_t, &
240  momx_t_tb => atmos_phy_tb_momx_t, &
241  momy_t_tb => atmos_phy_tb_momy_t, &
242  rhot_t_tb => atmos_phy_tb_rhot_t, &
243  rhoq_t_tb => atmos_phy_tb_rhoq_t
244  use mod_atmos_phy_sf_vars, only: &
245  sflx_mw => atmos_phy_sf_sflx_mw, &
246  sflx_mu => atmos_phy_sf_sflx_mu, &
247  sflx_mv => atmos_phy_sf_sflx_mv, &
248  sflx_sh => atmos_phy_sf_sflx_sh, &
249  sflx_q => atmos_phy_sf_sflx_qtrc
250  implicit none
251 
252  logical, intent(in) :: update_flag
253 
254  ! eddy viscosity/diffusion flux
255  real(rp) :: qflx_momz(ka,ia,ja,3)
256  real(rp) :: qflx_momx(ka,ia,ja,3)
257  real(rp) :: qflx_momy(ka,ia,ja,3)
258  real(rp) :: qflx_rhot(ka,ia,ja,3)
259  real(rp) :: qflx_rhoq(ka,ia,ja,3,qa)
260 
261  real(rp) :: nu(ka,ia,ja) ! eddy viscosity
262  real(rp) :: ri(ka,ia,ja) ! Richardson number
263  real(rp) :: pr(ka,ia,ja) ! Prandtl number
264 
265  real(rp) :: n2(ka,ia,ja)
266 
267  real(rp) :: tend(ka,ia,ja)
268 
269  ! for monitor
270  real(rp) :: qflx_x(ka,ja)
271  real(rp) :: qflx_y(ka,ia)
272 
273  integer :: jjs, jje
274  integer :: iis, iie
275 
276  integer :: k, i, j, iq
277  !---------------------------------------------------------------------------
278 
279  if ( update_flag ) then
280 
281  call atmos_vars_get_diagnostic( "N2", n2 )
282 
283  select case( atmos_phy_tb_type )
284  case( 'SMAGORINSKY' )
285  call atmos_phy_tb_smg( qflx_momz, qflx_momx, qflx_momy, & ! [OUT]
286  qflx_rhot, qflx_rhoq, & ! [OUT]
287  momz_t_tb, momx_t_tb, momy_t_tb, & ! [OUT]
288  rhot_t_tb, rhoq_t_tb, & ! [OUT]
289  nu, ri, pr, & ! [OUT]
290  momz, momx, momy, pott, dens, qtrc, n2, & ! [IN]
291  fz, fdz, rcdz, rfdz, cdx, fdx, cdy, fdy, & ! [IN]
292  gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
293  dt_tb ) ! [IN]
294  case( 'D1980' )
295  momz_t_tb(:,:,:) = 0.0_rp
296  momx_t_tb(:,:,:) = 0.0_rp
297  momy_t_tb(:,:,:) = 0.0_rp
298  rhot_t_tb(:,:,:) = 0.0_rp
299  call atmos_phy_tb_d1980( qflx_momz, qflx_momx, qflx_momy, & ! [OUT]
300  qflx_rhot, qflx_rhoq, & ! [OUT]
301  rhoq_t_tb, & ! [OUT]
302  nu, ri, pr, & ! [OUT]
303  momz, momx, momy, rhot, dens, qtrc, n2, & ! [IN]
304  sflx_mw, sflx_mu, sflx_mv, & ! [IN]
305  sflx_sh, sflx_q, & ! [IN]
306  gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
307  dt_tb ) ! [IN]
308  case( 'DNS' )
309  momz_t_tb(:,:,:) = 0.0_rp
310  momx_t_tb(:,:,:) = 0.0_rp
311  momy_t_tb(:,:,:) = 0.0_rp
312  rhot_t_tb(:,:,:) = 0.0_rp
313  call atmos_phy_tb_dns( qflx_momz, qflx_momx, qflx_momy, & ! [OUT]
314  qflx_rhot, qflx_rhoq, & ! [OUT]
315  rhoq_t_tb, & ! [OUT]
316  nu, ri, pr, & ! [OUT]
317  momz, momx, momy, rhot, dens, qtrc, n2, & ! [IN]
318  sflx_mw, sflx_mu, sflx_mv, & ! [IN]
319  sflx_sh, sflx_q, & ! [IN]
320  gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
321  dt_tb ) ! [IN]
322  end select
323 
324 
325  do jjs = js, je, jblock
326  jje = jjs+jblock-1
327  do iis = is, ie, iblock
328  iie = iis+iblock-1
329 
330  call calc_tend_momz( tend(:,:,:), & ! (out)
331  qflx_momz, & ! (in)
332  gsqrt, j13g, j23g, j33g, mapf, & ! (in)
333  iis, iie, jjs, jje ) ! (in)
334  !$omp parallel do
335  do j = jjs, jje
336  do i = iis, iie
337  do k = ks, ke-1
338  momz_t_tb(k,i,j) = momz_t_tb(k,i,j) + tend(k,i,j)
339  end do
340  end do
341  end do
342 
343  call calc_tend_momx( tend(:,:,:), & ! (out)
344  qflx_momx, & ! (in)
345  gsqrt, j13g, j23g, j33g, mapf, & ! (in)
346  iis, iie, jjs, jje ) ! (in)
347  !$omp parallel do
348  do j = jjs, jje
349  do i = iis, iie
350  do k = ks, ke
351  momx_t_tb(k,i,j) = momx_t_tb(k,i,j) + tend(k,i,j)
352  end do
353  end do
354  end do
355 
356  call calc_tend_momy( tend(:,:,:), & ! (out)
357  qflx_momy, & ! (in)
358  gsqrt, j13g, j23g, j33g, mapf, & ! (in)
359  iis, iie, jjs, jje ) ! (in)
360  !$omp parallel do
361  do j = jjs, jje
362  do i = iis, iie
363  do k = ks, ke
364  momy_t_tb(k,i,j) = momy_t_tb(k,i,j) + tend(k,i,j)
365  end do
366  end do
367  end do
368 
369  call calc_tend_phi ( tend(:,:,:), & ! (out)
370  qflx_rhot, & ! (in)
371  gsqrt, j13g, j23g, j33g, mapf, & ! (in)
372  iis, iie, jjs, jje ) ! (in)
373  !$omp parallel do
374  do j = jjs, jje
375  do i = iis, iie
376  do k = ks, ke
377  rhot_t_tb(k,i,j) = rhot_t_tb(k,i,j) + tend(k,i,j)
378  end do
379  end do
380  end do
381 
382  do iq = 1, qa
383  if ( iq == i_tke .or. .not. tracer_advc(iq) ) cycle
384 
385  call calc_tend_phi( tend(:,:,:), & ! (out)
386  qflx_rhoq(:,:,:,:,iq), & ! (in)
387  gsqrt, j13g, j23g, j33g, mapf, & ! (in)
388  iis, iie, jjs, jje ) ! (in)
389 
390  !$omp parallel do
391  do j = jjs, jje
392  do i = iis, iie
393  do k = ks, ke
394  rhoq_t_tb(k,i,j,iq) = rhoq_t_tb(k,i,j,iq) + tend(k,i,j)
395  end do
396  end do
397  end do
398 
399  end do
400 
401  end do
402  end do
403 
404 
405  call file_history_in( nu(:,:,:), 'NU', 'eddy viscosity', 'm2/s' , fill_halo=.true. )
406  call file_history_in( ri(:,:,:), 'Ri', 'Richardson number', 'NIL' , fill_halo=.true. )
407  call file_history_in( pr(:,:,:), 'Pr', 'Prantle number', 'NIL' , fill_halo=.true. )
408 
409  call file_history_in( momz_t_tb(:,:,:), 'MOMZ_t_TB', 'MOMZ tendency (TB)', 'kg/m2/s2', fill_halo=.true. )
410  call file_history_in( momx_t_tb(:,:,:), 'MOMX_t_TB', 'MOMX tendency (TB)', 'kg/m2/s2', fill_halo=.true. )
411  call file_history_in( momy_t_tb(:,:,:), 'MOMY_t_TB', 'MOMY tendency (TB)', 'kg/m2/s2', fill_halo=.true. )
412  call file_history_in( rhot_t_tb(:,:,:), 'RHOT_t_TB', 'RHOT tendency (TB)', 'K.kg/m3/s', fill_halo=.true. )
413 
414  do iq = 1, qa
415  call file_history_in( rhoq_t_tb(:,:,:,iq), trim(tracer_name(iq))//'_t_TB', &
416  'tendency rho*'//trim(tracer_name(iq))//' in TB', 'kg/m3/s', fill_halo=.true. )
417  enddo
418 
419  call file_history_in( qflx_momz(:,:,:,zdir), 'SGS_ZFLX_MOMZ', 'SGS Z FLUX of MOMZ', 'kg/m/s2', &
420  fill_halo=.true.)
421  call file_history_in( qflx_momz(:,:,:,xdir), 'SGS_XFLX_MOMZ', 'SGS X FLUX of MOMZ', 'kg/m/s2', &
422  dim_type='ZHXHY', fill_halo=.true.)
423  call file_history_in( qflx_momz(:,:,:,ydir), 'SGS_YFLX_MOMZ', 'SGS Y FLUX of MOMZ', 'kg/m/s2', &
424  dim_type='ZHXYH', fill_halo=.true.)
425 
426  call file_history_in( qflx_momx(:,:,:,zdir), 'SGS_ZFLX_MOMX', 'SGS Z FLUX of MOMX', 'kg/m/s2', &
427  dim_type='ZHXHY', fill_halo=.true.)
428  call file_history_in( qflx_momx(:,:,:,xdir), 'SGS_XFLX_MOMX', 'SGS X FLUX of MOMX', 'kg/m/s2', &
429  fill_halo=.true.)
430  call file_history_in( qflx_momx(:,:,:,ydir), 'SGS_YFLX_MOMX', 'SGS Y FLUX of MOMX', 'kg/m/s2', &
431  dim_type='ZXHYH', fill_halo=.true.)
432 
433  call file_history_in( qflx_momy(:,:,:,zdir), 'SGS_ZFLX_MOMY', 'SGS Z FLUX of MOMY', 'kg/m/s2', &
434  dim_type='ZHXYH', fill_halo=.true.)
435  call file_history_in( qflx_momy(:,:,:,xdir), 'SGS_XFLX_MOMY', 'SGS X FLUX of MOMY', 'kg/m/s2', &
436  dim_type='ZXHYH', fill_halo=.true.)
437  call file_history_in( qflx_momy(:,:,:,ydir), 'SGS_YFLX_MOMY', 'SGS Y FLUX of MOMY', 'kg/m/s2', &
438  fill_halo=.true.)
439 
440  call file_history_in( qflx_rhot(:,:,:,zdir), 'SGS_ZFLX_RHOT', 'SGS Z FLUX of RHOT', 'K*kg/m2/s', &
441  dim_type='ZHXY', fill_halo=.true.)
442  call file_history_in( qflx_rhot(:,:,:,xdir), 'SGS_XFLX_RHOT', 'SGS X FLUX of RHOT', 'K*kg/m2/s', &
443  dim_type='ZXHY', fill_halo=.true.)
444  call file_history_in( qflx_rhot(:,:,:,ydir), 'SGS_YFLX_RHOT', 'SGS Y FLUX of RHOT', 'K*kg/m2/s', &
445  dim_type='ZXYH', fill_halo=.true.)
446 
447 
448  do iq = 1, qa
449  if ( iq == i_tke .or. .not. tracer_advc(iq) ) cycle
450 
451  call file_history_in( qflx_rhoq(:,:,:,zdir,iq), &
452  'SGS_ZFLX_'//trim(tracer_name(iq)), 'SGS Z FLUX of '//trim(tracer_name(iq)), 'kg/m2/s', &
453  dim_type='ZHXY', fill_halo=.true.)
454  call file_history_in( qflx_rhoq(:,:,:,xdir,iq), &
455  'SGS_XFLX_'//trim(tracer_name(iq)), 'SGS X FLUX of '//trim(tracer_name(iq)), 'kg/m2/s', &
456  dim_type='ZXHY', fill_halo=.true.)
457  call file_history_in( qflx_rhoq(:,:,:,ydir,iq), &
458  'SGS_YFLX_'//trim(tracer_name(iq)), 'SGS Y FLUX of '//trim(tracer_name(iq)), 'kg/m2/s', &
459  dim_type='ZXYH', fill_halo=.true.)
460  end do
461 
462  if ( statistics_checktotal ) then
463  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
464  momz_t_tb(:,:,:), 'MOMZ_t_TB', &
467  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
468  momx_t_tb(:,:,:), 'MOMX_t_TB', &
471  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
472  momy_t_tb(:,:,:), 'MOMY_t_TB', &
475  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
476  rhot_t_tb(:,:,:), 'RHOT_t_TB', &
479  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
480  nu(:,:,:), 'Nu', &
483  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
484  ri(:,:,:), 'Ri', &
487  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
488  pr(:,:,:), 'Pr', &
491 
492  do iq = 1, qa
493  call statistics_total( ka, ks, ke, ia, is, ie, ja, js, je, &
494  rhoq_t_tb(:,:,:,iq), trim(tracer_name(iq))//'_t_TB', &
497  enddo
498  endif
499 
500  endif
501 
502  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
503  !$omp shared(JS,JE,IS,IE,KS,KE,MOMZ_t,MOMZ_t_TB,MOMX_t,MOMX_t_TB,MOMY_t,MOMY_t_TB,RHOT_t,RHOT_t_TB)
504  do j = js, je
505  do i = is, ie
506  do k = ks, ke
507  momz_t(k,i,j) = momz_t(k,i,j) + momz_t_tb(k,i,j)
508  momx_t(k,i,j) = momx_t(k,i,j) + momx_t_tb(k,i,j)
509  momy_t(k,i,j) = momy_t(k,i,j) + momy_t_tb(k,i,j)
510  rhot_t(k,i,j) = rhot_t(k,i,j) + rhot_t_tb(k,i,j)
511  enddo
512  enddo
513  enddo
514 
515  do iq = 1, qa
516 
517  if ( .not. ( iq == i_tke .or. tracer_advc(iq) ) ) cycle
518 
519  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
520  do j = js, je
521  do i = is, ie
522  do k = ks, ke
523  rhoq_t(k,i,j,iq) = rhoq_t(k,i,j,iq) + rhoq_t_tb(k,i,j,iq)
524  enddo
525  enddo
526  enddo
527 
528  enddo
529 
530  if ( monit_west > 0 ) then
531  qflx_x(:,:) = 0.0_rp
532  if ( .not. prc_has_w ) then
533  do iq = 1, qa
534  if ( tracer_advc(iq) .and. tracer_mass(iq) == 1.0_rp ) then
535  do j = js, je
536  do k = ks, ke
537  qflx_x(k,j) = qflx_x(k,j) + qflx_rhoq(k,is-1,j,xdir,iq)
538  end do
539  end do
540  end if
541  end do
542  end if
543  call monitor_put( monit_west, qflx_x(:,:) )
544  end if
545  if ( monit_east > 0 ) then
546  qflx_x(:,:) = 0.0_rp
547  if ( .not. prc_has_e ) then
548  do iq = 1, qa
549  if ( tracer_advc(iq) .and. tracer_mass(iq) == 1.0_rp ) then
550  do j = js, je
551  do k = ks, ke
552  qflx_x(k,j) = qflx_x(k,j) + qflx_rhoq(k,ie,j,xdir,iq)
553  end do
554  end do
555  end if
556  end do
557  end if
558  call monitor_put( monit_east, qflx_x(:,:) )
559  end if
560  if ( monit_south > 0 ) then
561  qflx_y(:,:) = 0.0_rp
562  if ( .not. prc_has_s ) then
563  do iq = 1, qa
564  if ( tracer_advc(iq) .and. tracer_mass(iq) == 1.0_rp ) then
565  do i = is, ie
566  do k = ks, ke
567  qflx_y(k,i) = qflx_y(k,i) + qflx_rhoq(k,i,js-1,ydir,iq)
568  end do
569  end do
570  end if
571  end do
572  end if
573  call monitor_put( monit_south, qflx_y(:,:) )
574  end if
575  if ( monit_north > 0 ) then
576  qflx_y(:,:) = 0.0_rp
577  if ( .not. prc_has_n ) then
578  do iq = 1, qa
579  if ( tracer_advc(iq) .and. tracer_mass(iq) == 1.0_rp ) then
580  do i = is, ie
581  do k = ks, ke
582  qflx_y(k,i) = qflx_y(k,i) + qflx_rhoq(k,i,je,ydir,iq)
583  end do
584  end do
585  end if
586  end do
587  end if
588  call monitor_put( monit_north, qflx_y(:,:) )
589  end if
590 
591  return
592  end subroutine atmos_phy_tb_driver_calc_tendency
593 
594 end module mod_atmos_phy_tb_driver
mod_atmos_vars::momz_av
real(rp), dimension(:,:,:), pointer, public momz_av
Definition: mod_atmos_vars.F90:90
scale_statistics
module Statistics
Definition: scale_statistics.F90:11
scale_time::time_dtsec_atmos_phy_tb
real(dp), public time_dtsec_atmos_phy_tb
time interval of physics(turbulence ) [sec]
Definition: scale_time.F90:41
scale_atmos_grid_cartesc_index::ke
integer, public ke
end point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:52
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mu
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mu
Definition: mod_atmos_phy_sf_vars.F90:77
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
mod_atmos_admin::atmos_sw_phy_tb
logical, public atmos_sw_phy_tb
Definition: mod_atmos_admin.F90:57
scale_atmos_phy_tb_dns::atmos_phy_tb_dns
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)
Definition: scale_atmos_phy_tb_dns.F90:130
scale_tracer::qa
integer, public qa
Definition: scale_tracer.F90:34
scale_atmos_phy_tb_dns
module ATMOSPHERE / Physics Turbulence
Definition: scale_atmos_phy_tb_dns.F90:12
mod_atmos_vars::momx_av
real(rp), dimension(:,:,:), pointer, public momx_av
Definition: mod_atmos_vars.F90:91
scale_atmos_grid_cartesc::atmos_grid_cartesc_cdz
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cdz
z-length of control volume [m]
Definition: scale_atmos_grid_cartesC.F90:42
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvolzxv
real(rp), public atmos_grid_cartesc_real_totvolzxv
total volume (zxv, local) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:90
mod_atmos_vars::rhoq_tp
real(rp), dimension(:,:,:,:), allocatable, public rhoq_tp
Definition: mod_atmos_vars.F90:120
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_cz
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_cz
geopotential height [m] (zxy)
Definition: scale_atmos_grid_cartesC_real.F90:38
mod_atmos_vars::pott
real(rp), dimension(:,:,:), allocatable, target, public pott
Definition: mod_atmos_vars.F90:132
scale_tracer::tracer_mass
real(rp), dimension(qa_max), public tracer_mass
Definition: scale_tracer.F90:46
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_volwxy
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_volwxy
control volume (wxy) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:84
mod_atmos_phy_tb_vars::atmos_phy_tb_momz_t
real(rp), dimension(:,:,:), allocatable, public atmos_phy_tb_momz_t
Definition: mod_atmos_phy_tb_vars.F90:55
scale_prc_cartesc::prc_has_s
logical, public prc_has_s
Definition: scale_prc_cartesC.F90:50
mod_atmos_vars::qtrc_av
real(rp), dimension(:,:,:,:), pointer, public qtrc_av
Definition: mod_atmos_vars.F90:94
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
mod_atmos_admin
module ATMOS admin
Definition: mod_atmos_admin.F90:11
scale_atmos_grid_cartesc_index::jblock
integer, public jblock
block size for cache blocking: y
Definition: scale_atmos_grid_cartesC_index.F90:41
scale_tracer::tracer_advc
logical, dimension(qa_max), public tracer_advc
Definition: scale_tracer.F90:45
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvolwxy
real(rp), public atmos_grid_cartesc_real_totvolwxy
total volume (wxy, local) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:88
mod_atmos_phy_tb_vars::atmos_phy_tb_momx_t
real(rp), dimension(:,:,:), allocatable, public atmos_phy_tb_momx_t
Definition: mod_atmos_phy_tb_vars.F90:56
scale_atmos_grid_cartesc_metric
module Atmosphere Grid CartesianC metirc
Definition: scale_atmos_grid_cartesC_metric.F90:12
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:1594
mod_atmos_vars::rhot_av
real(rp), dimension(:,:,:), pointer, public rhot_av
Definition: mod_atmos_vars.F90:93
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:1653
scale_atmos_grid_cartesc::atmos_grid_cartesc_cdy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cdy
y-length of control volume [m]
Definition: scale_atmos_grid_cartesC.F90:61
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mv
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mv
Definition: mod_atmos_phy_sf_vars.F90:78
scale_atmos_phy_tb_dns::atmos_phy_tb_dns_setup
subroutine, public atmos_phy_tb_dns_setup(CDZ, CDX, CDY, CZ)
Definition: scale_atmos_phy_tb_dns.F90:88
scale_atmos_grid_cartesc::atmos_grid_cartesc_fdy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fdy
y-length of grid(j+1) to grid(j) [m]
Definition: scale_atmos_grid_cartesC.F90:63
scale_atmos_phy_tb_dns::atmos_phy_tb_dns_config
subroutine, public atmos_phy_tb_dns_config(TYPE_TB, I_TKE_out)
Config.
Definition: scale_atmos_phy_tb_dns.F90:64
mod_atmos_phy_sf_vars
module ATMOSPHERIC Surface Variables
Definition: mod_atmos_phy_sf_vars.F90:12
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_mw
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_mw
Definition: mod_atmos_phy_sf_vars.F90:76
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:80
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:44
mod_atmos_vars::rhot
real(rp), dimension(:,:,:), allocatable, target, public rhot
Definition: mod_atmos_vars.F90:79
scale_atmos_grid_cartesc_real
module Atmosphere GRID CartesC Real(real space)
Definition: scale_atmos_grid_cartesC_real.F90:11
mod_atmos_vars::qtrc
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc
Definition: mod_atmos_vars.F90:80
scale_file_history
module file_history
Definition: scale_file_history.F90:15
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_volzuy
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_volzuy
control volume (zuy) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:85
scale_prc_cartesc::prc_has_n
logical, public prc_has_n
Definition: scale_prc_cartesC.F90:48
scale_prc_cartesc::prc_has_e
logical, public prc_has_e
Definition: scale_prc_cartesC.F90:49
mod_atmos_phy_tb_driver::atmos_phy_tb_driver_calc_tendency
subroutine, public atmos_phy_tb_driver_calc_tendency(update_flag)
calclate tendency
Definition: mod_atmos_phy_tb_driver.F90:172
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:172
scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_mapf
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_mapf
map factor
Definition: scale_atmos_grid_cartesC_metric.F90:34
scale_atmos_grid_cartesc_index::i_xy
integer, public i_xy
Definition: scale_atmos_grid_cartesC_index.F90:98
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
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_io
module STDIO
Definition: scale_io.F90:10
mod_atmos_vars::dens
real(rp), dimension(:,:,:), allocatable, target, public dens
Definition: mod_atmos_vars.F90:75
scale_tracer::k
real(rp), public k
Definition: scale_tracer.F90:44
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_volzxv
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_volzxv
control volume (zxv) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:86
scale_atmos_grid_cartesc_index
module atmosphere / grid / cartesC index
Definition: scale_atmos_grid_cartesC_index.F90:12
scale_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_vol
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_vol
control volume (zxy) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:83
mod_atmos_vars::momz
real(rp), dimension(:,:,:), allocatable, target, public momz
Definition: mod_atmos_vars.F90:76
scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_gsqrt
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_gsqrt
transformation metrics from Z to Xi, {G}^1/2
Definition: scale_atmos_grid_cartesC_metric.F90:37
scale_atmos_grid_cartesc_index::zdir
integer, parameter, public zdir
Definition: scale_atmos_grid_cartesC_index.F90:32
mod_atmos_phy_tb_vars
module Atmosphere / Physics Turbulence
Definition: mod_atmos_phy_tb_vars.F90:12
mod_atmos_vars::momy_av
real(rp), dimension(:,:,:), pointer, public momy_av
Definition: mod_atmos_vars.F90:92
scale_prc_cartesc
module process / cartesC
Definition: scale_prc_cartesC.F90:11
mod_atmos_phy_tb_driver::atmos_phy_tb_driver_tracer_setup
subroutine, public atmos_phy_tb_driver_tracer_setup
Tracer setup.
Definition: mod_atmos_phy_tb_driver.F90:55
scale_tracer::tracer_name
character(len=h_short), dimension(qa_max), public tracer_name
Definition: scale_tracer.F90:38
scale_prof
module profiler
Definition: scale_prof.F90:11
mod_atmos_vars::momz_tp
real(rp), dimension(:,:,:), allocatable, public momz_tp
Definition: mod_atmos_vars.F90:115
mod_atmos_vars::momx
real(rp), dimension(:,:,:), allocatable, target, public momx
Definition: mod_atmos_vars.F90:77
scale_atmos_grid_cartesc_index::is
integer, public is
start point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:53
mod_atmos_phy_tb_driver
module ATMOSPHERE / Physics Turbulence
Definition: mod_atmos_phy_tb_driver.F90:12
scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_j23g
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_j23g
(2,3) element of Jacobian matrix * {G}^1/2
Definition: scale_atmos_grid_cartesC_metric.F90:39
scale_monitor::monitor_reg
subroutine, public monitor_reg(name, desc, unit, itemid, ndims, dim_type, is_tendency)
Search existing item, or matching check between requested and registered item.
Definition: scale_monitor.F90:241
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvolzuy
real(rp), public atmos_grid_cartesc_real_totvolzuy
total volume (zuy, local) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:89
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
mod_atmos_vars::momy
real(rp), dimension(:,:,:), allocatable, target, public momy
Definition: mod_atmos_vars.F90:78
mod_atmos_phy_tb_driver::atmos_phy_tb_driver_setup
subroutine, public atmos_phy_tb_driver_setup
Setup.
Definition: mod_atmos_phy_tb_driver.F90:96
scale_time
module TIME
Definition: scale_time.F90:11
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_qtrc
real(rp), dimension(:,:,:), allocatable, target, public atmos_phy_sf_sflx_qtrc
Definition: mod_atmos_phy_sf_vars.F90:84
scale_atmos_grid_cartesc::atmos_grid_cartesc_fdx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fdx
x-length of grid(i+1) to grid(i) [m]
Definition: scale_atmos_grid_cartesC.F90:62
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:1713
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:114
mod_atmos_vars::dens_av
real(rp), dimension(:,:,:), pointer, public dens_av
Definition: mod_atmos_vars.F90:89
scale_atmos_grid_cartesc_index::ks
integer, public ks
start point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:51
mod_atmos_phy_tb_vars::i_tke
integer, public i_tke
Definition: mod_atmos_phy_tb_vars.F90:61
scale_statistics::statistics_checktotal
logical, public statistics_checktotal
calc&report variable totals to logfile?
Definition: scale_statistics.F90:64
mod_atmos_phy_tb_vars::atmos_phy_tb_momy_t
real(rp), dimension(:,:,:), allocatable, public atmos_phy_tb_momy_t
Definition: mod_atmos_phy_tb_vars.F90:57
scale_atmos_phy_tb_common
module ATMOSPHERE / Physics Turbulence
Definition: scale_atmos_phy_tb_common.F90:12
mod_atmos_vars
module ATMOSPHERIC Variables
Definition: mod_atmos_vars.F90:12
mod_atmos_vars::rhot_tp
real(rp), dimension(:,:,:), allocatable, public rhot_tp
Definition: mod_atmos_vars.F90:118
scale_atmos_phy_tb_smg::atmos_phy_tb_smg_setup
subroutine, public atmos_phy_tb_smg_setup(FZ, CZ, CDX, CDY, MAPF, horizontal)
Setup.
Definition: scale_atmos_phy_tb_smg.F90:100
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
mod_atmos_vars::momy_tp
real(rp), dimension(:,:,:), allocatable, public momy_tp
Definition: mod_atmos_vars.F90:124
mod_atmos_phy_sf_vars::atmos_phy_sf_sflx_sh
real(rp), dimension(:,:), allocatable, public atmos_phy_sf_sflx_sh
Definition: mod_atmos_phy_sf_vars.F90:79
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:45
mod_atmos_admin::atmos_sw_phy_bl
logical, public atmos_sw_phy_bl
Definition: mod_atmos_admin.F90:58
mod_atmos_vars::momx_tp
real(rp), dimension(:,:,:), allocatable, public momx_tp
Definition: mod_atmos_vars.F90:123
scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_j33g
real(rp), public atmos_grid_cartesc_metric_j33g
(3,3) element of Jacobian matrix * {G}^1/2
Definition: scale_atmos_grid_cartesC_metric.F90:40
scale_atmos_grid_cartesc::atmos_grid_cartesc_cdx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cdx
x-length of control volume [m]
Definition: scale_atmos_grid_cartesC.F90:60
scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_j13g
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_j13g
(1,3) element of Jacobian matrix * {G}^1/2
Definition: scale_atmos_grid_cartesC_metric.F90:38
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_fz
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_fz
geopotential height [m] (wxy)
Definition: scale_atmos_grid_cartesC_real.F90:42
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_totvol
real(rp), public atmos_grid_cartesc_real_totvol
total volume (zxy, local) [m3]
Definition: scale_atmos_grid_cartesC_real.F90:87
scale_atmos_grid_cartesc
module atmosphere / grid / cartesC
Definition: scale_atmos_grid_cartesC.F90:12
mod_atmos_phy_tb_vars::atmos_phy_tb_rhoq_t
real(rp), dimension(:,:,:,:), allocatable, public atmos_phy_tb_rhoq_t
Definition: mod_atmos_phy_tb_vars.F90:59
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: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:1534
scale_atmos_phy_tb_smg
module ATMOSPHERE / Physics Turbulence
Definition: scale_atmos_phy_tb_smg.F90:21
scale_prc_cartesc::prc_has_w
logical, public prc_has_w
Definition: scale_prc_cartesC.F90:47
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56
mod_atmos_phy_tb_vars::atmos_phy_tb_rhot_t
real(rp), dimension(:,:,:), allocatable, public atmos_phy_tb_rhot_t
Definition: mod_atmos_phy_tb_vars.F90:58
mod_atmos_admin::atmos_phy_tb_type
character(len=h_short), public atmos_phy_tb_type
Definition: mod_atmos_admin.F90:41
scale_atmos_phy_tb_smg::atmos_phy_tb_smg
subroutine, public atmos_phy_tb_smg(qflx_sgs_momz, qflx_sgs_momx, qflx_sgs_momy, qflx_sgs_rhot, qflx_sgs_rhoq, MOMZ_t, MOMX_t, MOMY_t, RHOT_t, RHOQ_t, Nu, Ri, Pr, MOMZ, MOMX, MOMY, POTT, DENS, QTRC, N2, FZ, FDZ, RCDZ, RFDZ, CDX, FDX, CDY, FDY, GSQRT, J13G, J23G, J33G, MAPF, dt)
Definition: scale_atmos_phy_tb_smg.F90:242
scale_monitor
module MONITOR
Definition: scale_monitor.F90:12