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