SCALE-RM
scale_atmos_dyn_tstep_short_fvm_hevi.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
11 
12 #ifdef PROFILE_FAPP
13 #define PROFILE_START(name) call fapp_start(name, 1, 1)
14 #define PROFILE_STOP(name) call fapp_stop (name, 1, 1)
15 #else
16 #define PROFILE_START(name)
17 #define PROFILE_STOP(name)
18 #endif
19 
20 #include "scalelib.h"
22  !-----------------------------------------------------------------------------
23  !
24  !++ used modules
25  !
26  use scale_precision
27  use scale_io
28  use scale_prof
30  use scale_index
31  use scale_tracer
32 #if defined DEBUG || defined QUICKDEBUG
33  use scale_debug, only: &
34  check
35  use scale_const, only: &
36  undef => const_undef, &
37  iundef => const_undef2
38 #endif
39  !-----------------------------------------------------------------------------
40  implicit none
41  private
42  !-----------------------------------------------------------------------------
43  !
44  !++ Public procedure
45  !
49 
50  !-----------------------------------------------------------------------------
51  !
52  !++ Public parameters & variables
53  !
54  !-----------------------------------------------------------------------------
55  !
56  !++ Private procedure
57  !
58 #if 1
59 #define F2H(k,p,idx) (CDZ(k+p-1)*GSQRT(k+p-1,i,j,idx)/(CDZ(k)*GSQRT(k,i,j,idx)+CDZ(k+1)*GSQRT(k+1,i,j,idx)))
60 #else
61 # define F2H(k,p,idx) 0.5_RP
62 #endif
63 
64  !-----------------------------------------------------------------------------
65  !
66  !++ Private parameters & variables
67  !
68  integer, private, parameter :: NB = 1
69  integer, private, parameter :: VA_FVM_HEVI = 0
70  integer :: IFS_OFF
71  integer :: JFS_OFF
72 
73  !-----------------------------------------------------------------------------
74 contains
75  !-----------------------------------------------------------------------------
78  ATMOS_DYN_TYPE, &
79  VA_out, &
80  VAR_NAME, &
81  VAR_DESC, &
82  VAR_UNIT )
83  use scale_prc, only: &
84  prc_abort
85  implicit none
86 
87  character(len=*), intent(in) :: atmos_dyn_type
88  integer, intent(out) :: va_out
89  character(len=H_SHORT), intent(out) :: var_name(:)
90  character(len=H_MID), intent(out) :: var_desc(:)
91  character(len=H_SHORT), intent(out) :: var_unit(:)
92  !---------------------------------------------------------------------------
93 
94  if ( atmos_dyn_type /= 'FVM-HEVI' .AND. atmos_dyn_type /= 'HEVI' ) then
95  log_error("ATMOS_DYN_Tstep_short_fvm_hevi_regist",*) 'ATMOS_DYN_TYPE is not FVM-HEVI. Check!'
96  call prc_abort
97  endif
98 
99  va_out = va_fvm_hevi
100  var_name(:) = ""
101  var_desc(:) = ""
102  var_unit(:) = ""
103 
104  log_newline
105  log_info("ATMOS_DYN_Tstep_short_fvm_hevi_regist",*) 'Register additional prognostic variables (HEVI)'
106  if ( va_out < 1 ) then
107  log_info_cont(*) '=> nothing.'
108  endif
109 
110  return
112 
113  !-----------------------------------------------------------------------------
116  implicit none
117  !---------------------------------------------------------------------------
118 
119  log_info("ATMOS_DYN_Tstep_short_fvm_hevi_setup",*) 'HEVI Setup'
120 
121  return
123 
124  !-----------------------------------------------------------------------------
125  subroutine atmos_dyn_tstep_short_fvm_hevi( &
126  DENS_RK, MOMZ_RK, MOMX_RK, MOMY_RK, RHOT_RK, &
127  PROG_RK, &
128  mflx_hi, tflx_hi, &
129  DENS0, MOMZ0, MOMX0, MOMY0, RHOT0, &
130  DENS, MOMZ, MOMX, MOMY, RHOT, &
131  DENS_t, MOMZ_t, MOMX_t, MOMY_t, RHOT_t, &
132  PROG0, PROG, &
133  DPRES0, RT2P, CORIOLI, &
134  num_diff, wdamp_coef, divdmp_coef, DDIV, &
135  FLAG_FCT_MOMENTUM, FLAG_FCT_T, &
136  FLAG_FCT_ALONG_STREAM, &
137  CDZ, FDZ, FDX, FDY, &
138  RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, &
139  PHI, GSQRT, J13G, J23G, J33G, MAPF, &
140  REF_dens, REF_rhot, &
141  BND_W, BND_E, BND_S, BND_N, TwoD, &
142  dtrk, last )
144  use scale_const, only: &
145 #ifdef dry
146  rdry => const_rdry, &
147  cvdry => const_cvdry, &
148  cpdry => const_cpdry, &
149 #endif
150  eps => const_eps, &
151  grav => const_grav, &
152  p00 => const_pre00
153  use scale_atmos_dyn_fvm_fct, only: &
155  use scale_atmos_dyn_fvm_flux, only: &
175 #ifdef HIST_TEND
176  use scale_file_history, only: &
177  file_history_in
178 #endif
179  use scale_matrix, only: &
180  matrix_solver_tridiagonal, &
182  implicit none
183 
184  real(rp), intent(out) :: dens_rk(ka,ia,ja) ! prognostic variables
185  real(rp), intent(out) :: momz_rk(ka,ia,ja) !
186  real(rp), intent(out) :: momx_rk(ka,ia,ja) !
187  real(rp), intent(out) :: momy_rk(ka,ia,ja) !
188  real(rp), intent(out) :: rhot_rk(ka,ia,ja) !
189 
190  real(rp), intent(out) :: prog_rk(ka,ia,ja,va) !
191 
192  real(rp), intent(inout) :: mflx_hi(ka,ia,ja,3) ! rho * vel(x,y,z)
193  real(rp), intent(out) :: tflx_hi(ka,ia,ja,3) ! rho * theta * vel(x,y,z)
194 
195  real(rp), intent(in),target :: dens0(ka,ia,ja) ! prognostic variables
196  real(rp), intent(in),target :: momz0(ka,ia,ja) ! at previous dynamical time step
197  real(rp), intent(in),target :: momx0(ka,ia,ja) !
198  real(rp), intent(in),target :: momy0(ka,ia,ja) !
199  real(rp), intent(in),target :: rhot0(ka,ia,ja) !
200 
201  real(rp), intent(in) :: dens(ka,ia,ja) ! prognostic variables
202  real(rp), intent(in) :: momz(ka,ia,ja) ! at previous RK step
203  real(rp), intent(in) :: momx(ka,ia,ja) !
204  real(rp), intent(in) :: momy(ka,ia,ja) !
205  real(rp), intent(in) :: rhot(ka,ia,ja) !
206 
207  real(rp), intent(in) :: dens_t(ka,ia,ja)
208  real(rp), intent(in) :: momz_t(ka,ia,ja)
209  real(rp), intent(in) :: momx_t(ka,ia,ja)
210  real(rp), intent(in) :: momy_t(ka,ia,ja)
211  real(rp), intent(in) :: rhot_t(ka,ia,ja)
212 
213  real(rp), intent(in) :: prog0(ka,ia,ja,va)
214  real(rp), intent(in) :: prog (ka,ia,ja,va)
215 
216  real(rp), intent(in) :: dpres0(ka,ia,ja)
217  real(rp), intent(in) :: rt2p(ka,ia,ja)
218  real(rp), intent(in) :: corioli(ia,ja)
219  real(rp), intent(in) :: num_diff(ka,ia,ja,5,3)
220  real(rp), intent(in) :: wdamp_coef(ka)
221  real(rp), intent(in) :: divdmp_coef
222  real(rp), intent(in) :: ddiv(ka,ia,ja)
223 
224  logical, intent(in) :: flag_fct_momentum
225  logical, intent(in) :: flag_fct_t
226  logical, intent(in) :: flag_fct_along_stream
227 
228  real(rp), intent(in) :: cdz(ka)
229  real(rp), intent(in) :: fdz(ka-1)
230  real(rp), intent(in) :: fdx(ia-1)
231  real(rp), intent(in) :: fdy(ja-1)
232  real(rp), intent(in) :: rcdz(ka)
233  real(rp), intent(in) :: rcdx(ia)
234  real(rp), intent(in) :: rcdy(ja)
235  real(rp), intent(in) :: rfdz(ka-1)
236  real(rp), intent(in) :: rfdx(ia-1)
237  real(rp), intent(in) :: rfdy(ja-1)
238 
239  real(rp), intent(in) :: phi (ka,ia,ja)
240  real(rp), intent(in) :: gsqrt (ka,ia,ja,7)
241  real(rp), intent(in) :: j13g (ka,ia,ja,7)
242  real(rp), intent(in) :: j23g (ka,ia,ja,7)
243  real(rp), intent(in) :: j33g
244  real(rp), intent(in) :: mapf (ia,ja,2,4)
245  real(rp), intent(in) :: ref_dens(ka,ia,ja)
246  real(rp), intent(in) :: ref_rhot(ka,ia,ja)
247 
248  logical, intent(in) :: bnd_w
249  logical, intent(in) :: bnd_e
250  logical, intent(in) :: bnd_s
251  logical, intent(in) :: bnd_n
252  logical, intent(in) :: twod
253 
254  real(rp), intent(in) :: dtrk
255  logical, intent(in) :: last
256 
257 
258  ! diagnostic variables (work space)
259  real(rp) :: pott(ka,ia,ja) ! potential temperature [K]
260  real(rp) :: dpres(ka,ia,ja) ! pressure deviation from reference pressure
261 
262  real(rp) :: qflx_hi (ka,ia,ja,3)
263  real(rp) :: qflx_j13(ka,ia,ja)
264  real(rp) :: qflx_j23(ka,ia,ja)
265 
266  real(rp) :: advch ! horizontal advection
267  real(rp) :: advcv ! vertical advection
268  real(rp) :: wdmp ! Raileight damping
269  real(rp) :: div ! divergence damping
270  real(rp) :: pg ! pressure gradient force
271  real(rp) :: cf ! colioris force
272 #ifdef HIST_TEND
273  real(rp) :: advch_t(ka,ia,ja,5)
274  real(rp) :: advcv_t(ka,ia,ja,5)
275  real(rp) :: wdmp_t(ka,ia,ja)
276  real(rp) :: ddiv_t(ka,ia,ja,3)
277  real(rp) :: pg_t(ka,ia,ja,3)
278  real(rp) :: cf_t(ka,ia,ja,2)
279  logical :: lhist
280 #endif
281 
282  ! for implicit solver
283  real(rp) :: a(ka)
284  real(rp) :: b
285  real(rp) :: sr(ka,ia,ja)
286  real(rp) :: sw(ka,ia,ja)
287  real(rp) :: st(ka,ia,ja)
288 
289  real(rp) :: pt(ka,lsize)
290  real(rp) :: ci(ks:ke-1,lsize)
291  real(rp) :: co(ks:ke-1,lsize)
292  real(rp) :: f1(ks:ke-1,lsize)
293  real(rp) :: f2(ks:ke-1,lsize)
294  real(rp) :: f3(ks:ke-1,lsize)
295 
296 #ifdef _OPENACC
297  real(rp) :: work(kmax-1,4) ! for CR
298 #endif
299 
300  integer :: iis, iie, jjs, jje
301  integer :: k, i, j, ii
302  integer :: iss, iee
303 #if LSIZE == 1
304  integer, parameter :: l = 1
305 #else
306  integer :: l
307 #endif
308 
309 
310 #ifdef DEBUG
311  pott(:,:,:) = undef
312  dpres(:,:,:) = undef
313 
314  qflx_hi(:,:,:,:) = undef
315  qflx_j13(:,:,:) = undef
316  qflx_j23(:,:,:) = undef
317 #endif
318 
319 #if defined DEBUG || defined QUICKDEBUG
320  dens_rk( 1:ks-1,:,:) = undef
321  dens_rk(ke+1:ka ,:,:) = undef
322  momz_rk( 1:ks-1,:,:) = undef
323  momz_rk(ke+1:ka ,:,:) = undef
324  momx_rk( 1:ks-1,:,:) = undef
325  momx_rk(ke+1:ka ,:,:) = undef
326  momy_rk( 1:ks-1,:,:) = undef
327  momy_rk(ke+1:ka ,:,:) = undef
328  rhot_rk( 1:ks-1,:,:) = undef
329  rhot_rk(ke+1:ka ,:,:) = undef
330  prog_rk( 1:ks-1,:,:,:) = undef
331  prog_rk(ke+1:ka ,:,:,:) = undef
332 #endif
333 
334 #ifdef HIST_TEND
335  lhist = last
336 #endif
337 
338 #ifdef PROFILE_FIPP
339  call fipp_start()
340 #endif
341 
342  ifs_off = 1
343  jfs_off = 1
344  if ( bnd_w ) ifs_off = 0
345  if ( bnd_s ) jfs_off = 0
346 
347 
348  !$acc data &
349  !$acc copy(mflx_hi) &
350  !$acc copyout(DENS_RK,MOMZ_RK,MOMX_RK,MOMY_RK,RHOT_RK,PROG_RK, &
351  !$acc tflx_hi) &
352 #ifdef HIST_TEND
353  !$acc copyout(advch_t,advcv_t,wdmp_t,ddiv_t,pg_t,cf_t) &
354 #endif
355  !$acc copyin(DENS0,MOMZ0,MOMX0,MOMY0,RHOT0, &
356  !$acc DENS,MOMZ,MOMX,MOMY,RHOT,DENS_t,MOMZ_t,MOMX_t,MOMY_t,RHOT_t, &
357  !$acc PROG0,PROG, &
358  !$acc DPRES0,RT2P,CORIOLI,num_diff,wdamp_coef,divdmp_coef,DDIV, &
359  !$acc FLAG_FCT_MOMENTUM,FLAG_FCT_T,FLAG_FCT_ALONG_STREAM, &
360  !$acc CDZ,FDZ,FDX,FDY,RCDZ,RCDX,RCDY,RFDZ,RFDX,RFDY, &
361  !$acc PHI,GSQRT,J13G,J23G,J33G,MAPF,REF_dens,REF_rhot, &
362  !$acc BND_W,BND_E,BND_S,BND_N,TwoD,dtrk,last) &
363  !$acc create(POTT,DPRES, &
364  !$acc qflx_hi,qflx_J13,qflx_J23, &
365  !$acc Sr,Sw,St)
366 
367  do jjs = js, je, jblock
368  jje = jjs+jblock-1
369  do iis = is, ie, iblock
370  iie = iis+iblock-1
371 
372  profile_start("hevi_pres")
373  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
374  !$omp shared(JJS,JJE,IIS,IIE,IA,KS,KE,DPRES0,RT2P,RHOT,REF_rhot,DPRES,DENS,PHI)
375  !$acc kernels
376  do j = jjs, jje+1
377  do i = iis, min(iie+1,ia)
378  do k = ks, ke
379 #ifdef DEBUG
380  call check( __line__, dpres0(k,i,j) )
381  call check( __line__, rt2p(k,i,j) )
382  call check( __line__, rhot(k,i,j) )
383  call check( __line__, ref_rhot(k,i,j) )
384 #endif
385  dpres(k,i,j) = dpres0(k,i,j) + rt2p(k,i,j) * ( rhot(k,i,j) - ref_rhot(k,i,j) )
386  enddo
387  dpres(ks-1,i,j) = dpres0(ks-1,i,j) - dens(ks,i,j) * ( phi(ks-1,i,j) - phi(ks+1,i,j) )
388  dpres(ke+1,i,j) = dpres0(ke+1,i,j) - dens(ke,i,j) * ( phi(ke+1,i,j) - phi(ke-1,i,j) )
389  enddo
390  enddo
391  !$acc end kernels
392  profile_stop("hevi_pres")
393 
394  !##### continuity equation #####
395 
396  profile_start("hevi_mflx_z")
397  ! at (x, y, w)
398  if ( twod ) then
399  !$omp parallel do default(none) private(j,k) OMP_SCHEDULE_ &
400  !$omp shared(JJS,JJE,IS,KS,KE,GSQRT,I_XYW,MOMY,J23G,mflx_hi,MAPF,I_XY,num_diff)
401  !$acc kernels
402  do j = jjs, jje
403  mflx_hi(ks-1,is,j,zdir) = 0.0_rp
404  do k = ks, ke-1
405 #ifdef DEBUG
406  call check( __line__, momy(k+1,is,j) )
407  call check( __line__, momy(k+1,is,j-1) )
408  call check( __line__, momy(k ,is,j) )
409  call check( __line__, momy(k ,is,j-1) )
410 #endif
411  mflx_hi(k,is,j,zdir) = j23g(k,is,j,i_xyw) * 0.25_rp &
412  * ( momy(k+1,is,j)+momy(k+1,is,j-1) &
413  + momy(k ,is,j)+momy(k ,is,j-1) ) &
414  + gsqrt(k,is,j,i_xyw) / mapf(is,j,2,i_xy) * num_diff(k,is,j,i_dens,zdir)
415  enddo
416  mflx_hi(ke,is,j,zdir) = 0.0_rp
417  enddo
418  !$acc end kernels
419  else
420  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
421  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,GSQRT,I_XYW,MOMX,MOMY,J13G,J23G,mflx_hi,MAPF,I_XY,num_diff)
422  !$acc kernels
423  do j = jjs, jje
424  do i = iis-1, iie
425  mflx_hi(ks-1,i,j,zdir) = 0.0_rp
426  do k = ks, ke-1
427 #ifdef DEBUG
428  call check( __line__, momx(k+1,i ,j) )
429  call check( __line__, momx(k+1,i-1,j) )
430  call check( __line__, momx(k ,i ,j) )
431  call check( __line__, momx(k ,i+1,j) )
432  call check( __line__, momy(k+1,i,j) )
433  call check( __line__, momy(k+1,i,j-1) )
434  call check( __line__, momy(k ,i,j) )
435  call check( __line__, momy(k ,i,j-1) )
436 #endif
437  mflx_hi(k,i,j,zdir) = j13g(k,i,j,i_xyw) * 0.25_rp / mapf(i,j,2,i_xy) &
438  * ( momx(k+1,i,j)+momx(k+1,i-1,j) &
439  + momx(k ,i,j)+momx(k ,i-1,j) ) &
440  + j23g(k,i,j,i_xyw) * 0.25_rp / mapf(i,j,1,i_xy) &
441  * ( momy(k+1,i,j)+momy(k+1,i,j-1) &
442  + momy(k ,i,j)+momy(k ,i,j-1) ) &
443  + gsqrt(k,i,j,i_xyw) / ( mapf(i,j,1,i_xy)*mapf(i,j,2,i_xy) ) * num_diff(k,i,j,i_dens,zdir)
444  enddo
445  mflx_hi(ke ,i,j,zdir) = 0.0_rp
446  enddo
447  enddo
448  !$acc end kernels
449  end if
450 #ifdef DEBUG
451  k = iundef; i = iundef; j = iundef
452 #endif
453  profile_stop("hevi_mflx_z")
454 
455  if ( .not. twod ) then
456  profile_start("hevi_mflx_x")
457  iss = max(iis-1,is-ifs_off)
458  iee = min(iie,ieh)
459  ! at (u, y, z)
460  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
461  !$omp shared(JJS,JJE,iss,iee,KS,KE,GSQRT,I_UYZ,MOMX,num_diff,mflx_hi,MAPF,I_UY)
462  !$acc kernels
463  do j = jjs, jje
464  do i = iss, iee
465  do k = ks, ke
466 #ifdef DEBUG
467  call check( __line__, gsqrt(k,i,j,i_uyz) )
468  call check( __line__, momx(k,i,j) )
469  call check( __line__, num_diff(k,i,j,i_dens,xdir) )
470 #endif
471  mflx_hi(k,i,j,xdir) = gsqrt(k,i,j,i_uyz) / mapf(i,j,2,i_uy) &
472  * ( momx(k,i,j) + num_diff(k,i,j,i_dens,xdir) )
473  enddo
474  enddo
475  enddo
476  !$acc end kernels
477 #ifdef DEBUG
478  k = iundef; i = iundef; j = iundef
479 #endif
480  profile_stop("hevi_mflx_x")
481  end if
482 
483  ! at (x, v, z)
484  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
485  !$omp shared(JJS,JS,JFS_OFF,JJE,JEH,IIS,IIE,KS,KE,GSQRT,I_XVZ,MOMY,num_diff,mflx_hi,MAPF,I_XV)
486  !$acc kernels
487  do j = max(jjs-1,js-jfs_off), min(jje,jeh)
488  do i = iis, iie
489  do k = ks, ke
490 #ifdef DEBUG
491  call check( __line__, gsqrt(k,i,j,i_xvz) )
492  call check( __line__, momy(k,i,j) )
493  call check( __line__, num_diff(k,i,j,i_dens,ydir) )
494 #endif
495  mflx_hi(k,i,j,ydir) = gsqrt(k,i,j,i_xvz) / mapf(i,j,1,i_xv) &
496  * ( momy(k,i,j) + num_diff(k,i,j,i_dens,ydir) )
497  enddo
498  enddo
499  enddo
500  !$acc end kernels
501 #ifdef DEBUG
502  k = iundef; i = iundef; j = iundef
503 #endif
504 
505  !--- update density
506  profile_start("hevi_sr")
507  if ( twod ) then
508  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
509  !$omp private(j,k,advcv,advch) &
510 #ifdef HIST_TEND
511  !$omp shared(advch_t,lhist) &
512 #endif
513  !$omp shared(JJS,JJE,IS,KS,KE) &
514  !$omp shared(DENS0,Sr,mflx_hi,DENS_t) &
515  !$omp shared(RCDZ,RCDY,MAPF,GSQRT,I_XY,I_XYZ)
516  !$acc kernels
517  do j = jjs, jje
518  do k = ks, ke
519 #ifdef DEBUG
520  call check( __line__, dens0(k,is,j) )
521  call check( __line__, mflx_hi(k ,is,j ,ydir) )
522  call check( __line__, mflx_hi(k ,is,j-1,ydir) )
523  call check( __line__, dens_t(k,is,j) )
524 #endif
525  advcv = - ( mflx_hi(k,is,j,zdir)-mflx_hi(k-1,is,j, zdir) ) * rcdz(k)
526  advch = - ( mflx_hi(k,is,j,ydir)-mflx_hi(k ,is,j-1,ydir) ) * rcdy(j)
527  sr(k,is,j) = ( advcv + advch ) * mapf(is,j,2,i_xy) / gsqrt(k,is,j,i_xyz) + dens_t(k,is,j)
528 #ifdef HIST_TEND
529  if ( lhist ) then
530  advch_t(k,is,j,i_dens) = ( advch + advcv ) * mapf(is,j,2,i_xy) / gsqrt(k,is,j,i_xyz)
531  end if
532 #endif
533  enddo
534  enddo
535  !$acc end kernels
536  else
537  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
538  !$omp private(i,j,k,advcv,advch) &
539 #ifdef HIST_TEND
540  !$omp shared(advch_t,lhist) &
541 #endif
542  !$omp shared(JJS,JJE,IIS,IIE,KS,KE) &
543  !$omp shared(DENS0,Sr,mflx_hi,DENS_t) &
544  !$omp shared(RCDZ,RCDX,RCDY,MAPF,GSQRT,I_XY,I_XYZ)
545  !$acc kernels
546  do j = jjs, jje
547  do i = iis, iie
548  do k = ks, ke
549 #ifdef DEBUG
550  call check( __line__, dens0(k,i,j) )
551  call check( __line__, mflx_hi(k ,i ,j ,xdir) )
552  call check( __line__, mflx_hi(k ,i-1,j ,xdir) )
553  call check( __line__, mflx_hi(k ,i ,j ,ydir) )
554  call check( __line__, mflx_hi(k ,i ,j-1,ydir) )
555  call check( __line__, dens_t(k,i,j) )
556 #endif
557  advcv = - ( mflx_hi(k,i,j,zdir)-mflx_hi(k-1,i ,j, zdir) ) * rcdz(k)
558  advch = - ( ( mflx_hi(k,i,j,xdir)-mflx_hi(k ,i-1,j, xdir) ) * rcdx(i) &
559  + ( mflx_hi(k,i,j,ydir)-mflx_hi(k ,i, j-1,ydir) ) * rcdy(j) )
560  sr(k,i,j) = ( advcv + advch ) * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz) + dens_t(k,i,j)
561 #ifdef HIST_TEND
562  if ( lhist ) then
563  advch_t(k,i,j,i_dens) = ( advch + advcv ) * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz)
564  end if
565 #endif
566  enddo
567  enddo
568  enddo
569  !$acc end kernels
570  end if
571 #ifdef DEBUG
572  k = iundef; i = iundef; j = iundef
573 #endif
574  profile_stop("hevi_sr")
575 
576  !##### momentum equation (z) #####
577 
578  ! at (x, y, z)
579  ! not that z-index is added by -1
580  profile_start("hevi_momz_qflxhi_z")
581  call atmos_dyn_fvm_fluxz_xyw( qflx_hi(:,:,:,zdir), & ! (out)
582  momz, momz, dens, & ! (in)
583  gsqrt(:,:,:,i_xyz), j33g, & ! (in)
584  num_diff(:,:,:,i_momz,zdir), & ! (in)
585  cdz, fdz, dtrk, &
586  iis, iie, jjs, jje ) ! (in)
587  profile_stop("hevi_momz_qflxhi_z")
588 
589  profile_start("hevi_momz_qflxj")
590  if ( .not. twod ) &
591  call atmos_dyn_fvm_fluxj13_xyw( qflx_j13, & ! (out)
592  momx, momz, dens, & ! (in)
593  gsqrt(:,:,:,i_xyz), j13g(:,:,:,i_xyz), mapf(:,:,:,i_xy), & ! (in)
594  cdz, twod, &
595  iis, iie, jjs, jje ) ! (in)
596  call atmos_dyn_fvm_fluxj23_xyw( qflx_j23, & ! (out)
597  momy, momz, dens, & ! (in)
598  gsqrt(:,:,:,i_xyz), j23g(:,:,:,i_xyz), mapf(:,:,:,i_xy), & ! (in)
599  cdz, twod, &
600  iis, iie, jjs, jje ) ! (in)
601  profile_stop("hevi_momz_qflxj")
602 
603  ! at (u, y, w)
604  if ( .not. twod ) then
605  profile_start("hevi_momz_qflxhi_x")
606  call atmos_dyn_fvm_fluxx_xyw( qflx_hi(:,:,:,xdir), & ! (out)
607  momx, momz, dens, & ! (in)
608  gsqrt(:,:,:,i_uyw), mapf(:,:,:,i_uy), & ! (in)
609  num_diff(:,:,:,i_momz,xdir), & ! (in)
610  cdz, twod, & ! (in)
611  iis, iie, jjs, jje ) ! (in)
612  profile_stop("hevi_momz_qflxhi_x")
613  end if
614 
615  ! at (x, v, w)
616  profile_start("hevi_momz_qflxhi_y")
617  call atmos_dyn_fvm_fluxy_xyw( qflx_hi(:,:,:,ydir), & ! (out)
618  momy, momz, dens, & ! (in)
619  gsqrt(:,:,:,i_xvw), mapf(:,:,:,i_xv), & ! (in)
620  num_diff(:,:,:,i_momz,ydir), & ! (in)
621  cdz, twod, & ! (in)
622  iis, iie, jjs, jje ) ! (in)
623  profile_stop("hevi_momz_qflxhi_y")
624 
625  !--- update momentum(z)
626  profile_start("hevi_sw")
627  if ( twod ) then
628  !$omp parallel do default(none) OMP_SCHEDULE_ &
629  !$omp private(j,k,advcv,advch,cf,wdmp,div) &
630 #ifdef HIST_TEND
631  !$omp shared(lhist,advcv_t,advch_t,wdmp_t,ddiv_t) &
632 #endif
633  !$omp shared(JJS,JJE,IS,KS,KE) &
634  !$omp shared(qflx_hi,qflx_J23,DDIV,MOMZ0,MOMZ_t,Sw) &
635  !$omp shared(RFDZ,RCDY,FDZ,dtrk,wdamp_coef,divdmp_coef) &
636  !$omp shared(MAPF,GSQRT,I_XY,I_XYW)
637  !$acc kernels
638  do j = jjs, jje
639 !OCL NORECURRENCE
640  do k = ks, ke-1
641 #ifdef DEBUG
642  call check( __line__, qflx_hi(k ,is,j ,zdir) )
643  call check( __line__, qflx_hi(k-1,is,j ,zdir) )
644  call check( __line__, qflx_j23(k ,is,j) )
645  call check( __line__, qflx_j23(k-1,is,j) )
646  call check( __line__, qflx_hi(k ,is,j ,ydir) )
647  call check( __line__, qflx_hi(k ,is,j-1,ydir) )
648  call check( __line__, ddiv(k ,is,j) )
649  call check( __line__, ddiv(k+1,is,j) )
650  call check( __line__, momz0(k,is,j) )
651  call check( __line__, momz_t(k,is,j) )
652 #endif
653  advcv = - ( qflx_hi(k,is,j,zdir) - qflx_hi(k-1,is,j ,zdir) &
654  + qflx_j23(k,is,j) - qflx_j23(k-1,is,j ) ) * rfdz(k)
655  advch = - ( qflx_hi(k,is,j,ydir) - qflx_hi(k, is,j-1,ydir) ) * rcdy(j) &
656  * mapf(is,j,2,i_xy)
657  wdmp = - wdamp_coef(k) * momz0(k,is,j)
658  div = divdmp_coef / dtrk * ( ddiv(k+1,is,j)-ddiv(k,is,j) ) * fdz(k) ! divergence damping
659  sw(k,is,j) = ( advcv + advch ) / gsqrt(k,is,j,i_xyw) &
660  + wdmp + div + momz_t(k,is,j)
661 #ifdef HIST_TEND
662  if ( lhist ) then
663  advcv_t(k,is,j,i_momz) = advcv / gsqrt(k,is,j,i_xyw)
664  advch_t(k,is,j,i_momz) = advch / gsqrt(k,is,j,i_xyw)
665  wdmp_t(k,is,j) = wdmp
666  ddiv_t(k,is,j,1) = div
667  endif
668 #endif
669  enddo
670  enddo
671  !$acc end kernels
672  else
673  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
674  !$omp private(i,j,k,advcv,advch,cf,wdmp,div) &
675 #ifdef HIST_TEND
676  !$omp shared(lhist,advcv_t,advch_t,wdmp_t,ddiv_t) &
677 #endif
678  !$omp shared(JJS,JJE,IIS,IIE,KS,KE) &
679  !$omp shared(qflx_hi,qflx_J13,qflx_J23,DDIV,MOMZ0,MOMZ_t,Sw) &
680  !$omp shared(RFDZ,RCDX,RCDY,FDZ,dtrk,wdamp_coef,divdmp_coef) &
681  !$omp shared(MAPF,GSQRT,I_XY,I_XYW)
682  !$acc kernels
683  do j = jjs, jje
684  do i = iis, iie
685 !OCL NORECURRENCE
686  do k = ks, ke-1
687 #ifdef DEBUG
688  call check( __line__, qflx_hi(k ,i ,j ,zdir) )
689  call check( __line__, qflx_hi(k-1,i ,j ,zdir) )
690  call check( __line__, qflx_j13(k ,i ,j) )
691  call check( __line__, qflx_j13(k-1,i ,j) )
692  call check( __line__, qflx_j23(k ,i ,j) )
693  call check( __line__, qflx_j23(k-1,i ,j) )
694  call check( __line__, qflx_hi(k ,i ,j ,xdir) )
695  call check( __line__, qflx_hi(k ,i-1,j ,xdir) )
696  call check( __line__, qflx_hi(k ,i ,j ,ydir) )
697  call check( __line__, qflx_hi(k ,i ,j-1,ydir) )
698  call check( __line__, ddiv(k ,i,j) )
699  call check( __line__, ddiv(k+1,i,j) )
700  call check( __line__, momz0(k,i,j) )
701  call check( __line__, momz_t(k,i,j) )
702 #endif
703  advcv = - ( qflx_hi(k,i,j,zdir) - qflx_hi(k-1,i ,j ,zdir) &
704  + qflx_j13(k,i,j) - qflx_j13(k-1,i ,j ) &
705  + qflx_j23(k,i,j) - qflx_j23(k-1,i ,j ) ) * rfdz(k)
706  advch = - ( ( qflx_hi(k,i,j,xdir) - qflx_hi(k,i-1,j ,xdir) ) * rcdx(i) &
707  + ( qflx_hi(k,i,j,ydir) - qflx_hi(k,i ,j-1,ydir) ) * rcdy(j) ) &
708  * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy)
709  wdmp = - wdamp_coef(k) * momz0(k,i,j)
710  div = divdmp_coef / dtrk * ( ddiv(k+1,i,j)-ddiv(k,i,j) ) * fdz(k) ! divergence damping
711  sw(k,i,j) = ( advcv + advch ) / gsqrt(k,i,j,i_xyw) &
712  + wdmp + div + momz_t(k,i,j)
713 #ifdef HIST_TEND
714  if ( lhist ) then
715  advcv_t(k,i,j,i_momz) = advcv / gsqrt(k,i,j,i_xyw)
716  advch_t(k,i,j,i_momz) = advch / gsqrt(k,i,j,i_xyw)
717  wdmp_t(k,i,j) = wdmp
718  ddiv_t(k,i,j,1) = div
719  endif
720 #endif
721  enddo
722  enddo
723  enddo
724  !$acc end kernels
725  end if
726  profile_stop("hevi_sw")
727 #ifdef DEBUG
728  k = iundef; i = iundef; j = iundef
729 #endif
730 
731 
732  !##### Thermodynamic Equation #####
733 
734  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
735  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,JHALO,IHALO,RHOT,DENS,POTT)
736  !$acc kernels
737  do j = jjs-jhalo, jje+jhalo
738  do i = iis-ihalo, iie+ihalo
739  do k = ks, ke
740 #ifdef DEBUG
741  call check( __line__, rhot(k,i,j) )
742  call check( __line__, dens(k,i,j) )
743 #endif
744  pott(k,i,j) = rhot(k,i,j) / dens(k,i,j)
745  enddo
746  enddo
747  enddo
748  !$acc end kernels
749 #ifdef DEBUG
750  k = iundef; i = iundef; j = iundef
751 #endif
752 
753  ! at (x, y, w)
754  call atmos_dyn_fvm_fluxz_xyz( tflx_hi(:,:,:,zdir), & ! (out)
755  mflx_hi(:,:,:,zdir), pott, gsqrt(:,:,:,i_xyw), & ! (in)
756  num_diff(:,:,:,i_rhot,zdir), & ! (in)
757  cdz, & ! (in)
758  iis, iie, jjs, jje ) ! (in)
759 
760  ! at (u, y, z)
761  if ( .not. twod ) &
762  call atmos_dyn_fvm_fluxx_xyz( tflx_hi(:,:,:,xdir), & ! (out)
763  mflx_hi(:,:,:,xdir), pott, gsqrt(:,:,:,i_uyz), & ! (in)
764  num_diff(:,:,:,i_rhot,xdir), & ! (in)
765  cdz, & ! (in)
766  iis, iie, jjs, jje ) ! (in)
767 
768  ! at (x, v, z)
769  call atmos_dyn_fvm_fluxy_xyz( tflx_hi(:,:,:,ydir), & ! (out)
770  mflx_hi(:,:,:,ydir), pott, gsqrt(:,:,:,i_xvz), & ! (in)
771  num_diff(:,:,:,i_rhot,ydir), & ! (in)
772  cdz, & ! (in)
773  iis, iie, jjs, jje ) ! (in)
774 
775 
776  profile_start("hevi_st")
777  if ( twod ) then
778  !$omp parallel do default(none) OMP_SCHEDULE_ &
779  !$omp private(j,k,advcv,advch) &
780 #ifdef HIST_TEND
781  !$omp shared(lhist,advcv_t,advch_t) &
782 #endif
783  !$omp shared(JJS,JJE,IS,KS,KE) &
784  !$omp shared(tflx_hi,RHOT_t,St,RCDZ,RCDY) &
785  !$omp shared(MAPF,GSQRT,I_XY,I_XYZ)
786  !$acc kernels
787  do j = jjs, jje
788  do k = ks, ke
789 #ifdef DEBUG
790  call check( __line__, tflx_hi(k ,is,j ,zdir) )
791  call check( __line__, tflx_hi(k-1,is,j ,zdir) )
792  call check( __line__, tflx_hi(k ,is,j ,ydir) )
793  call check( __line__, tflx_hi(k ,is,j-1,ydir) )
794  call check( __line__, rhot_t(k,is,j) )
795 #endif
796  advcv = - ( tflx_hi(k,is,j,zdir) - tflx_hi(k-1,is,j ,zdir) ) * rcdz(k)
797  advch = - ( tflx_hi(k,is,j,ydir) - tflx_hi(k ,is,j-1,ydir) ) * rcdy(j)
798  st(k,is,j) = ( advcv + advch ) * mapf(is,j,2,i_xy) / gsqrt(k,is,j,i_xyz) + rhot_t(k,is,j)
799 #ifdef HIST_TEND
800  if ( lhist ) then
801  advch_t(k,is,j,i_rhot) = ( advcv + advch ) * mapf(is,j,2,i_xy) / gsqrt(k,is,j,i_xyz)
802  endif
803 #endif
804  enddo
805  enddo
806  !$acc end kernels
807  else
808  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
809  !$omp private(i,j,k,advcv,advch) &
810 #ifdef HIST_TEND
811  !$omp shared(lhist,advcv_t,advch_t) &
812 #endif
813  !$omp shared(JJS,JJE,IIS,IIE,KS,KE) &
814  !$omp shared(tflx_hi,RHOT_t,St,RCDZ,RCDX,RCDY) &
815  !$omp shared(MAPF,GSQRT,I_XY,I_XYZ)
816  !$acc kernels
817  do j = jjs, jje
818  do i = iis, iie
819  do k = ks, ke
820 #ifdef DEBUG
821  call check( __line__, tflx_hi(k ,i ,j ,zdir) )
822  call check( __line__, tflx_hi(k-1,i ,j ,zdir) )
823  call check( __line__, tflx_hi(k ,i ,j ,xdir) )
824  call check( __line__, tflx_hi(k ,i-1,j ,xdir) )
825  call check( __line__, tflx_hi(k ,i ,j ,ydir) )
826  call check( __line__, tflx_hi(k ,i ,j-1,ydir) )
827  call check( __line__, rhot_t(k,i,j) )
828 #endif
829  advcv = - ( tflx_hi(k,i,j,zdir) - tflx_hi(k-1,i ,j ,zdir) ) * rcdz(k)
830  advch = - ( ( tflx_hi(k,i,j,xdir) - tflx_hi(k ,i-1,j ,xdir) ) * rcdx(i) &
831  + ( tflx_hi(k,i,j,ydir) - tflx_hi(k ,i ,j-1,ydir) ) * rcdy(j) )
832  st(k,i,j) = ( advcv + advch ) * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz) + rhot_t(k,i,j)
833 #ifdef HIST_TEND
834  if ( lhist ) then
835  advch_t(k,i,j,i_rhot) = ( advcv + advch ) * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz)
836  endif
837 #endif
838  enddo
839  enddo
840  enddo
841  !$acc end kernels
842  end if
843 #ifdef DEBUG
844  k = iundef; i = iundef; j = iundef
845 #endif
846  profile_stop("hevi_st")
847 
848  ! implicit solver
849 
850  profile_start("hevi_solver")
851 
852  call prof_rapstart("DYN_HEVI", 3)
853 
854 !OCL INDEPENDENT
855 !OCL PREFETCH_SEQUENTIAL(SOFT)
856 #ifndef __GFORTRAN__
857  !$omp parallel do default(none) OMP_SCHEDULE_ &
858  !$omp private(k,i,j,ii,l,A,B,Ci,Co,F1,F2,F3,PT,pg,advcv) &
859 #ifdef HIST_TEND
860  !$omp shared(lhist,pg_t,advcv_t) &
861 #endif
862 #ifdef DEBUG
863  !$omp shared(RHOT) &
864 #endif
865 #if defined DEBUG || defined QUICKDEBUG
866  !$omp shared(UNDEF) &
867 #endif
868  !$omp shared(JJS,JJE,IIS,IIE,KA,KMAX,KS,KE) &
869  !$omp shared(mflx_hi,tflx_hi,MOMZ_RK,MOMZ0) &
870  !$omp shared(DENS_RK,RHOT_RK,DENS0,RHOT0,DENS,MOMZ,POTT,DPRES) &
871  !$omp shared(GRAV,dtrk,REF_dens,Sr,Sw,St,RT2P) &
872  !$omp shared(ATMOS_DYN_FVM_flux_valueW_Z) &
873  !$omp shared(MAPF,GSQRT,J33G,I_XY,I_XYZ,I_XYW,CDZ,RCDZ,RFDZ)
874 #else
875  !$omp parallel do default(shared) private(i,j,k,ii,l) OMP_SCHEDULE_ &
876  !$omp private(A,B,Ci,Co,F1,F2,F3,PT,pg,advcv)
877 #endif
878  !$acc kernels
879  do j = jjs, jje
880 #if LSIZE == 1
881  !$acc loop private(F1,F2,F3,PT,Ci,Co,A,work)
882  do i = iis, iie
883 #else
884  do ii = iis, iie, lsize
885 
886 #if defined DEBUG || defined QUICKDEBUG
887  pt(:,:) = undef
888  ci(:,:) = undef
889  co(:,:) = undef
890 
891  f1(:,:) = undef
892  f2(:,:) = undef
893  f3(:,:) = 0.0_rp
894 #endif
895 
896  do l = 1, lsize
897  i = ii + l - 1
898  if ( i > iie ) exit
899 #endif
900 
901  call atmos_dyn_fvm_flux_valuew_z( pt(:,l), & ! (out)
902  momz(:,i,j), pott(:,i,j), gsqrt(:,i,j,i_xyz), & ! (in)
903  cdz )
904 
905  do k = ks, ke
906  a(k) = dtrk**2 * j33g * rcdz(k) * rt2p(k,i,j) * j33g / gsqrt(k,i,j,i_xyz)
907  enddo
908  b = grav * dtrk**2 * j33g / ( cdz(ks+1) + cdz(ks) )
909  f1(ks,l) = - ( pt(ks+1,l) * rfdz(ks) * a(ks+1) + b ) / gsqrt(ks,i,j,i_xyw)
910  f2(ks,l) = 1.0_rp + ( pt(ks ,l) * rfdz(ks) * ( a(ks+1)+a(ks) ) ) / gsqrt(ks,i,j,i_xyw)
911  do k = ks+1, ke-2
912  b = grav * dtrk**2 * j33g / ( cdz(k+1) + cdz(k) )
913  f1(k,l) = - ( pt(k+1,l) * rfdz(k) * a(k+1) + b ) / gsqrt(k,i,j,i_xyw)
914  f2(k,l) = 1.0_rp + ( pt(k ,l) * rfdz(k) * ( a(k+1)+a(k) ) ) / gsqrt(k,i,j,i_xyw)
915  f3(k,l) = - ( pt(k-1,l) * rfdz(k) * a(k) - b ) / gsqrt(k,i,j,i_xyw)
916  enddo
917 
918  b = grav * dtrk**2 * j33g / ( cdz(ke) + cdz(ke-1) )
919  f2(ke-1,l) = 1.0_rp + ( pt(ke-1,l) * rfdz(ke-1) * ( a(ke)+a(ke-1) ) ) / gsqrt(ke-1,i,j,i_xyw)
920  f3(ke-1,l) = - ( pt(ke-2,l) * rfdz(ke-1) * a(ke-1) - b ) / gsqrt(ke-1,i,j,i_xyw)
921  do k = ks, ke-1
922  ! use not density at the half level but mean density between CZ(k) and CZ(k+1)
923  pg = - ( dpres(k+1,i,j) + rt2p(k+1,i,j)*dtrk*st(k+1,i,j) &
924  - dpres(k ,i,j) - rt2p(k ,i,j)*dtrk*st(k ,i,j) ) &
925  * rfdz(k) * j33g / gsqrt(k,i,j,i_xyw) &
926  - grav * 0.5_rp &
927  * ( ( dens(k+1,i,j) - ref_dens(k+1,i,j) + sr(k+1,i,j) * dtrk ) &
928  + ( dens(k ,i,j) - ref_dens(k ,i,j) + sr(k ,i,j) * dtrk ) )
929  ci(k,l) = momz(k,i,j) + dtrk * ( pg + sw(k,i,j) )
930 #ifdef HIST_TEND
931  if ( lhist ) pg_t(k,i,j,1) = pg
932 #endif
933  enddo
934 
935 #if LSIZE == 1
936  call matrix_solver_tridiagonal_1d_cr( kmax-1, 1, kmax-1, &
937 #ifdef _OPENACC
938  work(:,:), &
939 #endif
940  f1(:,1), f2(:,1), f3(:,1), & ! (in)
941  ci(:,1), & ! (in)
942  co(:,1) ) ! (out)
943 #else
944  end do
945 
946  call matrix_solver_tridiagonal( kmax-1, 1, kmax-1, &
947  f1(:,:), f2(:,:), f3(:,:), & ! (in)
948  ci(:,:), & ! (in)
949  co(:,:) ) ! (out)
950 
951  do l = 1, lsize
952  i = ii + l - 1
953  if ( i > iie ) exit
954 #endif
955 
956 !OCL NORECURRENCE
957  do k = ks, ke-1
958 #ifdef DEBUG_HEVI2HEVE
959  ! for debug (change to explicit integration)
960  co(k,l) = momz(k,i,j)
961  mflx_hi(k,i,j,zdir) = mflx_hi(k,i,j,zdir) &
962  + j33g * momz(k,i,j) / ( mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) )
963  tflx_hi(k,i,j,zdir) = tflx_hi(k,i,j,zdir) &
964  + j33g * momz(k,i,j) * pt(k,l) / ( mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) )
965  ! use not density at the half level but mean density between CZ(k) and CZ(k+1)
966  momz_rk(k,i,j) = momz0(k,i,j) &
967  + dtrk*( &
968  - j33g * ( dpres(k+1,i,j)-dpres(k,i,j) ) * rfdz(k) / gsqrt(k,i,j,i_xyw) &
969  - grav * 0.5_rp * ( (dens(k,i,j)-ref_dens(k,i,j)) + (dens(k+1,i,j)-ref_dens(k+1,i,j)) ) &
970  + sw(k,i,j) )
971 #else
972  ! z-flux
973  mflx_hi(k,i,j,zdir) = mflx_hi(k,i,j,zdir) &
974  + j33g * co(k,l) / ( mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) )
975  tflx_hi(k,i,j,zdir) = tflx_hi(k,i,j,zdir) &
976  + j33g * co(k,l) * pt(k,l) / ( mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) )
977  ! z-momentum
978  momz_rk(k,i,j) = momz0(k,i,j) &
979  + ( co(k,l) - momz(k,i,j) )
980 #endif
981  enddo
982  momz_rk(ks-1,i,j) = 0.0_rp
983  momz_rk(ke ,i,j) = 0.0_rp
984 
985  ! density and rho*theta
986  advcv = - co(ks,l) * j33g * rcdz(ks) / gsqrt(ks,i,j,i_xyz) ! Co(KS-1) = 0
987  dens_rk(ks,i,j) = dens0(ks,i,j) + dtrk * ( advcv + sr(ks,i,j) )
988 #ifdef HIST_TEND
989  if ( lhist ) advcv_t(ks,i,j,i_dens) = advcv
990 #endif
991  advcv = - co(ks,l) * pt(ks,l) * j33g * rcdz(ks) / gsqrt(ks,i,j,i_xyz) ! Co(KS-1) = 0
992  rhot_rk(ks,i,j) = rhot0(ks,i,j) + dtrk * ( advcv + st(ks,i,j) )
993 #ifdef HIST_TEND
994  if ( lhist ) advcv_t(ks,i,j,i_rhot) = advcv
995 #endif
996 !OCL NORECURRENCE
997  do k = ks+1, ke-1
998  advcv = - ( co(k,l) - co(k-1,l) ) &
999  * j33g * rcdz(k) / gsqrt(k,i,j,i_xyz)
1000  dens_rk(k,i,j) = dens0(k,i,j) + dtrk * ( advcv + sr(k,i,j) )
1001 #ifdef HIST_TEND
1002  if ( lhist ) advcv_t(k,i,j,i_dens) = advcv
1003 #endif
1004  advcv = - ( co(k,l) * pt(k,l) - co(k-1,l) * pt(k-1,l) ) &
1005  * j33g * rcdz(k) / gsqrt(k,i,j,i_xyz)
1006  rhot_rk(k,i,j) = rhot0(k,i,j) + dtrk * ( advcv + st(k,i,j) )
1007 #ifdef HIST_TEND
1008  if ( lhist ) advcv_t(k,i,j,i_rhot) = advcv
1009 #endif
1010  enddo
1011  advcv = co(ke-1,l) * j33g * rcdz(ke) / gsqrt(ke,i,j,i_xyz) ! Co(KE) = 0
1012  dens_rk(ke,i,j) = dens0(ke,i,j) + dtrk * ( advcv + sr(ke,i,j) )
1013 #ifdef HIST_TEND
1014  if ( lhist ) advcv_t(ke,i,j,i_dens) = advcv
1015 #endif
1016  advcv = co(ke-1,l) * pt(ke-1,l) * j33g * rcdz(ke) / gsqrt(ke,i,j,i_xyz) ! Co(KE) = 0
1017  rhot_rk(ke,i,j) = rhot0(ke,i,j) + dtrk * ( advcv + st(ke,i,j) )
1018 #ifdef HIST_TEND
1019  if ( lhist ) advcv_t(ke,i,j,i_rhot) = advcv
1020 #endif
1021 
1022 #ifdef DEBUG
1023  call check_equation( &
1024  co(:,l), &
1025  dens(:,i,j), momz(:,i,j), rhot(:,i,j), dpres(:,i,j), &
1026  ref_dens(:,i,j), &
1027  sr(:,i,j), sw(:,i,j), st(:,i,j), &
1028  j33g, gsqrt(:,i,j,:), &
1029  rt2p(:,i,j), &
1030  dtrk, i, j )
1031 #endif
1032 
1033 #if LSIZE > 1
1034  end do
1035 #endif
1036 
1037  enddo
1038  enddo
1039  !$acc end kernels
1040 #ifdef DEBUG
1041  k = iundef; i = iundef; j = iundef
1042 #endif
1043 
1044  call prof_rapend("DYN_HEVI", 3)
1045 
1046  profile_stop("hevi_solver")
1047 
1048 
1049  !##### momentum equation (x) #####
1050 
1051  profile_start("hevi_momx")
1052  ! at (u, y, w)
1053  call atmos_dyn_fvm_fluxz_uyz( qflx_hi(:,:,:,zdir), & ! (out)
1054  momz, momx, dens, & ! (in)
1055  gsqrt(:,:,:,i_uyw), j33g, & ! (in)
1056  num_diff(:,:,:,i_momx,zdir), & ! (in)
1057  cdz, twod, & ! (in)
1058  iis, iie, jjs, jje ) ! (in)
1059  if ( .not. twod ) &
1060  call atmos_dyn_fvm_fluxj13_uyz( qflx_j13, & ! (out)
1061  momx, momx, dens, & ! (in)
1062  gsqrt(:,:,:,i_uyz), j13g(:,:,:,i_uyw), mapf(:,:,:,i_uy), & ! (in)
1063  cdz, twod, & ! (in)
1064  iis, iie, jjs, jje ) ! (in)
1065  call atmos_dyn_fvm_fluxj23_uyz( qflx_j23, & ! (out)
1066  momy, momx, dens, & ! (in)
1067  gsqrt(:,:,:,i_uyz), j23g(:,:,:,i_uyw), mapf(:,:,:,i_uy), & ! (in)
1068  cdz, twod, & ! (in)
1069  iis, iie, jjs, jje ) ! (in)
1070 
1071  ! at (x, y, z)
1072  ! note that x-index is added by -1
1073  if ( .not. twod ) &
1074  call atmos_dyn_fvm_fluxx_uyz( qflx_hi(:,:,:,xdir), & ! (out)
1075  momx, momx, dens, & ! (in)
1076  gsqrt(:,:,:,i_xyz), mapf(:,:,:,i_xy), & ! (in)
1077  num_diff(:,:,:,i_momx,xdir), & ! (in)
1078  cdz, twod, & ! (in)
1079  iis, iie, jjs, jje ) ! (in)
1080 
1081  ! at (u, v, z)
1082  call atmos_dyn_fvm_fluxy_uyz( qflx_hi(:,:,:,ydir), & ! (out)
1083  momy, momx, dens, & ! (in)
1084  gsqrt(:,:,:,i_uvz), mapf(:,:,1,i_uv), & ! (in)
1085  num_diff(:,:,:,i_momx,ydir), & ! (in)
1086  cdz, twod, & ! (in)
1087  iis, iie, jjs, jje ) ! (in)
1088 
1089  !--- update momentum(x)
1090  if ( twod ) then
1091  !$omp parallel do default(none) OMP_SCHEDULE_ &
1092  !$omp private(j,k,advch,advcv,cf) &
1093 #ifdef HIST_TEND
1094  !$omp shared(lhist,advch_t,advcv_t,pg_t,cf_t,ddiv_t) &
1095 #endif
1096  !$omp shared(JJS,JJE,IS,KS,KE) &
1097  !$omp shared(MOMX_RK,DENS,MOMX,MOMY,MOMX0,MOMX_t) &
1098  !$omp shared(qflx_hi,qflx_J23) &
1099  !$omp shared(RCDZ,RCDY,CDZ) &
1100  !$omp shared(MAPF,GSQRT,I_XY,I_UY,I_XYZ,I_UYZ) &
1101  !$omp shared(dtrk,CORIOLI)
1102  !$acc kernels
1103  do j = jjs, jje
1104 !OCL NORECURRENCE
1105  do k = ks, ke
1106 #ifdef DEBUG
1107  call check( __line__, qflx_hi(k ,is,j ,zdir) )
1108  call check( __line__, qflx_hi(k-1,is,j ,zdir) )
1109  call check( __line__, qflx_hi(k ,is,j ,ydir) )
1110  call check( __line__, qflx_hi(k ,is,j-1,ydir) )
1111  call check( __line__, corioli(is,j) )
1112  call check( __line__, momy(k,is,j ) )
1113  call check( __line__, momy(k,is,j-1) )
1114  call check( __line__, momx0(k,is,j) )
1115 #endif
1116  advcv = - ( qflx_hi(k,is,j,zdir) - qflx_hi(k-1,is,j ,zdir) &
1117  + qflx_j23(k,is,j) - qflx_j23(k-1,is,j ) ) * rcdz(k)
1118  advch = - ( qflx_hi(k,is,j,ydir) - qflx_hi(k ,is,j-1,ydir) ) * rcdy(j) &
1119  * mapf(is,j,2,i_uy)
1120  cf = 0.5_rp * corioli(is,j) * ( momy(k,is,j) + momy(k,is,j-1) )
1121  momx_rk(k,is,j) = momx0(k,is,j) &
1122  + dtrk * ( ( advcv + advch ) / gsqrt(k,is,j,i_uyz) + cf + momx_t(k,is,j) )
1123 #ifdef HIST_TEND
1124  if ( lhist ) then
1125  advcv_t(k,is,j,i_momx) = advcv / gsqrt(k,is,j,i_uyz)
1126  advch_t(k,is,j,i_momx) = advch / gsqrt(k,is,j,i_uyz)
1127  pg_t(k,is,j,2) = 0.0_rp
1128  cf_t(k,is,j,1) = cf
1129  ddiv_t(k,is,j,2) = 0.0_rp
1130  endif
1131 #endif
1132  enddo
1133  enddo
1134  !$acc end kernels
1135  else
1136  iee = min(iie,ieh)
1137  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
1138  !$omp private(i,j,k,advch,advcv,pg,cf,div) &
1139 #ifdef HIST_TEND
1140  !$omp shared(lhist,advch_t,advcv_t,pg_t,cf_t,ddiv_t) &
1141 #endif
1142  !$omp shared(JJS,JJE,IIS,iee,KS,KE) &
1143  !$omp shared(MOMX_RK,DPRES,DENS,MOMX,MOMY,DDIV,MOMX0,MOMX_t) &
1144  !$omp shared(qflx_hi,qflx_J13,qflx_J23) &
1145  !$omp shared(RCDZ,RCDY,RFDX,CDZ,FDX) &
1146  !$omp shared(MAPF,GSQRT,J13G,I_XY,I_UY,I_UV,I_XYZ,I_UYW,I_UYZ) &
1147  !$omp shared(dtrk,CORIOLI,divdmp_coef)
1148  !$acc kernels
1149  do j = jjs, jje
1150  do i = iis, iee
1151 !OCL NORECURRENCE
1152  do k = ks, ke
1153 #ifdef DEBUG
1154  call check( __line__, qflx_hi(k ,i ,j ,zdir) )
1155  call check( __line__, qflx_hi(k-1,i ,j ,zdir) )
1156  call check( __line__, qflx_hi(k ,i ,j ,xdir) )
1157  call check( __line__, qflx_hi(k ,i-1,j ,xdir) )
1158  call check( __line__, qflx_hi(k ,i ,j ,ydir) )
1159  call check( __line__, qflx_hi(k ,i ,j-1,ydir) )
1160  call check( __line__, dpres(k,i+1,j) )
1161  call check( __line__, dpres(k,i ,j) )
1162  call check( __line__, corioli(i ,j) )
1163  call check( __line__, corioli(i+1,j) )
1164  call check( __line__, momy(k,i ,j ) )
1165  call check( __line__, momy(k,i+1,j ) )
1166  call check( __line__, momy(k,i ,j-1) )
1167  call check( __line__, momy(k,i+1,j-1) )
1168  call check( __line__, ddiv(k,i+1,j) )
1169  call check( __line__, ddiv(k,i ,j) )
1170  call check( __line__, momx0(k,i,j) )
1171 #endif
1172  advcv = - ( qflx_hi(k,i,j,zdir) - qflx_hi(k-1,i ,j ,zdir) &
1173  + qflx_j13(k,i,j) - qflx_j13(k-1,i ,j ) &
1174  + qflx_j23(k,i,j) - qflx_j23(k-1,i ,j ) ) * rcdz(k)
1175  advch = - ( ( qflx_hi(k,i,j,xdir) - qflx_hi(k ,i-1,j ,xdir) ) * rfdx(i) &
1176  + ( qflx_hi(k,i,j,ydir) - qflx_hi(k ,i ,j-1,ydir) ) * rcdy(j) ) &
1177  * mapf(i,j,1,i_uy) * mapf(i,j,2,i_uy)
1178  pg = ( ( gsqrt(k,i+1,j,i_xyz) * dpres(k,i+1,j) & ! [x,y,z]
1179  - gsqrt(k,i ,j,i_xyz) * dpres(k,i ,j) & ! [x,y,z]
1180  ) * rfdx(i) &
1181  + ( j13g(k ,i,j,i_uyw) &
1182  * 0.5_rp * ( f2h(k,1,i_uyz) * ( dpres(k+1,i+1,j)+dpres(k+1,i,j) ) &
1183  + f2h(k,2,i_uyz) * ( dpres(k ,i+1,j)+dpres(k ,i,j) ) ) & ! [x,y,z->u,y,w]
1184  - j13g(k-1,i,j,i_uyw) &
1185  * 0.5_rp * ( f2h(k,1,i_uyz) * ( dpres(k ,i+1,j)+dpres(k ,i,j) ) &
1186  + f2h(k,2,i_uyz) * ( dpres(k-1,i+1,j)+dpres(k-1,i,j) ) ) & ! [x,y,z->u,y,w]
1187  ) * rcdz(k) ) &
1188  * mapf(i,j,1,i_uy)
1189  cf = 0.125_rp * ( corioli(i+1,j )+corioli(i,j ) ) & ! [x,y,z->u,y,z]
1190  * ( momy(k,i+1,j )+momy(k,i,j ) &
1191  + momy(k,i+1,j-1)+momy(k,i,j-1) ) & ! [x,v,z->u,y,z]
1192  + 0.25_rp * mapf(i,j,1,i_uy) * mapf(i,j,2,i_uy) &
1193  * ( momy(k,i,j) + momy(k,i,j-1) + momy(k,i+1,j) + momy(k,i+1,j-1) ) &
1194  * ( ( momy(k,i,j) + momy(k,i,j-1) + momy(k,i+1,j) + momy(k,i+1,j-1) ) * 0.25_rp &
1195  * ( 1.0_rp/mapf(i+1,j,2,i_xy) - 1.0_rp/mapf(i,j,2,i_xy) ) * rfdx(i) &
1196  - momx(k,i,j) &
1197  * ( 1.0_rp/mapf(i,j,1,i_uv) - 1.0_rp/mapf(i,j-1,1,i_uv) ) * rcdy(j) ) &
1198  * 2.0_rp / ( dens(k,i+1,j) + dens(k,i,j) ) ! metric term
1199  div = divdmp_coef / dtrk * ( ddiv(k,i+1,j)/mapf(i+1,j,2,i_xy) - ddiv(k,i,j)/mapf(i,j,1,i_xy) ) &
1200  * mapf(i,j,1,i_uy) * mapf(i,j,2,i_uy) * fdx(i) ! divergence damping
1201  momx_rk(k,i,j) = momx0(k,i,j) &
1202  + dtrk * ( ( advcv + advch - pg ) / gsqrt(k,i,j,i_uyz) + cf + div + momx_t(k,i,j) )
1203 #ifdef HIST_TEND
1204  if ( lhist ) then
1205  advcv_t(k,i,j,i_momx) = advcv / gsqrt(k,i,j,i_uyz)
1206  advch_t(k,i,j,i_momx) = advch / gsqrt(k,i,j,i_uyz)
1207  pg_t(k,i,j,2) = - pg / gsqrt(k,i,j,i_uyz)
1208  cf_t(k,i,j,1) = cf
1209  ddiv_t(k,i,j,2) = div
1210  endif
1211 #endif
1212  enddo
1213  enddo
1214  enddo
1215  !$acc end kernels
1216  end if
1217  profile_stop("hevi_momx")
1218 #ifdef DEBUG
1219  k = iundef; i = iundef; j = iundef
1220 #endif
1221 
1222  !##### momentum equation (y) #####
1223  profile_start("hevi_momy")
1224  ! at (x, v, w)
1225  call atmos_dyn_fvm_fluxz_xvz( qflx_hi(:,:,:,zdir), & ! (out)
1226  momz, momy, dens, & ! (in)
1227  gsqrt(:,:,:,i_xvw), j33g, & ! (in)
1228  num_diff(:,:,:,i_momy,zdir), & ! (in)
1229  cdz, twod, & ! (in)
1230  iis, iie, jjs, jje ) ! (in)
1231  if ( .not. twod ) &
1232  call atmos_dyn_fvm_fluxj13_xvz( qflx_j13, & ! (out)
1233  momx, momy, dens, & ! (in)
1234  gsqrt(:,:,:,i_xvz), j13g(:,:,:,i_xvw), mapf(:,:,:,i_xv), & ! (in)
1235  cdz, twod, & ! (in)
1236  iis, iie, jjs, jje ) ! (in)
1237  call atmos_dyn_fvm_fluxj23_xvz( qflx_j23, & ! (out)
1238  momy, momy, dens, & ! (in)
1239  gsqrt(:,:,:,i_xvz), j23g(:,:,:,i_xvw), mapf(:,:,:,i_xv), & ! (in)
1240  cdz, twod, & ! (in)
1241  iis, iie, jjs, jje ) ! (in)
1242 
1243  ! at (u, v, z)
1244  if ( .not. twod ) &
1245  call atmos_dyn_fvm_fluxx_xvz( qflx_hi(:,:,:,xdir), & ! (out)
1246  momx, momy, dens, & ! (in)
1247  gsqrt(:,:,:,i_uvz), mapf(:,:,:,i_uv), & ! (in)
1248  num_diff(:,:,:,i_momy,xdir), & ! (in)
1249  cdz, twod, & ! (in)
1250  iis, iie, jjs, jje ) ! (in)
1251 
1252  ! at (x, y, z)
1253  ! note that y-index is added by -1
1254  call atmos_dyn_fvm_fluxy_xvz( qflx_hi(:,:,:,ydir), & ! (out)
1255  momy, momy, dens, & ! (in)
1256  gsqrt(:,:,:,i_xyz), mapf(:,:,:,i_xy), & ! (in)
1257  num_diff(:,:,:,i_momy,ydir), & ! (in
1258  cdz, twod, & ! (in)
1259  iis, iie, jjs, jje ) ! (in)
1260 
1261  !--- update momentum(y)
1262  if ( twod ) then
1263  i = is
1264  !$omp parallel do default(none) OMP_SCHEDULE_ &
1265  !$omp private(j,k,advch,advcv,pg,cf,div) &
1266 #ifdef HIST_TEND
1267  !$omp shared(lhist,advch_t,advcv_t,pg_t,cf_t,ddiv_t) &
1268 #endif
1269  !$omp shared(JJS,JJE,JEH,IS,i,KS,KE) &
1270  !$omp shared(MOMY_RK,DPRES,DENS,MOMX,MOMY,DDIV,MOMY0,MOMY_t) &
1271  !$omp shared(qflx_hi,qflx_J23) &
1272  !$omp shared(RCDZ,RFDY,CDZ,FDY) &
1273  !$omp shared(MAPF,GSQRT,J23G,I_XV,I_XYZ,I_XVW,I_XVZ) &
1274  !$omp shared(dtrk,CORIOLI,divdmp_coef)
1275  !$acc kernels
1276  do j = jjs, min(jje,jeh)
1277 !OCL NORECURRENCE
1278  do k = ks, ke
1279 #ifdef DEBUG
1280  call check( __line__, qflx_hi(k ,is,j ,zdir) )
1281  call check( __line__, qflx_hi(k-1,is,j ,zdir) )
1282  call check( __line__, qflx_hi(k ,is,j ,ydir) )
1283  call check( __line__, qflx_hi(k ,is,j-1,ydir) )
1284  call check( __line__, dpres(k,is,j ) )
1285  call check( __line__, dpres(k,is,j+1) )
1286  call check( __line__, corioli(is,j ) )
1287  call check( __line__, corioli(is,j+1) )
1288  call check( __line__, momx(k,is,j ) )
1289  call check( __line__, momx(k,is,j+1) )
1290  call check( __line__, ddiv(k,is,j+1) )
1291  call check( __line__, ddiv(k,is,j ) )
1292  call check( __line__, momy_t(k,is,j) )
1293  call check( __line__, momy0(k,is,j) )
1294 #endif
1295  advcv = - ( qflx_hi(k,is,j,zdir) - qflx_hi(k-1,is,j ,zdir) &
1296  + qflx_j23(k,is,j) - qflx_j23(k-1,is,j ) ) * rcdz(k)
1297  advch = - ( qflx_hi(k,is,j,ydir) - qflx_hi(k ,is,j-1,ydir) ) * rfdy(j) &
1298  * mapf(is,j,2,i_xv)
1299  pg = ( ( gsqrt(k,is,j+1,i_xyz) * dpres(k,is,j+1) & ! [x,y,z]
1300  - gsqrt(k,is,j ,i_xyz) * dpres(k,is,j ) & ! [x,y,z]
1301  ) * rfdy(j) &
1302  + ( j23g(k ,is,j,i_xvw) &
1303  * 0.5_rp * ( f2h(k ,1,i_xvz) * ( dpres(k+1,is,j+1)+dpres(k+1,is,j) ) &
1304  + f2h(k ,2,i_xvz) * ( dpres(k ,is,j+1)+dpres(k ,is,j) ) ) & ! [x,y,z->x,v,w]
1305  - j23g(k-1,is,j,i_xvw) &
1306  * 0.5_rp * ( f2h(k-1,1,i_xvz) * ( dpres(k ,is,j+1)+dpres(k ,is,j) ) &
1307  + f2h(k-1,2,i_xvz) * ( dpres(k-1,is,j+1)+dpres(k-1,is,j) ) ) & ! [x,y,z->x,v,w]
1308  ) * rcdz(k) ) &
1309  * mapf(is,j,2,i_xv)
1310  cf = - 0.25_rp * ( corioli( is,j+1)+corioli( is,j) ) & ! [x,y,z->x,v,z]
1311  * ( momx(k,is,j+1)+momx(k,is,j) )
1312  div = divdmp_coef / dtrk * ( ddiv(k,is,j+1) - ddiv(k,is,j) ) &
1313  * mapf(is,j,2,i_xv) * fdy(j) ! divergence damping
1314  momy_rk(k,is,j) = momy0(k,is,j) &
1315  + dtrk * ( ( advcv + advch - pg ) / gsqrt(k,is,j,i_xvz) + cf + div + momy_t(k,is,j) )
1316 #ifdef HIST_TEND
1317  if ( lhist ) then
1318  advcv_t(k,is,j,i_momy) = advcv / gsqrt(k,is,j,i_xvz)
1319  advch_t(k,is,j,i_momy) = advch / gsqrt(k,is,j,i_xvz)
1320  pg_t(k,is,j,3) = - pg / gsqrt(k,is,j,i_xvz)
1321  cf_t(k,is,j,2) = cf
1322  ddiv_t(k,is,j,3) = div
1323  endif
1324 #endif
1325  enddo
1326  enddo
1327  !$acc end kernels
1328  else
1329  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
1330  !$omp private(i,j,k,advch,advcv,pg,cf,div) &
1331 #ifdef HIST_TEND
1332  !$omp shared(lhist,advch_t,advcv_t,pg_t,cf_t,ddiv_t) &
1333 #endif
1334  !$omp shared(JJS,JJE,JEH,IIS,IIE,KS,KE) &
1335  !$omp shared(MOMY_RK,DPRES,DENS,MOMX,MOMY,DDIV,MOMY0,MOMY_t) &
1336  !$omp shared(qflx_hi,qflx_J13,qflx_J23) &
1337  !$omp shared(RCDZ,RCDX,RFDY,CDZ,FDY) &
1338  !$omp shared(MAPF,GSQRT,J23G,I_XY,I_XV,I_UV,I_XYZ,I_XVW,I_XVZ) &
1339  !$omp shared(dtrk,CORIOLI,divdmp_coef)
1340  !$acc kernels
1341  do j = jjs, min(jje,jeh)
1342  do i = iis, iie
1343 !OCL NORECURRENCE
1344  do k = ks, ke
1345 #ifdef DEBUG
1346  call check( __line__, qflx_hi(k ,i ,j ,zdir) )
1347  call check( __line__, qflx_hi(k-1,i ,j ,zdir) )
1348  call check( __line__, qflx_hi(k ,i ,j ,xdir) )
1349  call check( __line__, qflx_hi(k ,i-1,j ,xdir) )
1350  call check( __line__, qflx_hi(k ,i ,j ,ydir) )
1351  call check( __line__, qflx_hi(k ,i ,j-1,ydir) )
1352  call check( __line__, dpres(k,i,j ) )
1353  call check( __line__, dpres(k,i,j+1) )
1354  call check( __line__, corioli(i,j ) )
1355  call check( __line__, corioli(i,j+1) )
1356  call check( __line__, momx(k,i ,j ) )
1357  call check( __line__, momx(k,i ,j+1) )
1358  call check( __line__, momx(k,i-1,j ) )
1359  call check( __line__, momx(k,i-1,j+1) )
1360  call check( __line__, ddiv(k,i,j+1) )
1361  call check( __line__, ddiv(k,i,j ) )
1362  call check( __line__, momy_t(k,i,j) )
1363  call check( __line__, momy0(k,i,j) )
1364 #endif
1365  advcv = - ( qflx_hi(k,i,j,zdir) - qflx_hi(k-1,i ,j ,zdir) &
1366  + qflx_j13(k,i,j) - qflx_j13(k-1,i ,j ) &
1367  + qflx_j23(k,i,j) - qflx_j23(k-1,i ,j ) ) * rcdz(k)
1368  advch = - ( ( qflx_hi(k,i,j,xdir) - qflx_hi(k ,i-1,j ,xdir) ) * rcdx(i) &
1369  + ( qflx_hi(k,i,j,ydir) - qflx_hi(k ,i ,j-1,ydir) ) * rfdy(j) ) &
1370  * mapf(i,j,1,i_xv) * mapf(i,j,2,i_xv)
1371  pg = ( ( gsqrt(k,i,j+1,i_xyz) * dpres(k,i,j+1) & ! [x,y,z]
1372  - gsqrt(k,i,j ,i_xyz) * dpres(k,i,j ) & ! [x,y,z]
1373  ) * rfdy(j) &
1374  + ( j23g(k ,i,j,i_xvw) &
1375  * 0.5_rp * ( f2h(k ,1,i_xvz) * ( dpres(k+1,i,j+1)+dpres(k+1,i,j) ) &
1376  + f2h(k ,2,i_xvz) * ( dpres(k ,i,j+1)+dpres(k ,i,j) ) ) & ! [x,y,z->x,v,w]
1377  - j23g(k-1,i,j,i_xvw) &
1378  * 0.5_rp * ( f2h(k-1,1,i_xvz) * ( dpres(k ,i,j+1)+dpres(k ,i,j) ) &
1379  + f2h(k-1,2,i_xvz) * ( dpres(k-1,i,j+1)+dpres(k-1,i,j) ) ) & ! [x,y,z->x,v,w]
1380  ) * rcdz(k) ) &
1381  * mapf(i,j,2,i_xv)
1382  cf = - 0.125_rp * ( corioli(i ,j+1)+corioli(i ,j) ) & ! [x,y,z->x,v,z]
1383  * ( momx(k,i ,j+1)+momx(k,i ,j) &
1384  + momx(k,i-1,j+1)+momx(k,i-1,j) ) & ! [u,y,z->x,v,z]
1385  - 0.25_rp * mapf(i,j,1,i_xv) * mapf(i,j,2,i_xv) &
1386  * ( momx(k,i,j) + momx(k,i-1,j) + momx(k,i,j+1) + momx(k,i-1,j+1) ) &
1387  * ( momy(k,i,j) &
1388  * ( 1.0_rp/mapf(i,j,2,i_uv) - 1.0_rp/mapf(i-1,j,2,i_uv) ) * rcdx(i) &
1389  - 0.25_rp * ( momx(k,i,j)+momx(k,i-1,j)+momx(k,i,j+1)+momx(k,i-1,j+1) ) &
1390  * ( 1.0_rp/mapf(i,j+1,1,i_xy) - 1.0_rp/mapf(i,j,1,i_xy) ) * rfdy(j) ) &
1391  * 2.0_rp / ( dens(k,i,j+1) + dens(k,i,j) ) ! metoric term
1392  div = divdmp_coef / dtrk * ( ddiv(k,i,j+1)/mapf(i,j+1,1,i_xy) - ddiv(k,i,j)/mapf(i,j,1,i_xy) ) &
1393  * mapf(i,j,1,i_xv) * mapf(i,j,2,i_xv) * fdy(j) ! divergence damping
1394  momy_rk(k,i,j) = momy0(k,i,j) &
1395  + dtrk * ( ( advcv + advch - pg ) / gsqrt(k,i,j,i_xvz) + cf + div + momy_t(k,i,j) )
1396 #ifdef HIST_TEND
1397  if ( lhist ) then
1398  advcv_t(k,i,j,i_momy) = advcv / gsqrt(k,i,j,i_xvz)
1399  advch_t(k,i,j,i_momy) = advch / gsqrt(k,i,j,i_xvz)
1400  pg_t(k,i,j,3) = - pg / gsqrt(k,i,j,i_xvz)
1401  cf_t(k,i,j,2) = cf
1402  ddiv_t(k,i,j,3) = div
1403  endif
1404 #endif
1405  enddo
1406  enddo
1407  enddo
1408  !$acc end kernels
1409  end if
1410  profile_stop("hevi_momy")
1411 #ifdef DEBUG
1412  k = iundef; i = iundef; j = iundef
1413 #endif
1414 
1415  enddo
1416  enddo
1417 
1418 #ifdef PROFILE_FIPP
1419  call fipp_stop()
1420 #endif
1421 
1422 #ifdef HIST_TEND
1423  if ( lhist ) then
1424  call file_history_in(advcv_t(:,:,:,i_dens), 'DENS_t_advcv', 'tendency of density (vert. advection) (w/ HIST_TEND)', 'kg/m3/s' )
1425  call file_history_in(advcv_t(:,:,:,i_momz), 'MOMZ_t_advcv', 'tendency of momentum z (vert. advection) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZHXY')
1426  call file_history_in(advcv_t(:,:,:,i_momx), 'MOMX_t_advcv', 'tendency of momentum x (vert. advection) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZXHY')
1427  call file_history_in(advcv_t(:,:,:,i_momy), 'MOMY_t_advcv', 'tendency of momentum y (vert. advection) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZXYH')
1428  call file_history_in(advcv_t(:,:,:,i_rhot), 'RHOT_t_advcv', 'tendency of rho*theta (vert. advection) (w/ HIST_TEND)', 'K kg/m3/s' )
1429 
1430  call file_history_in(advch_t(:,:,:,i_dens), 'DENS_t_advch', 'tendency of density (horiz. advection) (w/ HIST_TEND)', 'kg/m3/s' )
1431  call file_history_in(advch_t(:,:,:,i_momz), 'MOMZ_t_advch', 'tendency of momentum z (horiz. advection) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZHXY')
1432  call file_history_in(advch_t(:,:,:,i_momx), 'MOMX_t_advch', 'tendency of momentum x (horiz. advection) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZXHY')
1433  call file_history_in(advch_t(:,:,:,i_momy), 'MOMY_t_advch', 'tendency of momentum y (horiz. advection) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZXYH')
1434  call file_history_in(advch_t(:,:,:,i_rhot), 'RHOT_t_advch', 'tendency of rho*theta (horiz. advection) (w/ HIST_TEND)', 'K kg/m3/s' )
1435 
1436  call file_history_in(pg_t(:,:,:,1), 'MOMZ_t_pg', 'tendency of momentum z (pressure gradient) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZHXY')
1437  if ( .not. twod ) &
1438  call file_history_in(pg_t(:,:,:,2), 'MOMX_t_pg', 'tendency of momentum x (pressure gradient) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZXHY')
1439  call file_history_in(pg_t(:,:,:,3), 'MOMY_t_pg', 'tendency of momentum y (pressure gradient) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZXYH')
1440 
1441  call file_history_in(wdmp_t(:,:,:), 'MOMZ_t_wdamp', 'tendency of momentum z (Rayleigh damping) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZHXY')
1442 
1443  call file_history_in(ddiv_t(:,:,:,1), 'MOMZ_t_ddiv', 'tendency of momentum z (divergence damping) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZHXY')
1444  if ( .not. twod ) &
1445  call file_history_in(ddiv_t(:,:,:,2), 'MOMX_t_ddiv', 'tendency of momentum x (divergence damping) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZXHY')
1446  call file_history_in(ddiv_t(:,:,:,3), 'MOMY_t_ddiv', 'tendency of momentum y (divergence damping) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZXYH')
1447 
1448  call file_history_in(cf_t(:,:,:,1), 'MOMX_t_cf', 'tendency of momentum x (coliolis force) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZXHY')
1449  call file_history_in(cf_t(:,:,:,2), 'MOMY_t_cf', 'tendency of momentum y (coliolis force) (w/ HIST_TEND)', 'kg/m2/s2', dim_type='ZXYH')
1450  endif
1451 #endif
1452 
1453  !$acc end data
1454 
1455  return
1456  end subroutine atmos_dyn_tstep_short_fvm_hevi
1457 
1458 #ifdef DEBUG
1459 !OCL SERIAL
1460  subroutine check_equation( &
1461  VECT, &
1462  DENS, MOMZ, RHOT, DPRES, &
1463  REF_dens, &
1464  Sr, Sw, St, &
1465  J33G, G, &
1466  RT2P, &
1467  dt, i, j )
1468  use scale_const, only: &
1469  eps => const_eps, &
1470  grav => const_grav
1471  use scale_prc, only: &
1472  prc_abort
1473  use scale_atmos_grid_cartesc, only: &
1474  rcdz => atmos_grid_cartesc_rcdz, &
1475  rfdz => atmos_grid_cartesc_rfdz
1476  implicit none
1477  real(rp), intent(in) :: vect(ks:ke-1)
1478  real(rp), intent(in) :: dens(ka)
1479  real(rp), intent(in) :: momz(ka)
1480  real(rp), intent(in) :: rhot(ka)
1481  real(rp), intent(in) :: dpres(ka)
1482  real(rp), intent(in) :: ref_dens(ka)
1483  real(rp), intent(in) :: sr(ka)
1484  real(rp), intent(in) :: sw(ka)
1485  real(rp), intent(in) :: st(ka)
1486  real(rp), intent(in) :: j33g
1487  real(rp), intent(in) :: g(ka,8)
1488  real(rp), intent(in) :: rt2p(ka)
1489  real(rp), intent(in) :: dt
1490  integer , intent(in) :: i
1491  integer , intent(in) :: j
1492 
1493  real(rp), parameter :: small = 1e-6_rp
1494 
1495  real(rp) :: momz_n(ka)
1496  real(rp) :: dens_n(ka)
1497  real(rp) :: rhot_n(ka)
1498  real(rp) :: dpres_n(ka)
1499 
1500  real(rp) :: pott(ka)
1501  real(rp) :: pt(ka)
1502 
1503  real(rp) :: error, lhs, rhs
1504  integer :: k
1505 
1506 
1507  do k = ks, ke-1
1508  momz_n(k) = vect(k)
1509  enddo
1510  momz_n(:ks-1) = 0.0_rp
1511  momz_n(ke:) = 0.0_rp
1512 
1513  ! density
1514  do k = ks+1, ke-1
1515  dens_n(k) = dens(k) &
1516  + dt * ( - j33g * ( momz_n(k) - momz_n(k-1) ) * rcdz(k) / g(k,i_xyz) + sr(k) )
1517  enddo
1518  dens_n(ks) = dens(ks) &
1519  + dt * ( - j33g * momz_n(ks) * rcdz(ks) / g(ks,i_xyz) + sr(ks) )
1520  dens_n(ke) = dens(ke) &
1521  + dt * ( j33g * momz_n(ke-1) * rcdz(ke) / g(ke,i_xyz) + sr(ke) )
1522 
1523  ! rho*theta
1524  do k = ks, ke
1525  pott(k) = rhot(k) / dens(k)
1526  enddo
1527  do k = ks+1, ke-2
1528  pt(k) = ( 7.0_rp * ( pott(k+1) + pott(k ) ) &
1529  - ( pott(k+2) + pott(k-1) ) ) / 12.0_rp
1530  enddo
1531  pt(ks-1) = 0.0_rp
1532  pt(ks ) = ( pott(ks+1) + pott(ks ) ) * 0.5_rp
1533  pt(ke-1) = ( pott(ke ) + pott(ke-1) ) * 0.5_rp
1534  pt(ke ) = 0.0_rp
1535  do k = ks+1, ke-1
1536  rhot_n(k) = rhot(k) &
1537  + dt * ( - j33g * ( momz_n(k)*pt(k) - momz_n(k-1)*pt(k-1) ) * rcdz(k) / g(k,i_xyz) &
1538  + st(k) )
1539  enddo
1540  rhot_n(ks) = rhot(ks) &
1541  + dt * ( - j33g * momz_n(ks)*pt(ks) * rcdz(ks) / g(ks,i_xyz) + st(ks) )
1542  rhot_n(ke) = rhot(ke) &
1543  + dt * ( j33g * momz_n(ke-1)*pt(ke-1) * rcdz(ke) / g(ke-1,i_xyz) + st(ke) )
1544 
1545 
1546  do k = ks, ke
1547  dpres_n(k) = dpres(k) + rt2p(k) * ( rhot_n(k) - rhot(k) )
1548  enddo
1549 
1550  do k = ks, ke
1551  lhs = ( dens_n(k) - dens(k) ) / dt
1552  rhs = - j33g * ( momz_n(k) - momz_n(k-1) ) * rcdz(k) / g(k,i_xyz) + sr(k)
1553  if ( abs(lhs) < small ) then
1554  error = rhs
1555  else
1556  error = ( lhs - rhs ) / lhs
1557  endif
1558  if ( abs(error) > small ) then
1559  log_error("check_equation",*)"DENS error", k, i, j, error, lhs, rhs
1560  log_error_cont(*)eps
1561  call prc_abort
1562  endif
1563  enddo
1564 
1565  do k = ks, ke-1
1566  lhs = ( momz_n(k) - momz(k) ) / dt
1567  rhs = - j33g * ( dpres_n(k+1) - dpres_n(k) ) * rfdz(k) / g(k,i_xyw) &
1568  - grav * ( ( dens_n(k+1) - ref_dens(k+1) ) + ( dens_n(k) - ref_dens(k) ) ) * 0.5_rp &
1569  + sw(k)
1570  if ( abs(lhs) < small ) then
1571  error = rhs
1572  else
1573  error = ( lhs - rhs ) / lhs
1574  endif
1575  if ( abs(error) > small ) then
1576  log_error("check_equation",*)"MOMZ error", k, i, j, error, lhs, rhs
1577  log_error_cont(*) momz_n(k), momz(k), dt
1578  log_error_cont(*) - j33g * ( dpres(k+1) - dpres(k) ) * rfdz(k) / g(k,i_xyw) &
1579  - grav * ( ( dens(k+1) - ref_dens(k+1) ) + ( dens(k) - ref_dens(k) ) ) * 0.5_rp &
1580  + sw(k)
1581  call prc_abort
1582  endif
1583  enddo
1584 
1585  do k = ks, ke
1586  lhs = ( rhot_n(k) - rhot(k) ) / dt
1587  rhs = - j33g * ( momz_n(k)*pt(k) - momz_n(k-1)*pt(k-1) ) * rcdz(k) / g(k,i_xyz) + st(k)
1588  if ( abs(lhs) < small ) then
1589  error = rhs
1590  else
1591  error = ( lhs - rhs ) / lhs
1592  endif
1593  if ( abs(error) > small ) then
1594  log_error("check_equation",*)"RHOT error", k, i, j, error, lhs, rhs
1595  call prc_abort
1596  endif
1597  enddo
1598 
1599  return
1600  end subroutine check_equation
1601 #endif
1602 
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::i_uy
integer, public i_uy
Definition: scale_atmos_grid_cartesC_index.F90:100
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_atmos_grid_cartesc_index::i_xv
integer, public i_xv
Definition: scale_atmos_grid_cartesC_index.F90:101
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxx_uyz
procedure(flux_mom), pointer, public atmos_dyn_fvm_fluxx_uyz
Definition: scale_atmos_dyn_fvm_flux.F90:214
scale_index::i_momz
integer, parameter, public i_momz
Definition: scale_index.F90:29
scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxy_xyz
procedure(flux_phi), pointer, public atmos_dyn_fvm_fluxy_xyz
Definition: scale_atmos_dyn_fvm_flux.F90:175
scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxj23_xyw
procedure(flux_j), pointer, public atmos_dyn_fvm_fluxj23_xyw
Definition: scale_atmos_dyn_fvm_flux.F90:203
scale_atmos_dyn_fvm_fct::atmos_dyn_fvm_fct
subroutine, public atmos_dyn_fvm_fct(qflx_anti, phi_in, DENS0, DENS, qflx_hi, qflx_lo, mflx_hi, rdz, rdx, rdy, GSQRT, MAPF, TwoD, dt, flag_vect)
Setup.
Definition: scale_atmos_dyn_fvm_fct.F90:72
scale_index::i_momx
integer, parameter, public i_momx
Definition: scale_index.F90:30
scale_index
module Index
Definition: scale_index.F90:11
scale_atmos_grid_cartesc_index::ihalo
integer, public ihalo
Definition: scale_atmos_grid_cartesC_index.F90:44
scale_const::const_undef2
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:40
scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxj23_uyz
procedure(flux_j), pointer, public atmos_dyn_fvm_fluxj23_uyz
Definition: scale_atmos_dyn_fvm_flux.F90:230
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_dyn_fvm_flux
module scale_atmos_dyn_fvm_flux
Definition: scale_atmos_dyn_fvm_flux.F90:13
scale_matrix::matrix_solver_tridiagonal_1d_cr
subroutine, public matrix_solver_tridiagonal_1d_cr( KA, KS, KE, ifdef _OPENACC
solve tridiagonal matrix with Cyclic Reduction method
Definition: scale_matrix.F90:186
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_index::i_uv
integer, public i_uv
Definition: scale_atmos_grid_cartesC_index.F90:102
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:35
scale_prof::prof_rapstart
subroutine, public prof_rapstart(rapname_base, level, disable_barrier)
Start raptime.
Definition: scale_prof.F90:174
scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxj13_xyw
procedure(flux_j), pointer, public atmos_dyn_fvm_fluxj13_xyw
Definition: scale_atmos_dyn_fvm_flux.F90:198
scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxx_xyz
procedure(flux_phi), pointer, public atmos_dyn_fvm_fluxx_xyz
Definition: scale_atmos_dyn_fvm_flux.F90:170
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_dyn_fvm_flux::atmos_dyn_fvm_fluxx_xvz
procedure(flux_mom), pointer, public atmos_dyn_fvm_fluxx_xvz
Definition: scale_atmos_dyn_fvm_flux.F90:241
scale_atmos_dyn_tstep_short_fvm_hevi::atmos_dyn_tstep_short_fvm_hevi
subroutine, public atmos_dyn_tstep_short_fvm_hevi(DENS_RK, MOMZ_RK, MOMX_RK, MOMY_RK, RHOT_RK, PROG_RK, mflx_hi, tflx_hi, DENS0, MOMZ0, MOMX0, MOMY0, RHOT0, DENS, MOMZ, MOMX, MOMY, RHOT, DENS_t, MOMZ_t, MOMX_t, MOMY_t, RHOT_t, PROG0, PROG, DPRES0, RT2P, CORIOLI, num_diff, wdamp_coef, divdmp_coef, DDIV, FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_ALONG_STREAM, CDZ, FDZ, FDX, FDY, RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, PHI, GSQRT, J13G, J23G, J33G, MAPF, REF_dens, REF_rhot, BND_W, BND_E, BND_S, BND_N, TwoD, dtrk, last)
Definition: scale_atmos_dyn_tstep_short_fvm_hevi.F90:143
scale_atmos_grid_cartesc_index::jeh
integer, public jeh
end point of inner domain: y, local (half level)
Definition: scale_atmos_grid_cartesC_index.F90:69
scale_index::i_dens
integer, parameter, public i_dens
Definition: scale_index.F90:28
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_file_history
module file_history
Definition: scale_file_history.F90:15
scale_index::i_momy
integer, parameter, public i_momy
Definition: scale_index.F90:31
scale_atmos_grid_cartesc_index::i_xyz
integer, public i_xyz
Definition: scale_atmos_grid_cartesC_index.F90:91
scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxy_xyw
procedure(flux_mom), pointer, public atmos_dyn_fvm_fluxy_xyw
Definition: scale_atmos_dyn_fvm_flux.F90:192
scale_index::va
integer, public va
Definition: scale_index.F90:35
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxx_xyw
procedure(flux_mom), pointer, public atmos_dyn_fvm_fluxx_xyw
Definition: scale_atmos_dyn_fvm_flux.F90:187
scale_atmos_grid_cartesc_index::i_xy
integer, public i_xy
Definition: scale_atmos_grid_cartesC_index.F90:99
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_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_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::kmax
integer, public kmax
Definition: scale_atmos_grid_cartesC_index.F90:36
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_atmos_grid_cartesc_index::ieh
integer, public ieh
end point of inner domain: x, local (half level)
Definition: scale_atmos_grid_cartesC_index.F90:68
scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxz_xvz
procedure(flux_z), pointer, public atmos_dyn_fvm_fluxz_xvz
Definition: scale_atmos_dyn_fvm_flux.F90:236
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_atmos_dyn_tstep_short_fvm_hevi::atmos_dyn_tstep_short_fvm_hevi_regist
subroutine, public atmos_dyn_tstep_short_fvm_hevi_regist(ATMOS_DYN_TYPE, VA_out, VAR_NAME, VAR_DESC, VAR_UNIT)
Register.
Definition: scale_atmos_dyn_tstep_short_fvm_hevi.F90:83
scale_index::i_rhot
integer, parameter, public i_rhot
Definition: scale_index.F90:32
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_atmos_grid_cartesc_index::jhalo
integer, public jhalo
Definition: scale_atmos_grid_cartesC_index.F90:45
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_flux_valuew_z
procedure(valuew), pointer, public atmos_dyn_fvm_flux_valuew_z
Definition: scale_atmos_dyn_fvm_flux.F90:160
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxz_uyz
procedure(flux_z), pointer, public atmos_dyn_fvm_fluxz_uyz
Definition: scale_atmos_dyn_fvm_flux.F90:209
scale_atmos_grid_cartesc_index::xdir
integer, parameter, public xdir
Definition: scale_atmos_grid_cartesC_index.F90:33
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_dyn_fvm_flux::atmos_dyn_fvm_fluxj23_xvz
procedure(flux_j), pointer, public atmos_dyn_fvm_fluxj23_xvz
Definition: scale_atmos_dyn_fvm_flux.F90:257
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_dyn_fvm_flux::atmos_dyn_fvm_fluxj13_xvz
procedure(flux_j), pointer, public atmos_dyn_fvm_fluxj13_xvz
Definition: scale_atmos_dyn_fvm_flux.F90:252
scale_atmos_dyn_tstep_short_fvm_hevi
module Atmosphere / Dynamics RK
Definition: scale_atmos_dyn_tstep_short_fvm_hevi.F90:21
scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxz_xyz
procedure(flux_phi), pointer, public atmos_dyn_fvm_fluxz_xyz
Definition: scale_atmos_dyn_fvm_flux.F90:165
scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxj13_uyz
procedure(flux_j), pointer, public atmos_dyn_fvm_fluxj13_uyz
Definition: scale_atmos_dyn_fvm_flux.F90:225
scale_prof::prof_rapend
subroutine, public prof_rapend(rapname_base, level, disable_barrier)
Save raptime.
Definition: scale_prof.F90:246
scale_atmos_grid_cartesc
module atmosphere / grid / cartesC
Definition: scale_atmos_grid_cartesC.F90:12
scale_atmos_grid_cartesc_index::i_uvz
integer, public i_uvz
Definition: scale_atmos_grid_cartesC_index.F90:97
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_atmos_dyn_fvm_fct
module Atmosphere / Dynamics common
Definition: scale_atmos_dyn_fvm_fct.F90:12
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_dyn_tstep_short_fvm_hevi::atmos_dyn_tstep_short_fvm_hevi_setup
subroutine, public atmos_dyn_tstep_short_fvm_hevi_setup
Setup.
Definition: scale_atmos_dyn_tstep_short_fvm_hevi.F90:116
scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxz_xyw
procedure(flux_wz), pointer, public atmos_dyn_fvm_fluxz_xyw
Definition: scale_atmos_dyn_fvm_flux.F90:182
scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxy_xvz
procedure(flux_mom), pointer, public atmos_dyn_fvm_fluxy_xvz
Definition: scale_atmos_dyn_fvm_flux.F90:246
scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_fluxy_uyz
procedure(flux_mom), pointer, public atmos_dyn_fvm_fluxy_uyz
Definition: scale_atmos_dyn_fvm_flux.F90:219
scale_atmos_grid_cartesc_index::i_xyw
integer, public i_xyw
Definition: scale_atmos_grid_cartesC_index.F90:92