SCALE-RM
scale_atmos_dyn_fvm_flux_cd8.F90
Go to the documentation of this file.
1 
2 !-------------------------------------------------------------------------------
11 !-------------------------------------------------------------------------------
12 ! Warning: This file was generated from atmos-rm/dynamics/scale_atmos_dyn_fvm_flux_udcd.F90.erb.
13 ! Do not edit this file.
14 !-------------------------------------------------------------------------------
15 #include "scalelib.h"
17  !-----------------------------------------------------------------------------
18  !
19  !++ used modules
20  !
21  use scale_precision
22  use scale_io
23  use scale_prof
25  use scale_index
26  use scale_tracer
27  use scale_prc
28 #ifdef DEBUG
29  use scale_debug, only: &
30  check
31  use scale_const, only: &
32  undef => const_undef, &
33  iundef => const_undef2
34 #endif
35  !-----------------------------------------------------------------------------
36  implicit none
37  private
38  !-----------------------------------------------------------------------------
39  !
40  !++ Public procedure
41  !
43 
47 
53 
59 
65 
66  !-----------------------------------------------------------------------------
67  !
68  !++ Public parameters & variables
69  !
70  !-----------------------------------------------------------------------------
71  !
72  !++ Private procedure
73  !
74 #if 1
75 #define F2H(k,p,q) (CDZ(k+p-1)*GSQRT(k+p-1,i,j)/(CDZ(k)*GSQRT(k,i,j)+CDZ(k+1)*GSQRT(k+1,i,j)))
76 #else
77 #define F2H(k,p,q) 0.5_RP
78 #endif
79  !-----------------------------------------------------------------------------
80  !
81  !++ Private parameters & variables
82  !
83 
84  real(RP), parameter :: F2 = 0.5_rp ! F2 is always used to calculate flux near boundary.
85 
86 
87  real(RP), parameter :: F41 = 7.0_rp/12.0_rp
88  real(RP), parameter :: F42 = -1.0_rp/12.0_rp
89 
90 
91  real(RP), parameter :: F61 = 37.0_rp/60.0_rp
92  real(RP), parameter :: F62 = -8.0_rp/60.0_rp
93  real(RP), parameter :: F63 = 1.0_rp/60.0_rp
94 
95 
96  real(RP), parameter :: F81 = 533.0_rp/840.0_rp
97  real(RP), parameter :: F82 = -139.0_rp/840.0_rp
98  real(RP), parameter :: F83 = 29.0_rp/840.0_rp
99  real(RP), parameter :: F84 = -3.0_rp/840.0_rp
100 
101 
102 
103 contains
104  !-----------------------------------------------------------------------------
106 !OCL SERIAL
107  subroutine atmos_dyn_fvm_flux_valuew_z_cd8( &
108  valW, &
109  mflx, val, GSQRT, &
110  CDZ )
111  !$acc routine vector
112  implicit none
113 
114  real(rp), intent(out) :: valw (ka)
115  real(rp), intent(in) :: mflx (ka)
116  real(rp), intent(in) :: val (ka)
117  real(rp), intent(in) :: gsqrt(ka)
118  real(rp), intent(in) :: cdz (ka)
119 
120  integer :: k
121  !---------------------------------------------------------------------------
122 
123  do k = ks+3, ke-4
124 #ifdef DEBUG
125  call check( __line__, mflx(k) )
126 
127  call check( __line__, val(k) )
128  call check( __line__, val(k+1) )
129 
130  call check( __line__, val(k-1) )
131  call check( __line__, val(k+2) )
132 
133  call check( __line__, val(k-2) )
134  call check( __line__, val(k+3) )
135 
136  call check( __line__, val(k-3) )
137  call check( __line__, val(k+4) )
138 
139 #endif
140  valw(k) = f81 * ( val(k+1)+val(k) ) &
141  + f82 * ( val(k+2)+val(k-1) ) &
142  + f83 * ( val(k+3)+val(k-2) ) &
143  + f84 * ( val(k+4)+val(k-3) )
144  enddo
145 #ifdef DEBUG
146  k = iundef
147 #endif
148 
149 #ifdef DEBUG
150 
151  call check( __line__, mflx(ks) )
152  call check( __line__, val(ks ) )
153  call check( __line__, val(ks+1) )
154  call check( __line__, mflx(ke-1) )
155  call check( __line__, val(ke ) )
156  call check( __line__, val(ke-1) )
157 
158  call check( __line__, mflx(ks+1) )
159  call check( __line__, val(ks+2 ) )
160  call check( __line__, val(ks+3) )
161  call check( __line__, mflx(ke-2) )
162  call check( __line__, val(ke-2 ) )
163  call check( __line__, val(ke-3) )
164 
165  call check( __line__, mflx(ks+2) )
166  call check( __line__, val(ks+4 ) )
167  call check( __line__, val(ks+5) )
168  call check( __line__, mflx(ke-3) )
169  call check( __line__, val(ke-4 ) )
170  call check( __line__, val(ke-5) )
171 
172 #endif
173 
174  valw(ks) = f2 * ( val(ks+1)+val(ks) )
175  valw(ke-1) = f2 * ( val(ke)+val(ke-1) )
176 
177  valw(ks+1) = f41 * ( val(ks+2)+val(ks+1) ) &
178  + f42 * ( val(ks+3)+val(ks) )
179  valw(ke-2) = f41 * ( val(ke-1)+val(ke-2) ) &
180  + f42 * ( val(ke)+val(ke-3) )
181 
182  valw(ks+2) = f61 * ( val(ks+3)+val(ks+2) ) &
183  + f62 * ( val(ks+4)+val(ks+1) ) &
184  + f63 * ( val(ks+5)+val(ks) )
185  valw(ke-3) = f61 * ( val(ke-2)+val(ke-3) ) &
186  + f62 * ( val(ke-1)+val(ke-4) ) &
187  + f63 * ( val(ke)+val(ke-5) )
188 
189 
190  return
191  end subroutine atmos_dyn_fvm_flux_valuew_z_cd8
192 
193  !-----------------------------------------------------------------------------
195  subroutine atmos_dyn_fvm_fluxz_xyz_cd8( &
196  flux, &
197  mflx, val, GSQRT, &
198  num_diff, &
199  CDZ, &
200  IIS, IIE, JJS, JJE )
201  use scale_const, only: &
202  eps => const_eps
203  implicit none
204 
205  real(rp), intent(inout) :: flux (ka,ia,ja)
206  real(rp), intent(in) :: mflx (ka,ia,ja)
207  real(rp), intent(in) :: val (ka,ia,ja)
208  real(rp), intent(in) :: gsqrt (ka,ia,ja)
209  real(rp), intent(in) :: num_diff(ka,ia,ja)
210  real(rp), intent(in) :: cdz (ka)
211  integer, intent(in) :: iis, iie, jjs, jje
212 
213  real(rp) :: vel
214  integer :: k, i, j
215  !---------------------------------------------------------------------------
216 
217  !$omp parallel default(none) private(i,j,k, vel) &
218  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mflx,val,flux,GSQRT,num_diff,EPS)
219 
220  !$acc data copy(flux) copyin(mflx, val, GSQRT, num_diff, CDZ)
221 
222  !$omp do OMP_SCHEDULE_ collapse(2)
223  !$acc kernels
224  do j = jjs, jje
225  do i = iis, iie
226  do k = ks+3, ke-4
227 #ifdef DEBUG
228  call check( __line__, mflx(k,i,j) )
229 
230  call check( __line__, val(k,i,j) )
231  call check( __line__, val(k+1,i,j) )
232 
233  call check( __line__, val(k-1,i,j) )
234  call check( __line__, val(k+2,i,j) )
235 
236  call check( __line__, val(k-2,i,j) )
237  call check( __line__, val(k+3,i,j) )
238 
239  call check( __line__, val(k-3,i,j) )
240  call check( __line__, val(k+4,i,j) )
241 
242 #endif
243  vel = mflx(k,i,j)
244  flux(k,i,j) = vel &
245  * ( f81 * ( val(k+1,i,j)+val(k,i,j) ) &
246  + f82 * ( val(k+2,i,j)+val(k-1,i,j) ) &
247  + f83 * ( val(k+3,i,j)+val(k-2,i,j) ) &
248  + f84 * ( val(k+4,i,j)+val(k-3,i,j) ) ) &
249  + gsqrt(k,i,j) * num_diff(k,i,j)
250  enddo
251  enddo
252  enddo
253  !$acc end kernels
254  !$omp end do nowait
255 #ifdef DEBUG
256  k = iundef; i = iundef; j = iundef
257 #endif
258 
259  !$omp do OMP_SCHEDULE_ collapse(2)
260  !$acc kernels
261  do j = jjs, jje
262  do i = iis, iie
263 #ifdef DEBUG
264 
265  call check( __line__, mflx(ks,i,j) )
266  call check( __line__, val(ks ,i,j) )
267  call check( __line__, val(ks+1,i,j) )
268  call check( __line__, mflx(ke-1,i,j) )
269  call check( __line__, val(ke ,i,j) )
270  call check( __line__, val(ke-1,i,j) )
271 
272  call check( __line__, mflx(ks+1,i,j) )
273  call check( __line__, val(ks+2 ,i,j) )
274  call check( __line__, val(ks+3,i,j) )
275  call check( __line__, mflx(ke-2,i,j) )
276  call check( __line__, val(ke-2 ,i,j) )
277  call check( __line__, val(ke-3,i,j) )
278 
279  call check( __line__, mflx(ks+2,i,j) )
280  call check( __line__, val(ks+4 ,i,j) )
281  call check( __line__, val(ks+5,i,j) )
282  call check( __line__, mflx(ke-3,i,j) )
283  call check( __line__, val(ke-4 ,i,j) )
284  call check( __line__, val(ke-5,i,j) )
285 
286 #endif
287  flux(ks-1,i,j) = 0.0_rp
288 
289  vel = mflx(ks,i,j)
290  flux(ks,i,j) = vel &
291  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) ) &
292  + gsqrt(ks,i,j) * num_diff(ks,i,j)
293  vel = mflx(ke-1,i,j)
294  flux(ke-1,i,j) = vel &
295  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) ) &
296  + gsqrt(ke-1,i,j) * num_diff(ke-1,i,j)
297 
298  vel = mflx(ks+1,i,j)
299  flux(ks+1,i,j) = vel &
300  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
301  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) ) &
302  + gsqrt(ks+1,i,j) * num_diff(ks+1,i,j)
303  vel = mflx(ke-2,i,j)
304  flux(ke-2,i,j) = vel &
305  * ( f41 * ( val(ke-1,i,j)+val(ke-2,i,j) ) &
306  + f42 * ( val(ke,i,j)+val(ke-3,i,j) ) ) &
307  + gsqrt(ke-2,i,j) * num_diff(ke-2,i,j)
308 
309  vel = mflx(ks+2,i,j)
310  flux(ks+2,i,j) = vel &
311  * ( f61 * ( val(ks+3,i,j)+val(ks+2,i,j) ) &
312  + f62 * ( val(ks+4,i,j)+val(ks+1,i,j) ) &
313  + f63 * ( val(ks+5,i,j)+val(ks,i,j) ) ) &
314  + gsqrt(ks+2,i,j) * num_diff(ks+2,i,j)
315  vel = mflx(ke-3,i,j)
316  flux(ke-3,i,j) = vel &
317  * ( f61 * ( val(ke-2,i,j)+val(ke-3,i,j) ) &
318  + f62 * ( val(ke-1,i,j)+val(ke-4,i,j) ) &
319  + f63 * ( val(ke,i,j)+val(ke-5,i,j) ) ) &
320  + gsqrt(ke-3,i,j) * num_diff(ke-3,i,j)
321 
322  flux(ke ,i,j) = 0.0_rp
323  enddo
324  enddo
325  !$acc end kernels
326  !$omp end do nowait
327 
328  !$acc end data
329 
330  !$omp end parallel
331 #ifdef DEBUG
332  k = iundef; i = iundef; j = iundef
333 #endif
334 
335  return
336  end subroutine atmos_dyn_fvm_fluxz_xyz_cd8
337 
338  !-----------------------------------------------------------------------------
340  subroutine atmos_dyn_fvm_fluxx_xyz_cd8( &
341  flux, &
342  mflx, val, GSQRT, &
343  num_diff, &
344  CDZ, &
345  IIS, IIE, JJS, JJE )
346  implicit none
347 
348  real(rp), intent(inout) :: flux (ka,ia,ja)
349  real(rp), intent(in) :: mflx (ka,ia,ja)
350  real(rp), intent(in) :: val (ka,ia,ja)
351  real(rp), intent(in) :: gsqrt (ka,ia,ja)
352  real(rp), intent(in) :: num_diff(ka,ia,ja)
353  real(rp), intent(in) :: cdz(ka)
354  integer, intent(in) :: iis, iie, jjs, jje
355 
356  real(rp) :: vel
357  integer :: k, i, j
358  !---------------------------------------------------------------------------
359 
360  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
361  !$omp private(vel) &
362  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mflx,val,flux,GSQRT,num_diff)
363  !$acc kernels
364  do j = jjs, jje
365  do i = iis-1, iie
366  do k = ks, ke
367 #ifdef DEBUG
368  call check( __line__, mflx(k,i,j) )
369 
370  call check( __line__, val(k,i,j) )
371  call check( __line__, val(k,i+1,j) )
372 
373  call check( __line__, val(k,i-1,j) )
374  call check( __line__, val(k,i+2,j) )
375 
376  call check( __line__, val(k,i-2,j) )
377  call check( __line__, val(k,i+3,j) )
378 
379  call check( __line__, val(k,i-3,j) )
380  call check( __line__, val(k,i+4,j) )
381 
382 #endif
383  vel = mflx(k,i,j)
384  flux(k,i,j) = vel &
385  * ( f81 * ( val(k,i+1,j)+val(k,i,j) ) &
386  + f82 * ( val(k,i+2,j)+val(k,i-1,j) ) &
387  + f83 * ( val(k,i+3,j)+val(k,i-2,j) ) &
388  + f84 * ( val(k,i+4,j)+val(k,i-3,j) ) ) &
389  + gsqrt(k,i,j) * num_diff(k,i,j)
390  enddo
391  enddo
392  enddo
393  !$acc end kernels
394 #ifdef DEBUG
395  k = iundef; i = iundef; j = iundef
396 #endif
397 
398  return
399  end subroutine atmos_dyn_fvm_fluxx_xyz_cd8
400 
401  !-----------------------------------------------------------------------------
403  subroutine atmos_dyn_fvm_fluxy_xyz_cd8( &
404  flux, &
405  mflx, val, GSQRT, &
406  num_diff, &
407  CDZ, &
408  IIS, IIE, JJS, JJE )
409  implicit none
410 
411  real(rp), intent(inout) :: flux (ka,ia,ja)
412  real(rp), intent(in) :: mflx (ka,ia,ja)
413  real(rp), intent(in) :: val (ka,ia,ja)
414  real(rp), intent(in) :: gsqrt (ka,ia,ja)
415  real(rp), intent(in) :: num_diff(ka,ia,ja)
416  real(rp), intent(in) :: cdz(ka)
417  integer, intent(in) :: iis, iie, jjs, jje
418 
419  real(rp) :: vel
420  integer :: k, i, j
421  !---------------------------------------------------------------------------
422 
423  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
424  !$omp private(vel) &
425  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mflx,val,flux,GSQRT,num_diff)
426  !$acc kernels
427  do j = jjs-1, jje
428  do i = iis, iie
429  do k = ks, ke
430 #ifdef DEBUG
431  call check( __line__, mflx(k,i,j) )
432 
433  call check( __line__, val(k,i,j) )
434  call check( __line__, val(k,i,j+1) )
435 
436  call check( __line__, val(k,i,j-1) )
437  call check( __line__, val(k,i,j+2) )
438 
439  call check( __line__, val(k,i,j-2) )
440  call check( __line__, val(k,i,j+3) )
441 
442  call check( __line__, val(k,i,j-3) )
443  call check( __line__, val(k,i,j+4) )
444 
445 #endif
446  vel = mflx(k,i,j)
447  flux(k,i,j) = vel &
448  * ( f81 * ( val(k,i,j+1)+val(k,i,j) ) &
449  + f82 * ( val(k,i,j+2)+val(k,i,j-1) ) &
450  + f83 * ( val(k,i,j+3)+val(k,i,j-2) ) &
451  + f84 * ( val(k,i,j+4)+val(k,i,j-3) ) ) &
452  + gsqrt(k,i,j) * num_diff(k,i,j)
453  enddo
454  enddo
455  enddo
456  !$acc end kernels
457 #ifdef DEBUG
458  k = iundef; i = iundef; j = iundef
459 #endif
460 
461  return
462  end subroutine atmos_dyn_fvm_fluxy_xyz_cd8
463 
464 
465  !-----------------------------------------------------------------------------
467  subroutine atmos_dyn_fvm_fluxz_xyw_cd8( &
468  flux, &
469  mom, val, DENS, &
470  GSQRT, J33G, &
471  num_diff, &
472  CDZ, FDZ, &
473  dtrk, &
474  IIS, IIE, JJS, JJE )
475  implicit none
476 
477  real(rp), intent(inout) :: flux (ka,ia,ja)
478  real(rp), intent(in) :: mom (ka,ia,ja)
479  real(rp), intent(in) :: val (ka,ia,ja)
480  real(rp), intent(in) :: dens (ka,ia,ja)
481  real(rp), intent(in) :: gsqrt (ka,ia,ja)
482  real(rp), intent(in) :: j33g
483  real(rp), intent(in) :: num_diff(ka,ia,ja)
484  real(rp), intent(in) :: cdz (ka)
485  real(rp), intent(in) :: fdz (ka-1)
486  real(rp), intent(in) :: dtrk
487  integer, intent(in) :: iis, iie, jjs, jje
488 
489  real(rp) :: vel
490  integer :: k, i, j
491  !---------------------------------------------------------------------------
492 
493  ! note than z-index is added by -1
494 
495  !$omp parallel default(none) private(i,j,k,vel) &
496  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,flux,J33G,GSQRT,num_diff,DENS,FDZ,dtrk)
497 
498  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, num_diff, CDZ, FDZ)
499 
500  !$omp do OMP_SCHEDULE_ collapse(2)
501  !$acc kernels
502  do j = jjs, jje
503  do i = iis, iie
504  do k = ks+4, ke-3
505 #ifdef DEBUG
506  call check( __line__, mom(k-1,i,j) )
507  call check( __line__, mom(k ,i,j) )
508 
509  call check( __line__, val(k-1,i,j) )
510  call check( __line__, val(k,i,j) )
511 
512  call check( __line__, val(k-2,i,j) )
513  call check( __line__, val(k+1,i,j) )
514 
515  call check( __line__, val(k-3,i,j) )
516  call check( __line__, val(k+2,i,j) )
517 
518  call check( __line__, val(k-4,i,j) )
519  call check( __line__, val(k+3,i,j) )
520 
521 #endif
522  vel = ( 0.5_rp * ( mom(k-1,i,j) &
523  + mom(k,i,j) ) ) &
524  / dens(k,i,j)
525  flux(k-1,i,j) = j33g * vel &
526  * ( f81 * ( val(k,i,j)+val(k-1,i,j) ) &
527  + f82 * ( val(k+1,i,j)+val(k-2,i,j) ) &
528  + f83 * ( val(k+2,i,j)+val(k-3,i,j) ) &
529  + f84 * ( val(k+3,i,j)+val(k-4,i,j) ) ) &
530  + gsqrt(k,i,j) * num_diff(k,i,j)
531  enddo
532  enddo
533  enddo
534  !$acc end kernels
535  !$omp end do nowait
536 #ifdef DEBUG
537  k = iundef; i = iundef; j = iundef
538 #endif
539 
540  !$omp do OMP_SCHEDULE_ collapse(2)
541  !$acc kernels
542  do j = jjs, jje
543  do i = iis, iie
544 #ifdef DEBUG
545 
546  call check( __line__, val(ks,i,j) )
547  call check( __line__, val(ks+1,i,j) )
548  call check( __line__, val(ks+2,i,j) )
549  call check( __line__, val(ks+3,i,j) )
550  call check( __line__, val(ks+4,i,j) )
551  call check( __line__, val(ks+5,i,j) )
552 
553 
554  call check( __line__, val(ke-5,i,j) )
555  call check( __line__, val(ke-4,i,j) )
556  call check( __line__, val(ke-3,i,j) )
557  call check( __line__, val(ke-2,i,j) )
558  call check( __line__, val(ke-1,i,j) )
559 
560 #endif
561  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS.
562  ! The flux at KS can be non-zero.
563  ! To reduce calculations, all the fluxes are set to zero.
564  flux(ks-1,i,j) = 0.0_rp ! k = KS
565 
566  vel = ( 0.5_rp * ( mom(ks,i,j) &
567  + mom(ks+1,i,j) ) ) &
568  / dens(ks+1,i,j)
569  flux(ks,i,j) = j33g * vel &
570  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) ) &
571  + gsqrt(ks+1,i,j) * num_diff(ks+1,i,j) ! k = KS+1
572 
573  vel = ( 0.5_rp * ( mom(ks+1,i,j) &
574  + mom(ks+2,i,j) ) ) &
575  / dens(ks+2,i,j)
576  flux(ks+1,i,j) = j33g * vel &
577  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
578  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) ) &
579  + gsqrt(ks+2,i,j) * num_diff(ks+2,i,j) ! k = KS+2
580 
581  vel = ( 0.5_rp * ( mom(ks+2,i,j) &
582  + mom(ks+3,i,j) ) ) &
583  / dens(ks+3,i,j)
584  flux(ks+2,i,j) = j33g * vel &
585  * ( f61 * ( val(ks+3,i,j)+val(ks+2,i,j) ) &
586  + f62 * ( val(ks+4,i,j)+val(ks+1,i,j) ) &
587  + f63 * ( val(ks+5,i,j)+val(ks,i,j) ) ) &
588  + gsqrt(ks+3,i,j) * num_diff(ks+3,i,j) ! k = KS+3
589 
590 
591 
592  vel = ( 0.5_rp * ( mom(ke-3,i,j) &
593  + mom(ke-2,i,j) ) ) &
594  / dens(ke-1,i,j)
595  flux(ke-3,i,j) = j33g * vel &
596  * ( f61 * ( val(ke-2,i,j)+val(ke-3,i,j) ) &
597  + f62 * ( val(ke-1,i,j)+val(ke-4,i,j) ) &
598  + f63 * ( val(ke,i,j)+val(ke-5,i,j) ) ) &
599  + gsqrt(ke-2,i,j) * num_diff(ke-2,i,j) ! k = KE-2
600 
601  vel = ( 0.5_rp * ( mom(ke-2,i,j) &
602  + mom(ke-1,i,j) ) ) &
603  / dens(ke-1,i,j)
604  flux(ke-2,i,j) = j33g * vel &
605  * ( f41 * ( val(ke-1,i,j)+val(ke-2,i,j) ) &
606  + f42 * ( val(ke,i,j)+val(ke-3,i,j) ) ) &
607  + gsqrt(ke-1,i,j) * num_diff(ke-1,i,j) ! k = KE-1
608 
609  flux(ke-1,i,j) = 0.0_rp ! k = KE
610  flux(ke ,i,j) = 0.0_rp ! k = KE+1
611  enddo
612  enddo
613  !$acc end kernels
614  !$omp end do nowait
615 
616  !$acc end data
617 
618  !$omp end parallel
619 
620  return
621  end subroutine atmos_dyn_fvm_fluxz_xyw_cd8
622 
623 
624  !-----------------------------------------------------------------------------
626  subroutine atmos_dyn_fvm_fluxj13_xyw_cd8( &
627  flux, &
628  mom, val, DENS, &
629  GSQRT, J13G, MAPF, &
630  CDZ, TwoD, &
631  IIS, IIE, JJS, JJE )
632  implicit none
633 
634  real(rp), intent(inout) :: flux (ka,ia,ja)
635  real(rp), intent(in) :: mom (ka,ia,ja)
636  real(rp), intent(in) :: val (ka,ia,ja)
637  real(rp), intent(in) :: dens (ka,ia,ja)
638  real(rp), intent(in) :: gsqrt (ka,ia,ja)
639  real(rp), intent(in) :: j13g (ka,ia,ja)
640  real(rp), intent(in) :: mapf ( ia,ja,2)
641  real(rp), intent(in) :: cdz (ka)
642  logical, intent(in) :: twod
643  integer, intent(in) :: iis, iie, jjs, jje
644 
645  real(rp) :: vel
646  integer :: k, i, j
647  !---------------------------------------------------------------------------
648 
649  !$omp parallel default(none) private(i,j,k,vel) &
650  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J13G,MAPF)
651 
652  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, J13G, MAPF, CDZ)
653 
654  !$omp do OMP_SCHEDULE_ collapse(2)
655  !$acc kernels
656  do j = jjs, jje
657  do i = iis, iie
658  do k = ks+3, ke-3
659  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i-1,j) ) ) &
660  / dens(k,i,j)
661  vel = vel * j13g(k,i,j)
662  flux(k-1,i,j) = vel / mapf(i,j,+2) &
663  * ( f81 * ( val(k,i,j)+val(k-1,i,j) ) &
664  + f82 * ( val(k+1,i,j)+val(k-2,i,j) ) &
665  + f83 * ( val(k+2,i,j)+val(k-3,i,j) ) &
666  + f84 * ( val(k+3,i,j)+val(k-4,i,j) ) )
667  enddo
668  enddo
669  enddo
670  !$acc end kernels
671  !$omp end do nowait
672 
673  !$omp do OMP_SCHEDULE_ collapse(2)
674  !$acc kernels
675  do j = jjs, jje
676  do i = iis, iie
677  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS.
678  ! The flux at KS can be non-zero.
679  ! To reduce calculations, all the fluxes are set to zero.
680  flux(ks-1,i,j) = 0.0_rp ! k = KS
681 
682  ! physically incorrect but for numerical stability
683  vel = ( ( 0.5_rp * ( mom(ks+1,i,j)+mom(ks+1,i-1,j) ) ) / dens(ks+1,i,j) &
684  + ( 0.5_rp * ( mom(ks,i,j)+mom(ks,i-1,j) ) ) / dens(ks ,i,j) ) * 0.5_rp
685 ! vel = ( 0.5_RP * ( mom(KS+1,i,j)+mom(KS+1,i-1,j) ) ) &
686 ! / DENS(KS+1,i,j)
687  vel = vel * j13g(ks+1,i,j)
688  flux(ks,i,j) = vel / mapf(i,j,+2) &
689  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) ) ! k = KS+1
690 
691  vel = ( 0.5_rp * ( mom(ks+2,i,j)+mom(ks+2,i-1,j) ) ) &
692  / dens(ks+2,i,j)
693  vel = vel * j13g(ks,i,j)
694  flux(ks+1,i,j) = vel / mapf(i,j,+2) &
695  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
696  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) ) ! k = KS+2
697 
698 
699  vel = ( 0.5_rp * ( mom(ke-1,i,j)+mom(ke-1,i-1,j) ) ) &
700  / dens(ke-1,i,j)
701  vel = vel * j13g(ke-1,i,j)
702  flux(ke-2,i,j) = vel / mapf(i,j,+2) &
703  * ( f2 * ( val(ke-1,i,j)+val(ke-2,i,j) ) ) ! k = KE-3
704 
705  vel = ( 0.5_rp * ( mom(ke-2,i,j)+mom(ke-2,i-1,j) ) ) &
706  / dens(ke-2,i,j)
707  vel = vel * j13g(ke-2,i,j)
708  flux(ke-3,i,j) = vel / mapf(i,j,+2) &
709  * ( f41 * ( val(ke-2,i,j)+val(ke-3,i,j) ) &
710  + f42 * ( val(ke-1,i,j)+val(ke-4,i,j) ) ) ! k = KE-4
711 
712  flux(ke-1,i,j) = 0.0_rp
713  enddo
714  enddo
715  !$acc end kernels
716  !$omp end do nowait
717 
718  !$acc end data
719 
720  !$omp end parallel
721 
722  return
723  end subroutine atmos_dyn_fvm_fluxj13_xyw_cd8
724 
725  !-----------------------------------------------------------------------------
727  subroutine atmos_dyn_fvm_fluxj23_xyw_cd8( &
728  flux, &
729  mom, val, DENS, &
730  GSQRT, J23G, MAPF, &
731  CDZ, TwoD, &
732  IIS, IIE, JJS, JJE )
733  implicit none
734 
735  real(rp), intent(inout) :: flux (ka,ia,ja)
736  real(rp), intent(in) :: mom (ka,ia,ja)
737  real(rp), intent(in) :: val (ka,ia,ja)
738  real(rp), intent(in) :: dens (ka,ia,ja)
739  real(rp), intent(in) :: gsqrt (ka,ia,ja)
740  real(rp), intent(in) :: j23g (ka,ia,ja)
741  real(rp), intent(in) :: mapf ( ia,ja,2)
742  real(rp), intent(in) :: cdz (ka)
743  logical, intent(in) :: twod
744  integer, intent(in) :: iis, iie, jjs, jje
745 
746  real(rp) :: vel
747  integer :: k, i, j
748  !---------------------------------------------------------------------------
749 
750  !$omp parallel default(none) private(i,j,k,vel) &
751  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J23G,MAPF)
752 
753  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, J23G, MAPF, CDZ)
754 
755  !$omp do OMP_SCHEDULE_ collapse(2)
756  !$acc kernels
757  do j = jjs, jje
758  do i = iis, iie
759  do k = ks+3, ke-3
760  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i,j-1) ) ) &
761  / dens(k,i,j)
762  vel = vel * j23g(k,i,j)
763  flux(k-1,i,j) = vel / mapf(i,j,+1) &
764  * ( f81 * ( val(k,i,j)+val(k-1,i,j) ) &
765  + f82 * ( val(k+1,i,j)+val(k-2,i,j) ) &
766  + f83 * ( val(k+2,i,j)+val(k-3,i,j) ) &
767  + f84 * ( val(k+3,i,j)+val(k-4,i,j) ) )
768  enddo
769  enddo
770  enddo
771  !$acc end kernels
772  !$omp end do nowait
773 
774  !$omp do OMP_SCHEDULE_ collapse(2)
775  !$acc kernels
776  do j = jjs, jje
777  do i = iis, iie
778  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS.
779  ! The flux at KS can be non-zero.
780  ! To reduce calculations, all the fluxes are set to zero.
781  flux(ks-1,i,j) = 0.0_rp ! k = KS
782 
783  ! physically incorrect but for numerical stability
784  vel = ( ( 0.5_rp * ( mom(ks+1,i,j)+mom(ks+1,i,j-1) ) ) / dens(ks+1,i,j) &
785  + ( 0.5_rp * ( mom(ks,i,j)+mom(ks,i,j-1) ) ) / dens(ks ,i,j) ) * 0.5_rp
786 ! vel = ( 0.5_RP * ( mom(KS+1,i,j)+mom(KS+1,i,j-1) ) ) &
787 ! / DENS(KS+1,i,j)
788  vel = vel * j23g(ks+1,i,j)
789  flux(ks,i,j) = vel / mapf(i,j,+1) &
790  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) ) ! k = KS+1
791 
792  vel = ( 0.5_rp * ( mom(ks+2,i,j)+mom(ks+2,i,j-1) ) ) &
793  / dens(ks+2,i,j)
794  vel = vel * j23g(ks,i,j)
795  flux(ks+1,i,j) = vel / mapf(i,j,+1) &
796  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
797  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) ) ! k = KS+2
798 
799 
800  vel = ( 0.5_rp * ( mom(ke-1,i,j)+mom(ke-1,i,j-1) ) ) &
801  / dens(ke-1,i,j)
802  vel = vel * j23g(ke-1,i,j)
803  flux(ke-2,i,j) = vel / mapf(i,j,+1) &
804  * ( f2 * ( val(ke-1,i,j)+val(ke-2,i,j) ) ) ! k = KE-3
805 
806  vel = ( 0.5_rp * ( mom(ke-2,i,j)+mom(ke-2,i,j-1) ) ) &
807  / dens(ke-2,i,j)
808  vel = vel * j23g(ke-2,i,j)
809  flux(ke-3,i,j) = vel / mapf(i,j,+1) &
810  * ( f41 * ( val(ke-2,i,j)+val(ke-3,i,j) ) &
811  + f42 * ( val(ke-1,i,j)+val(ke-4,i,j) ) ) ! k = KE-4
812 
813  flux(ke-1,i,j) = 0.0_rp
814  enddo
815  enddo
816  !$acc end kernels
817  !$omp end do nowait
818 
819  !$acc end data
820 
821  !$omp end parallel
822 
823  return
824  end subroutine atmos_dyn_fvm_fluxj23_xyw_cd8
825 
826 
827  !-----------------------------------------------------------------------------
829  subroutine atmos_dyn_fvm_fluxx_xyw_cd8( &
830  flux, &
831  mom, val, DENS, &
832  GSQRT, MAPF, &
833  num_diff, &
834  CDZ, TwoD, &
835  IIS, IIE, JJS, JJE )
836  implicit none
837 
838  real(rp), intent(inout) :: flux (ka,ia,ja)
839  real(rp), intent(in) :: mom (ka,ia,ja)
840  real(rp), intent(in) :: val (ka,ia,ja)
841  real(rp), intent(in) :: dens (ka,ia,ja)
842  real(rp), intent(in) :: gsqrt (ka,ia,ja)
843  real(rp), intent(in) :: mapf ( ia,ja,2)
844  real(rp), intent(in) :: num_diff(ka,ia,ja)
845  real(rp), intent(in) :: cdz (ka)
846  logical, intent(in) :: twod
847  integer, intent(in) :: iis, iie, jjs, jje
848 
849  real(rp) :: vel
850  integer :: k, i, j
851  !---------------------------------------------------------------------------
852 
853  !$omp parallel default(none) private(i,j,k,vel) &
854  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff) &
855  !$omp shared(CDZ)
856 
857  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, MAPF, num_diff, CDZ)
858 
859  !$omp do OMP_SCHEDULE_ collapse(2)
860  !$acc kernels
861  do j = jjs, jje
862  do i = iis-1, iie
863  do k = ks, ke-1
864 #ifdef DEBUG
865  call check( __line__, mom(k ,i,j) )
866  call check( __line__, mom(k+1,i,j) )
867 
868  call check( __line__, val(k,i,j) )
869  call check( __line__, val(k,i+1,j) )
870 
871  call check( __line__, val(k,i-1,j) )
872  call check( __line__, val(k,i+2,j) )
873 
874  call check( __line__, val(k,i-2,j) )
875  call check( __line__, val(k,i+3,j) )
876 
877  call check( __line__, val(k,i-3,j) )
878  call check( __line__, val(k,i+4,j) )
879 
880 #endif
881  vel = ( f2h(k,1,i_uyz) &
882  * mom(k+1,i,j) &
883  + f2h(k,2,i_uyz) &
884  * mom(k,i,j) ) &
885  / ( f2h(k,1,i_uyz) &
886  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i+1,j) ) &
887  + f2h(k,2,i_uyz) &
888  * 0.5_rp * ( dens(k,i,j)+dens(k,i+1,j) ) )
889  flux(k,i,j) = gsqrt(k,i,j) / mapf(i,j,+2) * vel &
890  * ( f81 * ( val(k,i+1,j)+val(k,i,j) ) &
891  + f82 * ( val(k,i+2,j)+val(k,i-1,j) ) &
892  + f83 * ( val(k,i+3,j)+val(k,i-2,j) ) &
893  + f84 * ( val(k,i+4,j)+val(k,i-3,j) ) ) &
894  + gsqrt(k,i,j) * num_diff(k,i,j)
895  enddo
896  enddo
897  enddo
898  !$acc end kernels
899  !$omp end do nowait
900 #ifdef DEBUG
901  k = iundef; i = iundef; j = iundef
902 #endif
903 
904  !$omp do OMP_SCHEDULE_ collapse(2)
905  !$acc kernels
906  do j = jjs, jje
907  do i = iis-1, iie
908  flux(ke,i,j) = 0.0_rp
909  enddo
910  enddo
911  !$acc end kernels
912  !$omp end do nowait
913 
914  !$acc end data
915 
916  !$omp end parallel
917 #ifdef DEBUG
918  k = iundef; i = iundef; j = iundef
919 #endif
920 
921  return
922  end subroutine atmos_dyn_fvm_fluxx_xyw_cd8
923 
924  !-----------------------------------------------------------------------------
926  subroutine atmos_dyn_fvm_fluxy_xyw_cd8( &
927  flux, &
928  mom, val, DENS, &
929  GSQRT, MAPF, &
930  num_diff, &
931  CDZ, TwoD, &
932  IIS, IIE, JJS, JJE )
933  implicit none
934 
935  real(rp), intent(inout) :: flux (ka,ia,ja)
936  real(rp), intent(in) :: mom (ka,ia,ja)
937  real(rp), intent(in) :: val (ka,ia,ja)
938  real(rp), intent(in) :: dens (ka,ia,ja)
939  real(rp), intent(in) :: gsqrt (ka,ia,ja)
940  real(rp), intent(in) :: mapf ( ia,ja,2)
941  real(rp), intent(in) :: num_diff(ka,ia,ja)
942  real(rp), intent(in) :: cdz (ka)
943  logical, intent(in) :: twod
944  integer, intent(in) :: iis, iie, jjs, jje
945 
946  real(rp) :: vel
947  integer :: k, i, j
948  !---------------------------------------------------------------------------
949 
950  !$omp parallel default(none) private(i,j,k,vel) &
951  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff) &
952  !$omp shared(CDZ)
953 
954  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, MAPF, num_diff, CDZ)
955 
956  !$omp do OMP_SCHEDULE_ collapse(2)
957  !$acc kernels
958  do j = jjs-1, jje
959  do i = iis, iie
960  do k = ks, ke-1
961 #ifdef DEBUG
962  call check( __line__, mom(k ,i,j) )
963  call check( __line__, mom(k+1,i,j) )
964 
965  call check( __line__, val(k,i,j) )
966  call check( __line__, val(k,i,j+1) )
967 
968  call check( __line__, val(k,i,j-1) )
969  call check( __line__, val(k,i,j+2) )
970 
971  call check( __line__, val(k,i,j-2) )
972  call check( __line__, val(k,i,j+3) )
973 
974  call check( __line__, val(k,i,j-3) )
975  call check( __line__, val(k,i,j+4) )
976 
977 #endif
978  vel = ( f2h(k,1,i_xvz) &
979  * mom(k+1,i,j) &
980  + f2h(k,2,i_xvz) &
981  * mom(k,i,j) ) &
982  / ( f2h(k,1,i_xvz) &
983  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i,j+1) ) &
984  + f2h(k,2,i_xvz) &
985  * 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) ) )
986  flux(k,i,j) = gsqrt(k,i,j) / mapf(i,j,+1) * vel &
987  * ( f81 * ( val(k,i,j+1)+val(k,i,j) ) &
988  + f82 * ( val(k,i,j+2)+val(k,i,j-1) ) &
989  + f83 * ( val(k,i,j+3)+val(k,i,j-2) ) &
990  + f84 * ( val(k,i,j+4)+val(k,i,j-3) ) ) &
991  + gsqrt(k,i,j) * num_diff(k,i,j)
992  enddo
993  enddo
994  enddo
995  !$acc end kernels
996  !$omp end do nowait
997 #ifdef DEBUG
998  k = iundef; i = iundef; j = iundef
999 #endif
1000 
1001  !$omp do OMP_SCHEDULE_ collapse(2)
1002  !$acc kernels
1003  do j = jjs-1, jje
1004  do i = iis, iie
1005  flux(ke,i,j) = 0.0_rp
1006  enddo
1007  enddo
1008  !$acc end kernels
1009  !$omp end do nowait
1010 
1011  !$acc end data
1012 
1013  !$omp end parallel
1014 #ifdef DEBUG
1015  k = iundef; i = iundef; j = iundef
1016 #endif
1017 
1018  return
1019  end subroutine atmos_dyn_fvm_fluxy_xyw_cd8
1020 
1021 
1022  !-----------------------------------------------------------------------------
1024  subroutine atmos_dyn_fvm_fluxz_uyz_cd8( &
1025  flux, &
1026  mom, val, DENS, &
1027  GSQRT, J33G, &
1028  num_diff, &
1029  CDZ, TwoD, &
1030  IIS, IIE, JJS, JJE )
1031  implicit none
1032 
1033  real(rp), intent(inout) :: flux (ka,ia,ja)
1034  real(rp), intent(in) :: mom (ka,ia,ja)
1035  real(rp), intent(in) :: val (ka,ia,ja)
1036  real(rp), intent(in) :: dens (ka,ia,ja)
1037  real(rp), intent(in) :: gsqrt (ka,ia,ja)
1038  real(rp), intent(in) :: j33g
1039  real(rp), intent(in) :: num_diff(ka,ia,ja)
1040  real(rp), intent(in) :: cdz (ka)
1041  logical, intent(in) :: twod
1042  integer, intent(in) :: iis, iie, jjs, jje
1043 
1044  real(rp) :: vel
1045  integer :: k, i, j
1046  !---------------------------------------------------------------------------
1047 
1048  !$omp parallel default(none) private(i,j,k,vel) &
1049  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J33G,GSQRT,num_diff) &
1050  !$omp shared(CDZ,TwoD)
1051 
1052  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, num_diff, CDZ)
1053 
1054 
1055  if ( twod ) then
1056 
1057  !$omp do OMP_SCHEDULE_
1058  !$acc kernels
1059  do j = jjs, jje
1060  do k = ks+3, ke-4
1061  i = iis
1062 #ifdef DEBUG
1063  call check( __line__, mom(k,i,j) )
1064 
1065  call check( __line__, val(k,i,j) )
1066  call check( __line__, val(k+1,i,j) )
1067 
1068  call check( __line__, val(k-1,i,j) )
1069  call check( __line__, val(k+2,i,j) )
1070 
1071  call check( __line__, val(k-2,i,j) )
1072  call check( __line__, val(k+3,i,j) )
1073 
1074  call check( __line__, val(k-3,i,j) )
1075  call check( __line__, val(k+4,i,j) )
1076 
1077 #endif
1078  vel = ( mom(k,i,j) ) &
1079  / ( f2h(k,1,i_xyz) &
1080  * dens(k+1,i,j) &
1081  + f2h(k,2,i_xyz) &
1082  * dens(k,i,j) )
1083  flux(k,i,j) = j33g * vel &
1084  * ( f81 * ( val(k+1,i,j)+val(k,i,j) ) &
1085  + f82 * ( val(k+2,i,j)+val(k-1,i,j) ) &
1086  + f83 * ( val(k+3,i,j)+val(k-2,i,j) ) &
1087  + f84 * ( val(k+4,i,j)+val(k-3,i,j) ) ) &
1088  + gsqrt(k,i,j) * num_diff(k,i,j)
1089  enddo
1090  enddo
1091  !$acc end kernels
1092  !$omp end do nowait
1093 #ifdef DEBUG
1094  k = iundef; i = iundef; j = iundef
1095 #endif
1096 
1097  !$omp do OMP_SCHEDULE_
1098  !$acc kernels
1099  do j = jjs, jje
1100  i = iis
1101 #ifdef DEBUG
1102 
1103  call check( __line__, mom(ks,i ,j) )
1104  call check( __line__, val(ks+1,i,j) )
1105  call check( __line__, val(ks,i,j) )
1106 
1107  call check( __line__, mom(ks+1,i ,j) )
1108  call check( __line__, val(ks+3,i,j) )
1109  call check( __line__, val(ks+2,i,j) )
1110 
1111  call check( __line__, mom(ks+2,i ,j) )
1112  call check( __line__, val(ks+5,i,j) )
1113  call check( __line__, val(ks+4,i,j) )
1114 
1115 #endif
1116  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
1117  ! The flux at KS-1 can be non-zero.
1118  ! To reduce calculations, all the fluxes are set to zero.
1119  flux(ks-1,i,j) = 0.0_rp
1120 
1121  vel = ( mom(ks,i,j) ) &
1122  / ( f2h(ks,1,i_xyz) &
1123  * dens(ks+1,i,j) &
1124  + f2h(ks,2,i_xyz) &
1125  * dens(ks,i,j) )
1126  flux(ks,i,j) = j33g * vel &
1127  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) ) &
1128  + gsqrt(ks,i,j) * num_diff(ks,i,j)
1129  vel = ( mom(ke-1,i,j) ) &
1130  / ( f2h(ke-1,1,i_xyz) &
1131  * dens(ke,i,j) &
1132  + f2h(ke-1,2,i_xyz) &
1133  * dens(ke-1,i,j) )
1134  flux(ke-1,i,j) = j33g * vel &
1135  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) ) &
1136  + gsqrt(ke-1,i,j) * num_diff(ke-1,i,j)
1137 
1138  vel = ( mom(ks+1,i,j) ) &
1139  / ( f2h(ks+1,1,i_xyz) &
1140  * dens(ks+2,i,j) &
1141  + f2h(ks+1,2,i_xyz) &
1142  * dens(ks+1,i,j) )
1143  flux(ks+1,i,j) = j33g * vel &
1144  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
1145  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) ) &
1146  + gsqrt(ks+1,i,j) * num_diff(ks+1,i,j)
1147  vel = ( mom(ke-2,i,j) ) &
1148  / ( f2h(ke-2,1,i_xyz) &
1149  * dens(ke-1,i,j) &
1150  + f2h(ke-2,2,i_xyz) &
1151  * dens(ke-2,i,j) )
1152  flux(ke-2,i,j) = j33g * vel &
1153  * ( f41 * ( val(ke-1,i,j)+val(ke-2,i,j) ) &
1154  + f42 * ( val(ke,i,j)+val(ke-3,i,j) ) ) &
1155  + gsqrt(ke-2,i,j) * num_diff(ke-2,i,j)
1156 
1157  vel = ( mom(ks+2,i,j) ) &
1158  / ( f2h(ks+2,1,i_xyz) &
1159  * dens(ks+3,i,j) &
1160  + f2h(ks+2,2,i_xyz) &
1161  * dens(ks+2,i,j) )
1162  flux(ks+2,i,j) = j33g * vel &
1163  * ( f61 * ( val(ks+3,i,j)+val(ks+2,i,j) ) &
1164  + f62 * ( val(ks+4,i,j)+val(ks+1,i,j) ) &
1165  + f63 * ( val(ks+5,i,j)+val(ks,i,j) ) ) &
1166  + gsqrt(ks+2,i,j) * num_diff(ks+2,i,j)
1167  vel = ( mom(ke-3,i,j) ) &
1168  / ( f2h(ke-3,1,i_xyz) &
1169  * dens(ke-2,i,j) &
1170  + f2h(ke-3,2,i_xyz) &
1171  * dens(ke-3,i,j) )
1172  flux(ke-3,i,j) = j33g * vel &
1173  * ( f61 * ( val(ke-2,i,j)+val(ke-3,i,j) ) &
1174  + f62 * ( val(ke-1,i,j)+val(ke-4,i,j) ) &
1175  + f63 * ( val(ke,i,j)+val(ke-5,i,j) ) ) &
1176  + gsqrt(ke-3,i,j) * num_diff(ke-3,i,j)
1177 
1178  flux(ke,i,j) = 0.0_rp
1179  enddo
1180  !$acc end kernels
1181  !$omp end do nowait
1182 
1183  else
1184 
1185 
1186  !$omp do OMP_SCHEDULE_ collapse(2)
1187  !$acc kernels
1188  do j = jjs, jje
1189  do i = iis, iie
1190  do k = ks+3, ke-4
1191 #ifdef DEBUG
1192  call check( __line__, mom(k,i,j) )
1193  call check( __line__, mom(k,i+1,j) )
1194 
1195  call check( __line__, val(k,i,j) )
1196  call check( __line__, val(k+1,i,j) )
1197 
1198  call check( __line__, val(k-1,i,j) )
1199  call check( __line__, val(k+2,i,j) )
1200 
1201  call check( __line__, val(k-2,i,j) )
1202  call check( __line__, val(k+3,i,j) )
1203 
1204  call check( __line__, val(k-3,i,j) )
1205  call check( __line__, val(k+4,i,j) )
1206 
1207 #endif
1208  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i+1,j) ) ) &
1209  / ( f2h(k,1,i_uyz) &
1210  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i+1,j) ) &
1211  + f2h(k,2,i_uyz) &
1212  * 0.5_rp * ( dens(k,i,j)+dens(k,i+1,j) ) )
1213  flux(k,i,j) = j33g * vel &
1214  * ( f81 * ( val(k+1,i,j)+val(k,i,j) ) &
1215  + f82 * ( val(k+2,i,j)+val(k-1,i,j) ) &
1216  + f83 * ( val(k+3,i,j)+val(k-2,i,j) ) &
1217  + f84 * ( val(k+4,i,j)+val(k-3,i,j) ) ) &
1218  + gsqrt(k,i,j) * num_diff(k,i,j)
1219  enddo
1220  enddo
1221  enddo
1222  !$acc end kernels
1223  !$omp end do nowait
1224 #ifdef DEBUG
1225  k = iundef; i = iundef; j = iundef
1226 #endif
1227 
1228  !$omp do OMP_SCHEDULE_ collapse(2)
1229  !$acc kernels
1230  do j = jjs, jje
1231  do i = iis, iie
1232 #ifdef DEBUG
1233 
1234  call check( __line__, mom(ks,i ,j) )
1235  call check( __line__, mom(ks,i+1,j) )
1236  call check( __line__, val(ks+1,i,j) )
1237  call check( __line__, val(ks,i,j) )
1238 
1239  call check( __line__, mom(ks+1,i ,j) )
1240  call check( __line__, mom(ks+1,i+1,j) )
1241  call check( __line__, val(ks+3,i,j) )
1242  call check( __line__, val(ks+2,i,j) )
1243 
1244  call check( __line__, mom(ks+2,i ,j) )
1245  call check( __line__, mom(ks+2,i+1,j) )
1246  call check( __line__, val(ks+5,i,j) )
1247  call check( __line__, val(ks+4,i,j) )
1248 
1249 #endif
1250  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
1251  ! The flux at KS-1 can be non-zero.
1252  ! To reduce calculations, all the fluxes are set to zero.
1253  flux(ks-1,i,j) = 0.0_rp
1254 
1255  vel = ( 0.5_rp * ( mom(ks,i,j)+mom(ks,i+1,j) ) ) &
1256  / ( f2h(ks,1,i_uyz) &
1257  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i+1,j) ) &
1258  + f2h(ks,2,i_uyz) &
1259  * 0.5_rp * ( dens(ks,i,j)+dens(ks,i+1,j) ) )
1260  flux(ks,i,j) = j33g * vel &
1261  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) ) &
1262  + gsqrt(ks,i,j) * num_diff(ks,i,j)
1263  vel = ( 0.5_rp * ( mom(ke-1,i,j)+mom(ke-1,i+1,j) ) ) &
1264  / ( f2h(ke-1,1,i_uyz) &
1265  * 0.5_rp * ( dens(ke,i,j)+dens(ke,i+1,j) ) &
1266  + f2h(ke-1,2,i_uyz) &
1267  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i+1,j) ) )
1268  flux(ke-1,i,j) = j33g * vel &
1269  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) ) &
1270  + gsqrt(ke-1,i,j) * num_diff(ke-1,i,j)
1271 
1272  vel = ( 0.5_rp * ( mom(ks+1,i,j)+mom(ks+1,i+1,j) ) ) &
1273  / ( f2h(ks+1,1,i_uyz) &
1274  * 0.5_rp * ( dens(ks+2,i,j)+dens(ks+2,i+1,j) ) &
1275  + f2h(ks+1,2,i_uyz) &
1276  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i+1,j) ) )
1277  flux(ks+1,i,j) = j33g * vel &
1278  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
1279  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) ) &
1280  + gsqrt(ks+1,i,j) * num_diff(ks+1,i,j)
1281  vel = ( 0.5_rp * ( mom(ke-2,i,j)+mom(ke-2,i+1,j) ) ) &
1282  / ( f2h(ke-2,1,i_uyz) &
1283  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i+1,j) ) &
1284  + f2h(ke-2,2,i_uyz) &
1285  * 0.5_rp * ( dens(ke-2,i,j)+dens(ke-2,i+1,j) ) )
1286  flux(ke-2,i,j) = j33g * vel &
1287  * ( f41 * ( val(ke-1,i,j)+val(ke-2,i,j) ) &
1288  + f42 * ( val(ke,i,j)+val(ke-3,i,j) ) ) &
1289  + gsqrt(ke-2,i,j) * num_diff(ke-2,i,j)
1290 
1291  vel = ( 0.5_rp * ( mom(ks+2,i,j)+mom(ks+2,i+1,j) ) ) &
1292  / ( f2h(ks+2,1,i_uyz) &
1293  * 0.5_rp * ( dens(ks+3,i,j)+dens(ks+3,i+1,j) ) &
1294  + f2h(ks+2,2,i_uyz) &
1295  * 0.5_rp * ( dens(ks+2,i,j)+dens(ks+2,i+1,j) ) )
1296  flux(ks+2,i,j) = j33g * vel &
1297  * ( f61 * ( val(ks+3,i,j)+val(ks+2,i,j) ) &
1298  + f62 * ( val(ks+4,i,j)+val(ks+1,i,j) ) &
1299  + f63 * ( val(ks+5,i,j)+val(ks,i,j) ) ) &
1300  + gsqrt(ks+2,i,j) * num_diff(ks+2,i,j)
1301  vel = ( 0.5_rp * ( mom(ke-3,i,j)+mom(ke-3,i+1,j) ) ) &
1302  / ( f2h(ke-3,1,i_uyz) &
1303  * 0.5_rp * ( dens(ke-2,i,j)+dens(ke-2,i+1,j) ) &
1304  + f2h(ke-3,2,i_uyz) &
1305  * 0.5_rp * ( dens(ke-3,i,j)+dens(ke-3,i+1,j) ) )
1306  flux(ke-3,i,j) = j33g * vel &
1307  * ( f61 * ( val(ke-2,i,j)+val(ke-3,i,j) ) &
1308  + f62 * ( val(ke-1,i,j)+val(ke-4,i,j) ) &
1309  + f63 * ( val(ke,i,j)+val(ke-5,i,j) ) ) &
1310  + gsqrt(ke-3,i,j) * num_diff(ke-3,i,j)
1311 
1312  flux(ke,i,j) = 0.0_rp
1313  enddo
1314  enddo
1315  !$acc end kernels
1316  !$omp end do nowait
1317 
1318  end if
1319 
1320 
1321  !$acc end data
1322 
1323  !$omp end parallel
1324 #ifdef DEBUG
1325  k = iundef; i = iundef; j = iundef
1326 #endif
1327 
1328  return
1329  end subroutine atmos_dyn_fvm_fluxz_uyz_cd8
1330 
1331  !-----------------------------------------------------------------------------
1333  subroutine atmos_dyn_fvm_fluxj13_uyz_cd8( &
1334  flux, &
1335  mom, val, DENS, &
1336  GSQRT, J13G, MAPF, &
1337  CDZ, TwoD, &
1338  IIS, IIE, JJS, JJE )
1339  implicit none
1340 
1341  real(rp), intent(inout) :: flux (ka,ia,ja)
1342  real(rp), intent(in) :: mom (ka,ia,ja)
1343  real(rp), intent(in) :: val (ka,ia,ja)
1344  real(rp), intent(in) :: dens (ka,ia,ja)
1345  real(rp), intent(in) :: gsqrt (ka,ia,ja)
1346  real(rp), intent(in) :: j13g (ka,ia,ja)
1347  real(rp), intent(in) :: mapf ( ia,ja,2)
1348  real(rp), intent(in) :: cdz (ka)
1349  logical, intent(in) :: twod
1350  integer, intent(in) :: iis, iie, jjs, jje
1351 
1352  real(rp) :: vel
1353  integer :: k, i, j
1354  !---------------------------------------------------------------------------
1355 
1356  !$omp parallel default(none) private(i,j,k,vel) &
1357  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J13G,MAPF) &
1358  !$omp shared(GSQRT,CDZ,TwoD)
1359 
1360  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, J13G, MAPF, CDZ)
1361 
1362 
1363 
1364  !$omp do OMP_SCHEDULE_ collapse(2)
1365  !$acc kernels
1366  do j = jjs, jje
1367  do i = iis, iie
1368  do k = ks+3, ke-4
1369  vel = ( f2h(k,1,i_uyz) &
1370  * mom(k+1,i,j) &
1371  + f2h(k,2,i_uyz) &
1372  * mom(k,i,j) ) &
1373  / ( f2h(k,1,i_uyz) &
1374  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i+1,j) ) &
1375  + f2h(k,2,i_uyz) &
1376  * 0.5_rp * ( dens(k,i,j)+dens(k,i+1,j) ) )
1377  vel = vel * j13g(k,i,j)
1378  flux(k,i,j) = vel / mapf(i,j,+2) &
1379  * ( f81 * ( val(k+1,i,j)+val(k,i,j) ) &
1380  + f82 * ( val(k+2,i,j)+val(k-1,i,j) ) &
1381  + f83 * ( val(k+3,i,j)+val(k-2,i,j) ) &
1382  + f84 * ( val(k+4,i,j)+val(k-3,i,j) ) )
1383  enddo
1384  enddo
1385  enddo
1386  !$acc end kernels
1387  !$omp end do nowait
1388 
1389  !$omp do OMP_SCHEDULE_ collapse(2)
1390  !$acc kernels
1391  do j = jjs, jje
1392  do i = iis, iie
1393  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
1394  ! The flux at KS-1 can be non-zero.
1395  ! To reduce calculations, all the fluxes are set to zero.
1396  flux(ks-1,i,j) = 0.0_rp
1397 
1398  vel = ( f2h(ks,1,i_uyz) &
1399  * mom(ks+1,i,j) &
1400  + f2h(ks,2,i_uyz) &
1401  * mom(ks,i,j) ) &
1402  / ( f2h(ks,1,i_uyz) &
1403  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i+1,j) ) &
1404  + f2h(ks,2,i_uyz) &
1405  * 0.5_rp * ( dens(ks,i,j)+dens(ks,i+1,j) ) )
1406  vel = vel * j13g(ks,i,j)
1407  flux(ks,i,j) = vel / mapf(i,j,+2) &
1408  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) )
1409 
1410  vel = ( f2h(ke-1,1,i_uyz) &
1411  * mom(ke,i,j) &
1412  + f2h(ke-1,2,i_uyz) &
1413  * mom(ke-1,i,j) ) &
1414  / ( f2h(ke-1,1,i_uyz) &
1415  * 0.5_rp * ( dens(ke,i,j)+dens(ke,i+1,j) ) &
1416  + f2h(ke-1,2,i_uyz) &
1417  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i+1,j) ) )
1418  vel = vel * j13g(ke-1,i,j)
1419  flux(ke-1,i,j) = vel / mapf(i,j,+2) &
1420  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) )
1421 
1422  vel = ( f2h(ks+1,1,i_uyz) &
1423  * mom(ks+2,i,j) &
1424  + f2h(ks+1,2,i_uyz) &
1425  * mom(ks+1,i,j) ) &
1426  / ( f2h(ks+1,1,i_uyz) &
1427  * 0.5_rp * ( dens(ks+2,i,j)+dens(ks+2,i+1,j) ) &
1428  + f2h(ks+1,2,i_uyz) &
1429  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i+1,j) ) )
1430  vel = vel * j13g(ks+1,i,j)
1431  flux(ks+1,i,j) = vel / mapf(i,j,+2) &
1432  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
1433  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) )
1434 
1435  vel = ( f2h(ke-2,1,i_uyz) &
1436  * mom(ke-1,i,j) &
1437  + f2h(ke-2,2,i_uyz) &
1438  * mom(ke-2,i,j) ) &
1439  / ( f2h(ke-2,1,i_uyz) &
1440  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i+1,j) ) &
1441  + f2h(ke-2,2,i_uyz) &
1442  * 0.5_rp * ( dens(ke-2,i,j)+dens(ke-2,i+1,j) ) )
1443  vel = vel * j13g(ke-2,i,j)
1444  flux(ke-2,i,j) = vel / mapf(i,j,+2) &
1445  * ( f41 * ( val(ke-1,i,j)+val(ke-2,i,j) ) &
1446  + f42 * ( val(ke,i,j)+val(ke-3,i,j) ) )
1447 
1448  vel = ( f2h(ks+2,1,i_uyz) &
1449  * mom(ks+3,i,j) &
1450  + f2h(ks+2,2,i_uyz) &
1451  * mom(ks+2,i,j) ) &
1452  / ( f2h(ks+2,1,i_uyz) &
1453  * 0.5_rp * ( dens(ks+3,i,j)+dens(ks+3,i+1,j) ) &
1454  + f2h(ks+2,2,i_uyz) &
1455  * 0.5_rp * ( dens(ks+2,i,j)+dens(ks+2,i+1,j) ) )
1456  vel = vel * j13g(ks+2,i,j)
1457  flux(ks+2,i,j) = vel / mapf(i,j,+2) &
1458  * ( f61 * ( val(ks+3,i,j)+val(ks+2,i,j) ) &
1459  + f62 * ( val(ks+4,i,j)+val(ks+1,i,j) ) &
1460  + f63 * ( val(ks+5,i,j)+val(ks,i,j) ) )
1461 
1462  vel = ( f2h(ke-3,1,i_uyz) &
1463  * mom(ke-2,i,j) &
1464  + f2h(ke-3,2,i_uyz) &
1465  * mom(ke-3,i,j) ) &
1466  / ( f2h(ke-3,1,i_uyz) &
1467  * 0.5_rp * ( dens(ke-2,i,j)+dens(ke-2,i+1,j) ) &
1468  + f2h(ke-3,2,i_uyz) &
1469  * 0.5_rp * ( dens(ke-3,i,j)+dens(ke-3,i+1,j) ) )
1470  vel = vel * j13g(ke-3,i,j)
1471  flux(ke-3,i,j) = vel / mapf(i,j,+2) &
1472  * ( f61 * ( val(ke-2,i,j)+val(ke-3,i,j) ) &
1473  + f62 * ( val(ke-1,i,j)+val(ke-4,i,j) ) &
1474  + f63 * ( val(ke,i,j)+val(ke-5,i,j) ) )
1475 
1476  flux(ke ,i,j) = 0.0_rp
1477  enddo
1478  enddo
1479  !$acc end kernels
1480  !$omp end do nowait
1481 
1482 
1483 
1484  !$acc end data
1485 
1486  !$omp end parallel
1487  return
1488  end subroutine atmos_dyn_fvm_fluxj13_uyz_cd8
1489 
1490  !-----------------------------------------------------------------------------
1492  subroutine atmos_dyn_fvm_fluxj23_uyz_cd8( &
1493  flux, &
1494  mom, val, DENS, &
1495  GSQRT, J23G, MAPF, &
1496  CDZ, TwoD, &
1497  IIS, IIE, JJS, JJE )
1498  implicit none
1499 
1500  real(rp), intent(inout) :: flux (ka,ia,ja)
1501  real(rp), intent(in) :: mom (ka,ia,ja)
1502  real(rp), intent(in) :: val (ka,ia,ja)
1503  real(rp), intent(in) :: dens (ka,ia,ja)
1504  real(rp), intent(in) :: gsqrt (ka,ia,ja)
1505  real(rp), intent(in) :: j23g (ka,ia,ja)
1506  real(rp), intent(in) :: mapf ( ia,ja,2)
1507  real(rp), intent(in) :: cdz (ka)
1508  logical, intent(in) :: twod
1509  integer, intent(in) :: iis, iie, jjs, jje
1510 
1511  real(rp) :: vel
1512  integer :: k, i, j
1513  !---------------------------------------------------------------------------
1514 
1515  !$omp parallel default(none) private(i,j,k,vel) &
1516  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J23G,MAPF) &
1517  !$omp shared(GSQRT,CDZ,TwoD)
1518 
1519  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, J23G, MAPF, CDZ)
1520 
1521 
1522  if ( twod ) then
1523 
1524  !$omp do OMP_SCHEDULE_
1525  !$acc kernels
1526  do j = jjs, jje
1527  do k = ks+3, ke-4
1528  i = iis
1529  vel = ( f2h(k,1,i_xyz) &
1530  * 0.5_rp * ( mom(k+1,i,j)+mom(k+1,i,j-1) ) &
1531  + f2h(k,2,i_xyz) &
1532  * 0.5_rp * ( mom(k,i,j)+mom(k,i,j-1) ) ) &
1533  / ( f2h(k,1,i_xyz) &
1534  * dens(k+1,i,j) &
1535  + f2h(k,2,i_xyz) &
1536  * dens(k,i,j) )
1537  vel = vel * j23g(k,i,j)
1538  flux(k,i,j) = vel * ( f81 * ( val(k+1,i,j)+val(k,i,j) ) &
1539  + f82 * ( val(k+2,i,j)+val(k-1,i,j) ) &
1540  + f83 * ( val(k+3,i,j)+val(k-2,i,j) ) &
1541  + f84 * ( val(k+4,i,j)+val(k-3,i,j) ) )
1542  enddo
1543  enddo
1544  !$acc end kernels
1545  !$omp end do nowait
1546 
1547  !$omp do OMP_SCHEDULE_
1548  !$acc kernels
1549  do j = jjs, jje
1550  i = iis
1551  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
1552  ! The flux at KS-1 can be non-zero.
1553  ! To reduce calculations, all the fluxes are set to zero.
1554  flux(ks-1,i,j) = 0.0_rp
1555 
1556  vel = ( f2h(ks,1,i_xyz) &
1557  * 0.5_rp * ( mom(ks+1,i,j)+mom(ks+1,i,j-1) ) &
1558  + f2h(ks,2,i_xyz) &
1559  * 0.5_rp * ( mom(ks,i,j)+mom(ks,i,j-1) ) ) &
1560  / ( f2h(ks,1,i_xyz) &
1561  * dens(ks+1,i,j) &
1562  + f2h(ks,2,i_xyz) &
1563  * dens(ks,i,j) )
1564  vel = vel * j23g(ks,i,j)
1565  flux(ks,i,j) = vel / mapf(i,j,+1) &
1566  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) )
1567 
1568  vel = ( f2h(ke-1,1,i_xyz) &
1569  * 0.5_rp * ( mom(ke,i,j)+mom(ke,i,j-1) ) &
1570  + f2h(ke-1,2,i_xyz) &
1571  * 0.5_rp * ( mom(ke-1,i,j)+mom(ke-1,i,j-1) ) ) &
1572  / ( f2h(ke-1,1,i_xyz) &
1573  * dens(ke,i,j) &
1574  + f2h(ke-1,2,i_xyz) &
1575  * dens(ke-1,i,j) )
1576  vel = vel * j23g(ke-1,i,j)
1577  flux(ke-1,i,j) = vel / mapf(i,j,+1) &
1578  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) )
1579 
1580  vel = ( f2h(ks+1,1,i_xyz) &
1581  * 0.5_rp * ( mom(ks+2,i,j)+mom(ks+2,i,j-1) ) &
1582  + f2h(ks+1,2,i_xyz) &
1583  * 0.5_rp * ( mom(ks+1,i,j)+mom(ks+1,i,j-1) ) ) &
1584  / ( f2h(ks+1,1,i_xyz) &
1585  * dens(ks+2,i,j) &
1586  + f2h(ks+1,2,i_xyz) &
1587  * dens(ks+1,i,j) )
1588  vel = vel * j23g(ks+1,i,j)
1589  flux(ks+1,i,j) = vel / mapf(i,j,+1) &
1590  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
1591  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) )
1592 
1593  vel = ( f2h(ke-2,1,i_xyz) &
1594  * 0.5_rp * ( mom(ke-1,i,j)+mom(ke-1,i,j-1) ) &
1595  + f2h(ke-2,2,i_xyz) &
1596  * 0.5_rp * ( mom(ke-2,i,j)+mom(ke-2,i,j-1) ) ) &
1597  / ( f2h(ke-2,1,i_xyz) &
1598  * dens(ke-1,i,j) &
1599  + f2h(ke-2,2,i_xyz) &
1600  * dens(ke-2,i,j) )
1601  vel = vel * j23g(ke-2,i,j)
1602  flux(ke-2,i,j) = vel / mapf(i,j,+1) &
1603  * ( f41 * ( val(ke-1,i,j)+val(ke-2,i,j) ) &
1604  + f42 * ( val(ke,i,j)+val(ke-3,i,j) ) )
1605 
1606  vel = ( f2h(ks+2,1,i_xyz) &
1607  * 0.5_rp * ( mom(ks+3,i,j)+mom(ks+3,i,j-1) ) &
1608  + f2h(ks+2,2,i_xyz) &
1609  * 0.5_rp * ( mom(ks+2,i,j)+mom(ks+2,i,j-1) ) ) &
1610  / ( f2h(ks+2,1,i_xyz) &
1611  * dens(ks+3,i,j) &
1612  + f2h(ks+2,2,i_xyz) &
1613  * dens(ks+2,i,j) )
1614  vel = vel * j23g(ks+2,i,j)
1615  flux(ks+2,i,j) = vel / mapf(i,j,+1) &
1616  * ( f61 * ( val(ks+3,i,j)+val(ks+2,i,j) ) &
1617  + f62 * ( val(ks+4,i,j)+val(ks+1,i,j) ) &
1618  + f63 * ( val(ks+5,i,j)+val(ks,i,j) ) )
1619 
1620  vel = ( f2h(ke-3,1,i_xyz) &
1621  * 0.5_rp * ( mom(ke-2,i,j)+mom(ke-2,i,j-1) ) &
1622  + f2h(ke-3,2,i_xyz) &
1623  * 0.5_rp * ( mom(ke-3,i,j)+mom(ke-3,i,j-1) ) ) &
1624  / ( f2h(ke-3,1,i_xyz) &
1625  * dens(ke-2,i,j) &
1626  + f2h(ke-3,2,i_xyz) &
1627  * dens(ke-3,i,j) )
1628  vel = vel * j23g(ke-3,i,j)
1629  flux(ke-3,i,j) = vel / mapf(i,j,+1) &
1630  * ( f61 * ( val(ke-2,i,j)+val(ke-3,i,j) ) &
1631  + f62 * ( val(ke-1,i,j)+val(ke-4,i,j) ) &
1632  + f63 * ( val(ke,i,j)+val(ke-5,i,j) ) )
1633 
1634  flux(ke ,i,j) = 0.0_rp
1635  enddo
1636  !$acc end kernels
1637  !$omp end do nowait
1638 
1639  else
1640 
1641 
1642  !$omp do OMP_SCHEDULE_ collapse(2)
1643  !$acc kernels
1644  do j = jjs, jje
1645  do i = iis, iie
1646  do k = ks+3, ke-4
1647  vel = ( f2h(k,1,i_uyz) &
1648  * 0.25_rp * ( mom(k+1,i,j)+mom(k+1,i+1,j)+mom(k+1,i,j-1)+mom(k+1,i+1,j-1) ) &
1649  + f2h(k,2,i_uyz) &
1650  * 0.25_rp * ( mom(k,i,j)+mom(k,i+1,j)+mom(k,i,j-1)+mom(k,i+1,j-1) ) ) &
1651  / ( f2h(k,1,i_uyz) &
1652  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i+1,j) ) &
1653  + f2h(k,2,i_uyz) &
1654  * 0.5_rp * ( dens(k,i,j)+dens(k,i+1,j) ) )
1655  vel = vel * j23g(k,i,j)
1656  flux(k,i,j) = vel / mapf(i,j,+1) &
1657  * ( f81 * ( val(k+1,i,j)+val(k,i,j) ) &
1658  + f82 * ( val(k+2,i,j)+val(k-1,i,j) ) &
1659  + f83 * ( val(k+3,i,j)+val(k-2,i,j) ) &
1660  + f84 * ( val(k+4,i,j)+val(k-3,i,j) ) )
1661  enddo
1662  enddo
1663  enddo
1664  !$acc end kernels
1665  !$omp end do nowait
1666 
1667  !$omp do OMP_SCHEDULE_ collapse(2)
1668  !$acc kernels
1669  do j = jjs, jje
1670  do i = iis, iie
1671  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
1672  ! The flux at KS-1 can be non-zero.
1673  ! To reduce calculations, all the fluxes are set to zero.
1674  flux(ks-1,i,j) = 0.0_rp
1675 
1676  vel = ( f2h(ks,1,i_uyz) &
1677  * 0.25_rp * ( mom(ks+1,i,j)+mom(ks+1,i+1,j)+mom(ks+1,i,j-1)+mom(ks+1,i+1,j-1) ) &
1678  + f2h(ks,2,i_uyz) &
1679  * 0.25_rp * ( mom(ks,i,j)+mom(ks,i+1,j)+mom(ks,i,j-1)+mom(ks,i+1,j-1) ) ) &
1680  / ( f2h(ks,1,i_uyz) &
1681  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i+1,j) ) &
1682  + f2h(ks,2,i_uyz) &
1683  * 0.5_rp * ( dens(ks,i,j)+dens(ks,i+1,j) ) )
1684  vel = vel * j23g(ks,i,j)
1685  flux(ks,i,j) = vel / mapf(i,j,+1) &
1686  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) )
1687 
1688  vel = ( f2h(ke-1,1,i_uyz) &
1689  * 0.25_rp * ( mom(ke,i,j)+mom(ke,i+1,j)+mom(ke,i,j-1)+mom(ke,i+1,j-1) ) &
1690  + f2h(ke-1,2,i_uyz) &
1691  * 0.25_rp * ( mom(ke-1,i,j)+mom(ke-1,i+1,j)+mom(ke-1,i,j-1)+mom(ke-1,i+1,j-1) ) ) &
1692  / ( f2h(ke-1,1,i_uyz) &
1693  * 0.5_rp * ( dens(ke,i,j)+dens(ke,i+1,j) ) &
1694  + f2h(ke-1,2,i_uyz) &
1695  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i+1,j) ) )
1696  vel = vel * j23g(ke-1,i,j)
1697  flux(ke-1,i,j) = vel / mapf(i,j,+1) &
1698  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) )
1699 
1700  vel = ( f2h(ks+1,1,i_uyz) &
1701  * 0.25_rp * ( mom(ks+2,i,j)+mom(ks+2,i+1,j)+mom(ks+2,i,j-1)+mom(ks+2,i+1,j-1) ) &
1702  + f2h(ks+1,2,i_uyz) &
1703  * 0.25_rp * ( mom(ks+1,i,j)+mom(ks+1,i+1,j)+mom(ks+1,i,j-1)+mom(ks+1,i+1,j-1) ) ) &
1704  / ( f2h(ks+1,1,i_uyz) &
1705  * 0.5_rp * ( dens(ks+2,i,j)+dens(ks+2,i+1,j) ) &
1706  + f2h(ks+1,2,i_uyz) &
1707  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i+1,j) ) )
1708  vel = vel * j23g(ks+1,i,j)
1709  flux(ks+1,i,j) = vel / mapf(i,j,+1) &
1710  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
1711  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) )
1712 
1713  vel = ( f2h(ke-2,1,i_uyz) &
1714  * 0.25_rp * ( mom(ke-1,i,j)+mom(ke-1,i+1,j)+mom(ke-1,i,j-1)+mom(ke-1,i+1,j-1) ) &
1715  + f2h(ke-2,2,i_uyz) &
1716  * 0.25_rp * ( mom(ke-2,i,j)+mom(ke-2,i+1,j)+mom(ke-2,i,j-1)+mom(ke-2,i+1,j-1) ) ) &
1717  / ( f2h(ke-2,1,i_uyz) &
1718  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i+1,j) ) &
1719  + f2h(ke-2,2,i_uyz) &
1720  * 0.5_rp * ( dens(ke-2,i,j)+dens(ke-2,i+1,j) ) )
1721  vel = vel * j23g(ke-2,i,j)
1722  flux(ke-2,i,j) = vel / mapf(i,j,+1) &
1723  * ( f41 * ( val(ke-1,i,j)+val(ke-2,i,j) ) &
1724  + f42 * ( val(ke,i,j)+val(ke-3,i,j) ) )
1725 
1726  vel = ( f2h(ks+2,1,i_uyz) &
1727  * 0.25_rp * ( mom(ks+3,i,j)+mom(ks+3,i+1,j)+mom(ks+3,i,j-1)+mom(ks+3,i+1,j-1) ) &
1728  + f2h(ks+2,2,i_uyz) &
1729  * 0.25_rp * ( mom(ks+2,i,j)+mom(ks+2,i+1,j)+mom(ks+2,i,j-1)+mom(ks+2,i+1,j-1) ) ) &
1730  / ( f2h(ks+2,1,i_uyz) &
1731  * 0.5_rp * ( dens(ks+3,i,j)+dens(ks+3,i+1,j) ) &
1732  + f2h(ks+2,2,i_uyz) &
1733  * 0.5_rp * ( dens(ks+2,i,j)+dens(ks+2,i+1,j) ) )
1734  vel = vel * j23g(ks+2,i,j)
1735  flux(ks+2,i,j) = vel / mapf(i,j,+1) &
1736  * ( f61 * ( val(ks+3,i,j)+val(ks+2,i,j) ) &
1737  + f62 * ( val(ks+4,i,j)+val(ks+1,i,j) ) &
1738  + f63 * ( val(ks+5,i,j)+val(ks,i,j) ) )
1739 
1740  vel = ( f2h(ke-3,1,i_uyz) &
1741  * 0.25_rp * ( mom(ke-2,i,j)+mom(ke-2,i+1,j)+mom(ke-2,i,j-1)+mom(ke-2,i+1,j-1) ) &
1742  + f2h(ke-3,2,i_uyz) &
1743  * 0.25_rp * ( mom(ke-3,i,j)+mom(ke-3,i+1,j)+mom(ke-3,i,j-1)+mom(ke-3,i+1,j-1) ) ) &
1744  / ( f2h(ke-3,1,i_uyz) &
1745  * 0.5_rp * ( dens(ke-2,i,j)+dens(ke-2,i+1,j) ) &
1746  + f2h(ke-3,2,i_uyz) &
1747  * 0.5_rp * ( dens(ke-3,i,j)+dens(ke-3,i+1,j) ) )
1748  vel = vel * j23g(ke-3,i,j)
1749  flux(ke-3,i,j) = vel / mapf(i,j,+1) &
1750  * ( f61 * ( val(ke-2,i,j)+val(ke-3,i,j) ) &
1751  + f62 * ( val(ke-1,i,j)+val(ke-4,i,j) ) &
1752  + f63 * ( val(ke,i,j)+val(ke-5,i,j) ) )
1753 
1754  flux(ke ,i,j) = 0.0_rp
1755  enddo
1756  enddo
1757  !$acc end kernels
1758  !$omp end do nowait
1759 
1760 
1761  end if
1762 
1763 
1764  !$acc end data
1765 
1766  !$omp end parallel
1767  return
1768  end subroutine atmos_dyn_fvm_fluxj23_uyz_cd8
1769 
1770  !-----------------------------------------------------------------------------
1772  subroutine atmos_dyn_fvm_fluxx_uyz_cd8( &
1773  flux, &
1774  mom, val, DENS, &
1775  GSQRT, MAPF, &
1776  num_diff, &
1777  CDZ, TwoD, &
1778  IIS, IIE, JJS, JJE )
1779  implicit none
1780 
1781  real(rp), intent(inout) :: flux (ka,ia,ja)
1782  real(rp), intent(in) :: mom (ka,ia,ja)
1783  real(rp), intent(in) :: val (ka,ia,ja)
1784  real(rp), intent(in) :: dens (ka,ia,ja)
1785  real(rp), intent(in) :: gsqrt (ka,ia,ja)
1786  real(rp), intent(in) :: mapf ( ia,ja,2)
1787  real(rp), intent(in) :: num_diff(ka,ia,ja)
1788  real(rp), intent(in) :: cdz (ka)
1789  logical, intent(in) :: twod
1790  integer, intent(in) :: iis, iie, jjs, jje
1791 
1792  real(rp) :: vel
1793  integer :: k, i, j
1794  !---------------------------------------------------------------------------
1795 
1796  ! note that x-index is added by -1
1797 
1798 
1799 
1800  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
1801  !$omp private(vel) &
1802  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff)
1803  !$acc kernels
1804  do j = jjs, jje
1805  do i = iis, iie+1
1806  do k = ks, ke
1807 #ifdef DEBUG
1808  call check( __line__, mom(k,i ,j) )
1809  call check( __line__, mom(k,i-1,j) )
1810 
1811  call check( __line__, val(k,i-1,j) )
1812  call check( __line__, val(k,i,j) )
1813 
1814  call check( __line__, val(k,i-2,j) )
1815  call check( __line__, val(k,i+1,j) )
1816 
1817  call check( __line__, val(k,i-3,j) )
1818  call check( __line__, val(k,i+2,j) )
1819 
1820  call check( __line__, val(k,i-4,j) )
1821  call check( __line__, val(k,i+3,j) )
1822 
1823 #endif
1824  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i-1,j) ) ) &
1825  / ( dens(k,i,j) )
1826  flux(k,i-1,j) = gsqrt(k,i,j) / mapf(i,j,+2) * vel &
1827  * ( f81 * ( val(k,i,j)+val(k,i-1,j) ) &
1828  + f82 * ( val(k,i+1,j)+val(k,i-2,j) ) &
1829  + f83 * ( val(k,i+2,j)+val(k,i-3,j) ) &
1830  + f84 * ( val(k,i+3,j)+val(k,i-4,j) ) ) &
1831  + gsqrt(k,i,j) * num_diff(k,i,j)
1832  enddo
1833  enddo
1834  enddo
1835  !$acc end kernels
1836 #ifdef DEBUG
1837  k = iundef; i = iundef; j = iundef
1838 #endif
1839 
1840 
1841 
1842  return
1843  end subroutine atmos_dyn_fvm_fluxx_uyz_cd8
1844 
1845  !-----------------------------------------------------------------------------
1847  subroutine atmos_dyn_fvm_fluxy_uyz_cd8( &
1848  flux, &
1849  mom, val, DENS, &
1850  GSQRT, MAPF, &
1851  num_diff, &
1852  CDZ, TwoD, &
1853  IIS, IIE, JJS, JJE )
1854  implicit none
1855 
1856  real(rp), intent(inout) :: flux (ka,ia,ja)
1857  real(rp), intent(in) :: mom (ka,ia,ja)
1858  real(rp), intent(in) :: val (ka,ia,ja)
1859  real(rp), intent(in) :: dens (ka,ia,ja)
1860  real(rp), intent(in) :: gsqrt (ka,ia,ja)
1861  real(rp), intent(in) :: mapf ( ia,ja,2)
1862  real(rp), intent(in) :: num_diff(ka,ia,ja)
1863  real(rp), intent(in) :: cdz (ka)
1864  logical, intent(in) :: twod
1865  integer, intent(in) :: iis, iie, jjs, jje
1866 
1867  real(rp) :: vel
1868  integer :: k, i, j
1869  !---------------------------------------------------------------------------
1870 
1871 
1872 
1873  if ( twod ) then
1874 
1875  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ &
1876  !$omp private(vel) &
1877  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff,TwoD)
1878  !$acc kernels
1879  do j = jjs-1, jje
1880  do k = ks, ke
1881  i = iis
1882 #ifdef DEBUG
1883  call check( __line__, mom(k,i ,j) )
1884 
1885  call check( __line__, val(k,i,j) )
1886  call check( __line__, val(k,i,j+1) )
1887 
1888  call check( __line__, val(k,i,j-1) )
1889  call check( __line__, val(k,i,j+2) )
1890 
1891  call check( __line__, val(k,i,j-2) )
1892  call check( __line__, val(k,i,j+3) )
1893 
1894  call check( __line__, val(k,i,j-3) )
1895  call check( __line__, val(k,i,j+4) )
1896 
1897 #endif
1898  vel = ( mom(k,i,j) ) &
1899  / ( 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) ) )
1900  flux(k,i,j) = gsqrt(k,i,j) / mapf(i,j,+1) * vel &
1901  * ( f81 * ( val(k,i,j+1)+val(k,i,j) ) &
1902  + f82 * ( val(k,i,j+2)+val(k,i,j-1) ) &
1903  + f83 * ( val(k,i,j+3)+val(k,i,j-2) ) &
1904  + f84 * ( val(k,i,j+4)+val(k,i,j-3) ) ) &
1905  + gsqrt(k,i,j) * num_diff(k,i,j)
1906  enddo
1907  enddo
1908  !$acc end kernels
1909 #ifdef DEBUG
1910  k = iundef; i = iundef; j = iundef
1911 #endif
1912 
1913  else
1914 
1915 
1916  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
1917  !$omp private(vel) &
1918  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff)
1919  !$acc kernels
1920  do j = jjs-1, jje
1921  do i = iis, iie
1922  do k = ks, ke
1923 #ifdef DEBUG
1924  call check( __line__, mom(k,i ,j) )
1925  call check( __line__, mom(k,i-1,j) )
1926 
1927  call check( __line__, val(k,i,j) )
1928  call check( __line__, val(k,i,j+1) )
1929 
1930  call check( __line__, val(k,i,j-1) )
1931  call check( __line__, val(k,i,j+2) )
1932 
1933  call check( __line__, val(k,i,j-2) )
1934  call check( __line__, val(k,i,j+3) )
1935 
1936  call check( __line__, val(k,i,j-3) )
1937  call check( __line__, val(k,i,j+4) )
1938 
1939 #endif
1940  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i+1,j) ) ) &
1941  / ( 0.25_rp * ( dens(k,i,j)+dens(k,i+1,j)+dens(k,i,j+1)+dens(k,i+1,j+1) ) )
1942  flux(k,i,j) = gsqrt(k,i,j) / mapf(i,j,+1) * vel &
1943  * ( f81 * ( val(k,i,j+1)+val(k,i,j) ) &
1944  + f82 * ( val(k,i,j+2)+val(k,i,j-1) ) &
1945  + f83 * ( val(k,i,j+3)+val(k,i,j-2) ) &
1946  + f84 * ( val(k,i,j+4)+val(k,i,j-3) ) ) &
1947  + gsqrt(k,i,j) * num_diff(k,i,j)
1948  enddo
1949  enddo
1950  enddo
1951  !$acc end kernels
1952 #ifdef DEBUG
1953  k = iundef; i = iundef; j = iundef
1954 #endif
1955 
1956 
1957  end if
1958 
1959 
1960  return
1961  end subroutine atmos_dyn_fvm_fluxy_uyz_cd8
1962 
1963 
1964 
1965  !-----------------------------------------------------------------------------
1967  subroutine atmos_dyn_fvm_fluxz_xvz_cd8( &
1968  flux, &
1969  mom, val, DENS, &
1970  GSQRT, J33G, &
1971  num_diff, &
1972  CDZ, TwoD, &
1973  IIS, IIE, JJS, JJE )
1974  implicit none
1975 
1976  real(rp), intent(inout) :: flux (ka,ia,ja)
1977  real(rp), intent(in) :: mom (ka,ia,ja)
1978  real(rp), intent(in) :: val (ka,ia,ja)
1979  real(rp), intent(in) :: dens (ka,ia,ja)
1980  real(rp), intent(in) :: gsqrt (ka,ia,ja)
1981  real(rp), intent(in) :: j33g
1982  real(rp), intent(in) :: num_diff(ka,ia,ja)
1983  real(rp), intent(in) :: cdz (ka)
1984  logical, intent(in) :: twod
1985  integer, intent(in) :: iis, iie, jjs, jje
1986 
1987  real(rp) :: vel
1988  integer :: k, i, j
1989  !---------------------------------------------------------------------------
1990 
1991  !$omp parallel default(none) private(i,j,k,vel) &
1992  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J33G,GSQRT,num_diff) &
1993  !$omp shared(CDZ,TwoD)
1994 
1995  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, num_diff, CDZ)
1996 
1997 
1998  !$omp do OMP_SCHEDULE_ collapse(2)
1999  !$acc kernels
2000  do j = jjs, jje
2001  do i = iis, iie
2002  do k = ks+3, ke-4
2003 #ifdef DEBUG
2004  call check( __line__, mom(k,i,j) )
2005  call check( __line__, mom(k,i,j+1) )
2006 
2007  call check( __line__, val(k,i,j) )
2008  call check( __line__, val(k+1,i,j) )
2009 
2010  call check( __line__, val(k-1,i,j) )
2011  call check( __line__, val(k+2,i,j) )
2012 
2013  call check( __line__, val(k-2,i,j) )
2014  call check( __line__, val(k+3,i,j) )
2015 
2016  call check( __line__, val(k-3,i,j) )
2017  call check( __line__, val(k+4,i,j) )
2018 
2019 #endif
2020  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i,j+1) ) ) &
2021  / ( f2h(k,1,i_xvz) &
2022  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i,j+1) ) &
2023  + f2h(k,2,i_xvz) &
2024  * 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) ) )
2025  flux(k,i,j) = j33g * vel &
2026  * ( f81 * ( val(k+1,i,j)+val(k,i,j) ) &
2027  + f82 * ( val(k+2,i,j)+val(k-1,i,j) ) &
2028  + f83 * ( val(k+3,i,j)+val(k-2,i,j) ) &
2029  + f84 * ( val(k+4,i,j)+val(k-3,i,j) ) ) &
2030  + gsqrt(k,i,j) * num_diff(k,i,j)
2031  enddo
2032  enddo
2033  enddo
2034  !$acc end kernels
2035  !$omp end do nowait
2036 #ifdef DEBUG
2037  k = iundef; i = iundef; j = iundef
2038 #endif
2039 
2040  !$omp do OMP_SCHEDULE_ collapse(2)
2041  !$acc kernels
2042  do j = jjs, jje
2043  do i = iis, iie
2044 #ifdef DEBUG
2045 
2046  call check( __line__, mom(ks,i ,j) )
2047  call check( __line__, mom(ks,i,j+1) )
2048  call check( __line__, val(ks+1,i,j) )
2049  call check( __line__, val(ks,i,j) )
2050 
2051  call check( __line__, mom(ks+1,i ,j) )
2052  call check( __line__, mom(ks+1,i,j+1) )
2053  call check( __line__, val(ks+3,i,j) )
2054  call check( __line__, val(ks+2,i,j) )
2055 
2056  call check( __line__, mom(ks+2,i ,j) )
2057  call check( __line__, mom(ks+2,i,j+1) )
2058  call check( __line__, val(ks+5,i,j) )
2059  call check( __line__, val(ks+4,i,j) )
2060 
2061 #endif
2062  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
2063  ! The flux at KS-1 can be non-zero.
2064  ! To reduce calculations, all the fluxes are set to zero.
2065  flux(ks-1,i,j) = 0.0_rp
2066 
2067  vel = ( 0.5_rp * ( mom(ks,i,j)+mom(ks,i,j+1) ) ) &
2068  / ( f2h(ks,1,i_xvz) &
2069  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i,j+1) ) &
2070  + f2h(ks,2,i_xvz) &
2071  * 0.5_rp * ( dens(ks,i,j)+dens(ks,i,j+1) ) )
2072  flux(ks,i,j) = j33g * vel &
2073  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) ) &
2074  + gsqrt(ks,i,j) * num_diff(ks,i,j)
2075  vel = ( 0.5_rp * ( mom(ke-1,i,j)+mom(ke-1,i,j+1) ) ) &
2076  / ( f2h(ke-1,1,i_xvz) &
2077  * 0.5_rp * ( dens(ke,i,j)+dens(ke,i,j+1) ) &
2078  + f2h(ke-1,2,i_xvz) &
2079  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i,j+1) ) )
2080  flux(ke-1,i,j) = j33g * vel &
2081  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) ) &
2082  + gsqrt(ke-1,i,j) * num_diff(ke-1,i,j)
2083 
2084  vel = ( 0.5_rp * ( mom(ks+1,i,j)+mom(ks+1,i,j+1) ) ) &
2085  / ( f2h(ks+1,1,i_xvz) &
2086  * 0.5_rp * ( dens(ks+2,i,j)+dens(ks+2,i,j+1) ) &
2087  + f2h(ks+1,2,i_xvz) &
2088  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i,j+1) ) )
2089  flux(ks+1,i,j) = j33g * vel &
2090  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
2091  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) ) &
2092  + gsqrt(ks+1,i,j) * num_diff(ks+1,i,j)
2093  vel = ( 0.5_rp * ( mom(ke-2,i,j)+mom(ke-2,i,j+1) ) ) &
2094  / ( f2h(ke-2,1,i_xvz) &
2095  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i,j+1) ) &
2096  + f2h(ke-2,2,i_xvz) &
2097  * 0.5_rp * ( dens(ke-2,i,j)+dens(ke-2,i,j+1) ) )
2098  flux(ke-2,i,j) = j33g * vel &
2099  * ( f41 * ( val(ke-1,i,j)+val(ke-2,i,j) ) &
2100  + f42 * ( val(ke,i,j)+val(ke-3,i,j) ) ) &
2101  + gsqrt(ke-2,i,j) * num_diff(ke-2,i,j)
2102 
2103  vel = ( 0.5_rp * ( mom(ks+2,i,j)+mom(ks+2,i,j+1) ) ) &
2104  / ( f2h(ks+2,1,i_xvz) &
2105  * 0.5_rp * ( dens(ks+3,i,j)+dens(ks+3,i,j+1) ) &
2106  + f2h(ks+2,2,i_xvz) &
2107  * 0.5_rp * ( dens(ks+2,i,j)+dens(ks+2,i,j+1) ) )
2108  flux(ks+2,i,j) = j33g * vel &
2109  * ( f61 * ( val(ks+3,i,j)+val(ks+2,i,j) ) &
2110  + f62 * ( val(ks+4,i,j)+val(ks+1,i,j) ) &
2111  + f63 * ( val(ks+5,i,j)+val(ks,i,j) ) ) &
2112  + gsqrt(ks+2,i,j) * num_diff(ks+2,i,j)
2113  vel = ( 0.5_rp * ( mom(ke-3,i,j)+mom(ke-3,i,j+1) ) ) &
2114  / ( f2h(ke-3,1,i_xvz) &
2115  * 0.5_rp * ( dens(ke-2,i,j)+dens(ke-2,i,j+1) ) &
2116  + f2h(ke-3,2,i_xvz) &
2117  * 0.5_rp * ( dens(ke-3,i,j)+dens(ke-3,i,j+1) ) )
2118  flux(ke-3,i,j) = j33g * vel &
2119  * ( f61 * ( val(ke-2,i,j)+val(ke-3,i,j) ) &
2120  + f62 * ( val(ke-1,i,j)+val(ke-4,i,j) ) &
2121  + f63 * ( val(ke,i,j)+val(ke-5,i,j) ) ) &
2122  + gsqrt(ke-3,i,j) * num_diff(ke-3,i,j)
2123 
2124  flux(ke,i,j) = 0.0_rp
2125  enddo
2126  enddo
2127  !$acc end kernels
2128  !$omp end do nowait
2129 
2130 
2131  !$acc end data
2132 
2133  !$omp end parallel
2134 #ifdef DEBUG
2135  k = iundef; i = iundef; j = iundef
2136 #endif
2137 
2138  return
2139  end subroutine atmos_dyn_fvm_fluxz_xvz_cd8
2140 
2141  !-----------------------------------------------------------------------------
2143  subroutine atmos_dyn_fvm_fluxj13_xvz_cd8( &
2144  flux, &
2145  mom, val, DENS, &
2146  GSQRT, J13G, MAPF, &
2147  CDZ, TwoD, &
2148  IIS, IIE, JJS, JJE )
2149  implicit none
2150 
2151  real(rp), intent(inout) :: flux (ka,ia,ja)
2152  real(rp), intent(in) :: mom (ka,ia,ja)
2153  real(rp), intent(in) :: val (ka,ia,ja)
2154  real(rp), intent(in) :: dens (ka,ia,ja)
2155  real(rp), intent(in) :: gsqrt (ka,ia,ja)
2156  real(rp), intent(in) :: j13g (ka,ia,ja)
2157  real(rp), intent(in) :: mapf ( ia,ja,2)
2158  real(rp), intent(in) :: cdz (ka)
2159  logical, intent(in) :: twod
2160  integer, intent(in) :: iis, iie, jjs, jje
2161 
2162  real(rp) :: vel
2163  integer :: k, i, j
2164  !---------------------------------------------------------------------------
2165 
2166  !$omp parallel default(none) private(i,j,k,vel) &
2167  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J13G,MAPF) &
2168  !$omp shared(GSQRT,CDZ,TwoD)
2169 
2170  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, J13G, MAPF, CDZ)
2171 
2172 
2173 
2174  !$omp do OMP_SCHEDULE_ collapse(2)
2175  !$acc kernels
2176  do j = jjs, jje
2177  do i = iis, iie
2178  do k = ks+3, ke-4
2179  vel = ( f2h(k,1,i_xvz) &
2180  * 0.25_rp * ( mom(k+1,i,j)+mom(k+1,i-1,j)+mom(k+1,i,j+1)+mom(k+1,i-1,j+1) ) &
2181  + f2h(k,2,i_xvz) &
2182  * 0.25_rp * ( mom(k,i,j)+mom(k,i-1,j)+mom(k,i,j+1)+mom(k,i-1,j+1) ) ) &
2183  / ( f2h(k,1,i_xvz) &
2184  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i,j+1) ) &
2185  + f2h(k,2,i_xvz) &
2186  * 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) ) )
2187  vel = vel * j13g(k,i,j)
2188  flux(k,i,j) = vel / mapf(i,j,+2) &
2189  * ( f81 * ( val(k+1,i,j)+val(k,i,j) ) &
2190  + f82 * ( val(k+2,i,j)+val(k-1,i,j) ) &
2191  + f83 * ( val(k+3,i,j)+val(k-2,i,j) ) &
2192  + f84 * ( val(k+4,i,j)+val(k-3,i,j) ) )
2193  enddo
2194  enddo
2195  enddo
2196  !$acc end kernels
2197  !$omp end do nowait
2198 
2199  !$omp do OMP_SCHEDULE_ collapse(2)
2200  !$acc kernels
2201  do j = jjs, jje
2202  do i = iis, iie
2203  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
2204  ! The flux at KS-1 can be non-zero.
2205  ! To reduce calculations, all the fluxes are set to zero.
2206  flux(ks-1,i,j) = 0.0_rp
2207 
2208  vel = ( f2h(ks,1,i_xvz) &
2209  * 0.25_rp * ( mom(ks+1,i,j)+mom(ks+1,i-1,j)+mom(ks+1,i,j+1)+mom(ks+1,i-1,j+1) ) &
2210  + f2h(ks,2,i_xvz) &
2211  * 0.25_rp * ( mom(ks,i,j)+mom(ks,i-1,j)+mom(ks,i,j+1)+mom(ks,i-1,j+1) ) ) &
2212  / ( f2h(ks,1,i_xvz) &
2213  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i,j+1) ) &
2214  + f2h(ks,2,i_xvz) &
2215  * 0.5_rp * ( dens(ks,i,j)+dens(ks,i,j+1) ) )
2216  vel = vel * j13g(ks,i,j)
2217  flux(ks,i,j) = vel / mapf(i,j,+2) &
2218  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) )
2219 
2220  vel = ( f2h(ke-1,1,i_xvz) &
2221  * 0.25_rp * ( mom(ke,i,j)+mom(ke,i-1,j)+mom(ke,i,j+1)+mom(ke,i-1,j+1) ) &
2222  + f2h(ke-1,2,i_xvz) &
2223  * 0.25_rp * ( mom(ke-1,i,j)+mom(ke-1,i-1,j)+mom(ke-1,i,j+1)+mom(ke-1,i-1,j+1) ) ) &
2224  / ( f2h(ke-1,1,i_xvz) &
2225  * 0.5_rp * ( dens(ke,i,j)+dens(ke,i,j+1) ) &
2226  + f2h(ke-1,2,i_xvz) &
2227  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i,j+1) ) )
2228  vel = vel * j13g(ke-1,i,j)
2229  flux(ke-1,i,j) = vel / mapf(i,j,+2) &
2230  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) )
2231 
2232  vel = ( f2h(ks+1,1,i_xvz) &
2233  * 0.25_rp * ( mom(ks+2,i,j)+mom(ks+2,i-1,j)+mom(ks+2,i,j+1)+mom(ks+2,i-1,j+1) ) &
2234  + f2h(ks+1,2,i_xvz) &
2235  * 0.25_rp * ( mom(ks+1,i,j)+mom(ks+1,i-1,j)+mom(ks+1,i,j+1)+mom(ks+1,i-1,j+1) ) ) &
2236  / ( f2h(ks+1,1,i_xvz) &
2237  * 0.5_rp * ( dens(ks+2,i,j)+dens(ks+2,i,j+1) ) &
2238  + f2h(ks+1,2,i_xvz) &
2239  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i,j+1) ) )
2240  vel = vel * j13g(ks+1,i,j)
2241  flux(ks+1,i,j) = vel / mapf(i,j,+2) &
2242  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
2243  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) )
2244 
2245  vel = ( f2h(ke-2,1,i_xvz) &
2246  * 0.25_rp * ( mom(ke-1,i,j)+mom(ke-1,i-1,j)+mom(ke-1,i,j+1)+mom(ke-1,i-1,j+1) ) &
2247  + f2h(ke-2,2,i_xvz) &
2248  * 0.25_rp * ( mom(ke-2,i,j)+mom(ke-2,i-1,j)+mom(ke-2,i,j+1)+mom(ke-2,i-1,j+1) ) ) &
2249  / ( f2h(ke-2,1,i_xvz) &
2250  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i,j+1) ) &
2251  + f2h(ke-2,2,i_xvz) &
2252  * 0.5_rp * ( dens(ke-2,i,j)+dens(ke-2,i,j+1) ) )
2253  vel = vel * j13g(ke-2,i,j)
2254  flux(ke-2,i,j) = vel / mapf(i,j,+2) &
2255  * ( f41 * ( val(ke-1,i,j)+val(ke-2,i,j) ) &
2256  + f42 * ( val(ke,i,j)+val(ke-3,i,j) ) )
2257 
2258  vel = ( f2h(ks+2,1,i_xvz) &
2259  * 0.25_rp * ( mom(ks+3,i,j)+mom(ks+3,i-1,j)+mom(ks+3,i,j+1)+mom(ks+3,i-1,j+1) ) &
2260  + f2h(ks+2,2,i_xvz) &
2261  * 0.25_rp * ( mom(ks+2,i,j)+mom(ks+2,i-1,j)+mom(ks+2,i,j+1)+mom(ks+2,i-1,j+1) ) ) &
2262  / ( f2h(ks+2,1,i_xvz) &
2263  * 0.5_rp * ( dens(ks+3,i,j)+dens(ks+3,i,j+1) ) &
2264  + f2h(ks+2,2,i_xvz) &
2265  * 0.5_rp * ( dens(ks+2,i,j)+dens(ks+2,i,j+1) ) )
2266  vel = vel * j13g(ks+2,i,j)
2267  flux(ks+2,i,j) = vel / mapf(i,j,+2) &
2268  * ( f61 * ( val(ks+3,i,j)+val(ks+2,i,j) ) &
2269  + f62 * ( val(ks+4,i,j)+val(ks+1,i,j) ) &
2270  + f63 * ( val(ks+5,i,j)+val(ks,i,j) ) )
2271 
2272  vel = ( f2h(ke-3,1,i_xvz) &
2273  * 0.25_rp * ( mom(ke-2,i,j)+mom(ke-2,i-1,j)+mom(ke-2,i,j+1)+mom(ke-2,i-1,j+1) ) &
2274  + f2h(ke-3,2,i_xvz) &
2275  * 0.25_rp * ( mom(ke-3,i,j)+mom(ke-3,i-1,j)+mom(ke-3,i,j+1)+mom(ke-3,i-1,j+1) ) ) &
2276  / ( f2h(ke-3,1,i_xvz) &
2277  * 0.5_rp * ( dens(ke-2,i,j)+dens(ke-2,i,j+1) ) &
2278  + f2h(ke-3,2,i_xvz) &
2279  * 0.5_rp * ( dens(ke-3,i,j)+dens(ke-3,i,j+1) ) )
2280  vel = vel * j13g(ke-3,i,j)
2281  flux(ke-3,i,j) = vel / mapf(i,j,+2) &
2282  * ( f61 * ( val(ke-2,i,j)+val(ke-3,i,j) ) &
2283  + f62 * ( val(ke-1,i,j)+val(ke-4,i,j) ) &
2284  + f63 * ( val(ke,i,j)+val(ke-5,i,j) ) )
2285 
2286  flux(ke ,i,j) = 0.0_rp
2287  enddo
2288  enddo
2289  !$acc end kernels
2290  !$omp end do nowait
2291 
2292 
2293 
2294  !$acc end data
2295 
2296  !$omp end parallel
2297  return
2298  end subroutine atmos_dyn_fvm_fluxj13_xvz_cd8
2299 
2300  !-----------------------------------------------------------------------------
2302  subroutine atmos_dyn_fvm_fluxj23_xvz_cd8( &
2303  flux, &
2304  mom, val, DENS, &
2305  GSQRT, J23G, MAPF, &
2306  CDZ, TwoD, &
2307  IIS, IIE, JJS, JJE )
2308  implicit none
2309 
2310  real(rp), intent(inout) :: flux (ka,ia,ja)
2311  real(rp), intent(in) :: mom (ka,ia,ja)
2312  real(rp), intent(in) :: val (ka,ia,ja)
2313  real(rp), intent(in) :: dens (ka,ia,ja)
2314  real(rp), intent(in) :: gsqrt (ka,ia,ja)
2315  real(rp), intent(in) :: j23g (ka,ia,ja)
2316  real(rp), intent(in) :: mapf ( ia,ja,2)
2317  real(rp), intent(in) :: cdz (ka)
2318  logical, intent(in) :: twod
2319  integer, intent(in) :: iis, iie, jjs, jje
2320 
2321  real(rp) :: vel
2322  integer :: k, i, j
2323  !---------------------------------------------------------------------------
2324 
2325  !$omp parallel default(none) private(i,j,k,vel) &
2326  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J23G,MAPF) &
2327  !$omp shared(GSQRT,CDZ,TwoD)
2328 
2329  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, J23G, MAPF, CDZ)
2330 
2331 
2332 
2333  !$omp do OMP_SCHEDULE_ collapse(2)
2334  !$acc kernels
2335  do j = jjs, jje
2336  do i = iis, iie
2337  do k = ks+3, ke-4
2338  vel = ( f2h(k,1,i_xvz) &
2339  * mom(k+1,i,j) &
2340  + f2h(k,2,i_xvz) &
2341  * mom(k,i,j) ) &
2342  / ( f2h(k,1,i_xvz) &
2343  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i,j+1) ) &
2344  + f2h(k,2,i_xvz) &
2345  * 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) ) )
2346  vel = vel * j23g(k,i,j)
2347  flux(k,i,j) = vel / mapf(i,j,+1) &
2348  * ( f81 * ( val(k+1,i,j)+val(k,i,j) ) &
2349  + f82 * ( val(k+2,i,j)+val(k-1,i,j) ) &
2350  + f83 * ( val(k+3,i,j)+val(k-2,i,j) ) &
2351  + f84 * ( val(k+4,i,j)+val(k-3,i,j) ) )
2352  enddo
2353  enddo
2354  enddo
2355  !$acc end kernels
2356  !$omp end do nowait
2357 
2358  !$omp do OMP_SCHEDULE_ collapse(2)
2359  !$acc kernels
2360  do j = jjs, jje
2361  do i = iis, iie
2362  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
2363  ! The flux at KS-1 can be non-zero.
2364  ! To reduce calculations, all the fluxes are set to zero.
2365  flux(ks-1,i,j) = 0.0_rp
2366 
2367  vel = ( f2h(ks,1,i_xvz) &
2368  * mom(ks+1,i,j) &
2369  + f2h(ks,2,i_xvz) &
2370  * mom(ks,i,j) ) &
2371  / ( f2h(ks,1,i_xvz) &
2372  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i,j+1) ) &
2373  + f2h(ks,2,i_xvz) &
2374  * 0.5_rp * ( dens(ks,i,j)+dens(ks,i,j+1) ) )
2375  vel = vel * j23g(ks,i,j)
2376  flux(ks,i,j) = vel / mapf(i,j,+1) &
2377  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) )
2378 
2379  vel = ( f2h(ke-1,1,i_xvz) &
2380  * mom(ke,i,j) &
2381  + f2h(ke-1,2,i_xvz) &
2382  * mom(ke-1,i,j) ) &
2383  / ( f2h(ke-1,1,i_xvz) &
2384  * 0.5_rp * ( dens(ke,i,j)+dens(ke,i,j+1) ) &
2385  + f2h(ke-1,2,i_xvz) &
2386  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i,j+1) ) )
2387  vel = vel * j23g(ke-1,i,j)
2388  flux(ke-1,i,j) = vel / mapf(i,j,+1) &
2389  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) )
2390 
2391  vel = ( f2h(ks+1,1,i_xvz) &
2392  * mom(ks+2,i,j) &
2393  + f2h(ks+1,2,i_xvz) &
2394  * mom(ks+1,i,j) ) &
2395  / ( f2h(ks+1,1,i_xvz) &
2396  * 0.5_rp * ( dens(ks+2,i,j)+dens(ks+2,i,j+1) ) &
2397  + f2h(ks+1,2,i_xvz) &
2398  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i,j+1) ) )
2399  vel = vel * j23g(ks+1,i,j)
2400  flux(ks+1,i,j) = vel / mapf(i,j,+1) &
2401  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
2402  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) )
2403 
2404  vel = ( f2h(ke-2,1,i_xvz) &
2405  * mom(ke-1,i,j) &
2406  + f2h(ke-2,2,i_xvz) &
2407  * mom(ke-2,i,j) ) &
2408  / ( f2h(ke-2,1,i_xvz) &
2409  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i,j+1) ) &
2410  + f2h(ke-2,2,i_xvz) &
2411  * 0.5_rp * ( dens(ke-2,i,j)+dens(ke-2,i,j+1) ) )
2412  vel = vel * j23g(ke-2,i,j)
2413  flux(ke-2,i,j) = vel / mapf(i,j,+1) &
2414  * ( f41 * ( val(ke-1,i,j)+val(ke-2,i,j) ) &
2415  + f42 * ( val(ke,i,j)+val(ke-3,i,j) ) )
2416 
2417  vel = ( f2h(ks+2,1,i_xvz) &
2418  * mom(ks+3,i,j) &
2419  + f2h(ks+2,2,i_xvz) &
2420  * mom(ks+2,i,j) ) &
2421  / ( f2h(ks+2,1,i_xvz) &
2422  * 0.5_rp * ( dens(ks+3,i,j)+dens(ks+3,i,j+1) ) &
2423  + f2h(ks+2,2,i_xvz) &
2424  * 0.5_rp * ( dens(ks+2,i,j)+dens(ks+2,i,j+1) ) )
2425  vel = vel * j23g(ks+2,i,j)
2426  flux(ks+2,i,j) = vel / mapf(i,j,+1) &
2427  * ( f61 * ( val(ks+3,i,j)+val(ks+2,i,j) ) &
2428  + f62 * ( val(ks+4,i,j)+val(ks+1,i,j) ) &
2429  + f63 * ( val(ks+5,i,j)+val(ks,i,j) ) )
2430 
2431  vel = ( f2h(ke-3,1,i_xvz) &
2432  * mom(ke-2,i,j) &
2433  + f2h(ke-3,2,i_xvz) &
2434  * mom(ke-3,i,j) ) &
2435  / ( f2h(ke-3,1,i_xvz) &
2436  * 0.5_rp * ( dens(ke-2,i,j)+dens(ke-2,i,j+1) ) &
2437  + f2h(ke-3,2,i_xvz) &
2438  * 0.5_rp * ( dens(ke-3,i,j)+dens(ke-3,i,j+1) ) )
2439  vel = vel * j23g(ke-3,i,j)
2440  flux(ke-3,i,j) = vel / mapf(i,j,+1) &
2441  * ( f61 * ( val(ke-2,i,j)+val(ke-3,i,j) ) &
2442  + f62 * ( val(ke-1,i,j)+val(ke-4,i,j) ) &
2443  + f63 * ( val(ke,i,j)+val(ke-5,i,j) ) )
2444 
2445  flux(ke ,i,j) = 0.0_rp
2446  enddo
2447  enddo
2448  !$acc end kernels
2449  !$omp end do nowait
2450 
2451 
2452 
2453  !$acc end data
2454 
2455  !$omp end parallel
2456  return
2457  end subroutine atmos_dyn_fvm_fluxj23_xvz_cd8
2458 
2459  !-----------------------------------------------------------------------------
2461  subroutine atmos_dyn_fvm_fluxx_xvz_cd8( &
2462  flux, &
2463  mom, val, DENS, &
2464  GSQRT, MAPF, &
2465  num_diff, &
2466  CDZ, TwoD, &
2467  IIS, IIE, JJS, JJE )
2468  implicit none
2469 
2470  real(rp), intent(inout) :: flux (ka,ia,ja)
2471  real(rp), intent(in) :: mom (ka,ia,ja)
2472  real(rp), intent(in) :: val (ka,ia,ja)
2473  real(rp), intent(in) :: dens (ka,ia,ja)
2474  real(rp), intent(in) :: gsqrt (ka,ia,ja)
2475  real(rp), intent(in) :: mapf ( ia,ja,2)
2476  real(rp), intent(in) :: num_diff(ka,ia,ja)
2477  real(rp), intent(in) :: cdz (ka)
2478  logical, intent(in) :: twod
2479  integer, intent(in) :: iis, iie, jjs, jje
2480 
2481  real(rp) :: vel
2482  integer :: k, i, j
2483  !---------------------------------------------------------------------------
2484 
2485 
2486 
2487  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
2488  !$omp private(vel) &
2489  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff)
2490  !$acc kernels
2491  do j = jjs, jje
2492  do i = iis-1, iie
2493  do k = ks, ke
2494 #ifdef DEBUG
2495  call check( __line__, mom(k,i ,j) )
2496  call check( __line__, mom(k,i,j-1) )
2497 
2498  call check( __line__, val(k,i,j) )
2499  call check( __line__, val(k,i+1,j) )
2500 
2501  call check( __line__, val(k,i-1,j) )
2502  call check( __line__, val(k,i+2,j) )
2503 
2504  call check( __line__, val(k,i-2,j) )
2505  call check( __line__, val(k,i+3,j) )
2506 
2507  call check( __line__, val(k,i-3,j) )
2508  call check( __line__, val(k,i+4,j) )
2509 
2510 #endif
2511  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i,j+1) ) ) &
2512  / ( 0.25_rp * ( dens(k,i,j)+dens(k,i+1,j)+dens(k,i,j+1)+dens(k,i+1,j+1) ) )
2513  flux(k,i,j) = gsqrt(k,i,j) / mapf(i,j,+2) * vel &
2514  * ( f81 * ( val(k,i+1,j)+val(k,i,j) ) &
2515  + f82 * ( val(k,i+2,j)+val(k,i-1,j) ) &
2516  + f83 * ( val(k,i+3,j)+val(k,i-2,j) ) &
2517  + f84 * ( val(k,i+4,j)+val(k,i-3,j) ) ) &
2518  + gsqrt(k,i,j) * num_diff(k,i,j)
2519  enddo
2520  enddo
2521  enddo
2522  !$acc end kernels
2523 #ifdef DEBUG
2524  k = iundef; i = iundef; j = iundef
2525 #endif
2526 
2527 
2528 
2529  return
2530  end subroutine atmos_dyn_fvm_fluxx_xvz_cd8
2531 
2532  !-----------------------------------------------------------------------------
2534  subroutine atmos_dyn_fvm_fluxy_xvz_cd8( &
2535  flux, &
2536  mom, val, DENS, &
2537  GSQRT, MAPF, &
2538  num_diff, &
2539  CDZ, TwoD, &
2540  IIS, IIE, JJS, JJE )
2541  implicit none
2542 
2543  real(rp), intent(inout) :: flux (ka,ia,ja)
2544  real(rp), intent(in) :: mom (ka,ia,ja)
2545  real(rp), intent(in) :: val (ka,ia,ja)
2546  real(rp), intent(in) :: dens (ka,ia,ja)
2547  real(rp), intent(in) :: gsqrt (ka,ia,ja)
2548  real(rp), intent(in) :: mapf ( ia,ja,2)
2549  real(rp), intent(in) :: num_diff(ka,ia,ja)
2550  real(rp), intent(in) :: cdz (ka)
2551  logical, intent(in) :: twod
2552  integer, intent(in) :: iis, iie, jjs, jje
2553 
2554  real(rp) :: vel
2555  integer :: k, i, j
2556  !---------------------------------------------------------------------------
2557 
2558  ! note that y-index is added by -1
2559 
2560 
2561 
2562  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
2563  !$omp private(vel) &
2564  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff)
2565  !$acc kernels
2566  do j = jjs, jje+1
2567  do i = iis, iie
2568  do k = ks, ke
2569 #ifdef DEBUG
2570  call check( __line__, mom(k,i ,j) )
2571  call check( __line__, mom(k,i,j-1) )
2572 
2573  call check( __line__, val(k,i,j-1) )
2574  call check( __line__, val(k,i,j) )
2575 
2576  call check( __line__, val(k,i,j-2) )
2577  call check( __line__, val(k,i,j+1) )
2578 
2579  call check( __line__, val(k,i,j-3) )
2580  call check( __line__, val(k,i,j+2) )
2581 
2582  call check( __line__, val(k,i,j-4) )
2583  call check( __line__, val(k,i,j+3) )
2584 
2585 #endif
2586  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i,j-1) ) ) &
2587  / ( dens(k,i,j) )
2588  flux(k,i,j-1) = gsqrt(k,i,j) / mapf(i,j,+1) * vel &
2589  * ( f81 * ( val(k,i,j)+val(k,i,j-1) ) &
2590  + f82 * ( val(k,i,j+1)+val(k,i,j-2) ) &
2591  + f83 * ( val(k,i,j+2)+val(k,i,j-3) ) &
2592  + f84 * ( val(k,i,j+3)+val(k,i,j-4) ) ) &
2593  + gsqrt(k,i,j) * num_diff(k,i,j)
2594  enddo
2595  enddo
2596  enddo
2597  !$acc end kernels
2598 #ifdef DEBUG
2599  k = iundef; i = iundef; j = iundef
2600 #endif
2601 
2602 
2603 
2604  return
2605  end subroutine atmos_dyn_fvm_fluxy_xvz_cd8
2606 
2607 
2608 
2609 
2610 
2611 
2612 
2614 
2615 !--
2616 ! vi:set readonly sw=4 ts=8
2617 !
2618 !Local Variables:
2619 !mode: f90
2620 !buffer-read-only: t
2621 !End:
2622 !
2623 !++
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_dyn_fvm_flux_cd8::atmos_dyn_fvm_fluxz_uyz_cd8
subroutine, public atmos_dyn_fvm_fluxz_uyz_cd8(flux, mom, val, DENS, GSQRT, J33G, num_diff, CDZ, TwoD, IIS, IIE, JJS, JJE)
calculation z-flux at UY
Definition: scale_atmos_dyn_fvm_flux_cd8.F90:1031
scale_atmos_dyn_fvm_flux_cd8::atmos_dyn_fvm_fluxj13_xvz_cd8
subroutine, public atmos_dyn_fvm_fluxj13_xvz_cd8(flux, mom, val, DENS, GSQRT, J13G, MAPF, CDZ, TwoD, IIS, IIE, JJS, JJE)
calculation J13-flux at XVZ
Definition: scale_atmos_dyn_fvm_flux_cd8.F90:2149
scale_atmos_dyn_fvm_flux_cd8::atmos_dyn_fvm_fluxx_xyz_cd8
subroutine, public atmos_dyn_fvm_fluxx_xyz_cd8(flux, mflx, val, GSQRT, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation X-flux at XYZ
Definition: scale_atmos_dyn_fvm_flux_cd8.F90:346
scale_index
module Index
Definition: scale_index.F90:11
scale_const::const_undef2
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:40
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_cd8::atmos_dyn_fvm_fluxj23_xyw_cd8
subroutine, public atmos_dyn_fvm_fluxj23_xyw_cd8(flux, mom, val, DENS, GSQRT, J23G, MAPF, CDZ, TwoD, IIS, IIE, JJS, JJE)
calculation J23-flux at XYW
Definition: scale_atmos_dyn_fvm_flux_cd8.F90:733
scale_atmos_dyn_fvm_flux_cd8::atmos_dyn_fvm_fluxx_uyz_cd8
subroutine, public atmos_dyn_fvm_fluxx_uyz_cd8(flux, mom, val, DENS, GSQRT, MAPF, num_diff, CDZ, TwoD, IIS, IIE, JJS, JJE)
calculation X-flux at UY
Definition: scale_atmos_dyn_fvm_flux_cd8.F90:1779
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:35
scale_atmos_dyn_fvm_flux_cd8::atmos_dyn_fvm_fluxz_xyz_cd8
subroutine, public atmos_dyn_fvm_fluxz_xyz_cd8(flux, mflx, val, GSQRT, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation z-flux at XYZ
Definition: scale_atmos_dyn_fvm_flux_cd8.F90:201
scale_atmos_dyn_fvm_flux_cd8::atmos_dyn_fvm_fluxz_xvz_cd8
subroutine, public atmos_dyn_fvm_fluxz_xvz_cd8(flux, mom, val, DENS, GSQRT, J33G, num_diff, CDZ, TwoD, IIS, IIE, JJS, JJE)
calculation z-flux at XV
Definition: scale_atmos_dyn_fvm_flux_cd8.F90:1974
scale_atmos_dyn_fvm_flux_cd8::atmos_dyn_fvm_fluxj13_uyz_cd8
subroutine, public atmos_dyn_fvm_fluxj13_uyz_cd8(flux, mom, val, DENS, GSQRT, J13G, MAPF, CDZ, TwoD, IIS, IIE, JJS, JJE)
calculation J13-flux at UYZ
Definition: scale_atmos_dyn_fvm_flux_cd8.F90:1339
scale_atmos_dyn_fvm_flux_cd8::atmos_dyn_fvm_fluxy_xyw_cd8
subroutine, public atmos_dyn_fvm_fluxy_xyw_cd8(flux, mom, val, DENS, GSQRT, MAPF, num_diff, CDZ, TwoD, IIS, IIE, JJS, JJE)
calculation Y-flux at XYW
Definition: scale_atmos_dyn_fvm_flux_cd8.F90:933
scale_atmos_dyn_fvm_flux_cd8::atmos_dyn_fvm_flux_valuew_z_cd8
subroutine, public atmos_dyn_fvm_flux_valuew_z_cd8(valW, mflx, val, GSQRT, CDZ)
value at XYW
Definition: scale_atmos_dyn_fvm_flux_cd8.F90:111
scale_atmos_grid_cartesc_index::i_xyz
integer, public i_xyz
Definition: scale_atmos_grid_cartesC_index.F90:91
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_atmos_dyn_fvm_flux_cd8::atmos_dyn_fvm_fluxx_xyw_cd8
subroutine, public atmos_dyn_fvm_fluxx_xyw_cd8(flux, mom, val, DENS, GSQRT, MAPF, num_diff, CDZ, TwoD, IIS, IIE, JJS, JJE)
calculation X-flux at XYW
Definition: scale_atmos_dyn_fvm_flux_cd8.F90:836
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_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::i_xvz
integer, public i_xvz
Definition: scale_atmos_grid_cartesC_index.F90:96
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_atmos_dyn_fvm_flux_cd8::atmos_dyn_fvm_fluxy_xvz_cd8
subroutine, public atmos_dyn_fvm_fluxy_xvz_cd8(flux, mom, val, DENS, GSQRT, MAPF, num_diff, CDZ, TwoD, IIS, IIE, JJS, JJE)
calculation Y-flux at XV
Definition: scale_atmos_dyn_fvm_flux_cd8.F90:2541
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
scale_atmos_dyn_fvm_flux_cd8::atmos_dyn_fvm_fluxx_xvz_cd8
subroutine, public atmos_dyn_fvm_fluxx_xvz_cd8(flux, mom, val, DENS, GSQRT, MAPF, num_diff, CDZ, TwoD, IIS, IIE, JJS, JJE)
calculation X-flux at XV
Definition: scale_atmos_dyn_fvm_flux_cd8.F90:2468
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
scale_atmos_dyn_fvm_flux_cd8::atmos_dyn_fvm_fluxj23_uyz_cd8
subroutine, public atmos_dyn_fvm_fluxj23_uyz_cd8(flux, mom, val, DENS, GSQRT, J23G, MAPF, CDZ, TwoD, IIS, IIE, JJS, JJE)
calculation J23-flux at UYZ
Definition: scale_atmos_dyn_fvm_flux_cd8.F90:1498
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_atmos_dyn_fvm_flux_cd8::atmos_dyn_fvm_fluxz_xyw_cd8
subroutine, public atmos_dyn_fvm_fluxz_xyw_cd8(flux, mom, val, DENS, GSQRT, J33G, num_diff, CDZ, FDZ, dtrk, IIS, IIE, JJS, JJE)
calculation z-flux at XYW
Definition: scale_atmos_dyn_fvm_flux_cd8.F90:475
scale_debug
module DEBUG
Definition: scale_debug.F90:11
scale_atmos_dyn_fvm_flux_cd8::atmos_dyn_fvm_fluxy_uyz_cd8
subroutine, public atmos_dyn_fvm_fluxy_uyz_cd8(flux, mom, val, DENS, GSQRT, MAPF, num_diff, CDZ, TwoD, IIS, IIE, JJS, JJE)
calculation Y-flux at UY
Definition: scale_atmos_dyn_fvm_flux_cd8.F90:1854
scale_atmos_dyn_fvm_flux_cd8::atmos_dyn_fvm_fluxy_xyz_cd8
subroutine, public atmos_dyn_fvm_fluxy_xyz_cd8(flux, mflx, val, GSQRT, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation Y-flux at XYZ
Definition: scale_atmos_dyn_fvm_flux_cd8.F90:409
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_atmos_dyn_fvm_flux_cd8::atmos_dyn_fvm_fluxj13_xyw_cd8
subroutine, public atmos_dyn_fvm_fluxj13_xyw_cd8(flux, mom, val, DENS, GSQRT, J13G, MAPF, CDZ, TwoD, IIS, IIE, JJS, JJE)
calculation J13-flux at XYW
Definition: scale_atmos_dyn_fvm_flux_cd8.F90:632
scale_atmos_dyn_fvm_flux_cd8
module scale_atmos_dyn_fvm_flux_cd8
Definition: scale_atmos_dyn_fvm_flux_cd8.F90:16
scale_atmos_dyn_fvm_flux_cd8::atmos_dyn_fvm_fluxj23_xvz_cd8
subroutine, public atmos_dyn_fvm_fluxj23_xvz_cd8(flux, mom, val, DENS, GSQRT, J23G, MAPF, CDZ, TwoD, IIS, IIE, JJS, JJE)
calculation J23-flux at XVZ
Definition: scale_atmos_dyn_fvm_flux_cd8.F90:2308