SCALE-RM
scale_atmos_dyn_fvm_flux_cd6.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 
97 
98 contains
99  !-----------------------------------------------------------------------------
101 !OCL SERIAL
102  subroutine atmos_dyn_fvm_flux_valuew_z_cd6( &
103  valW, &
104  mflx, val, GSQRT, &
105  CDZ )
106  !$acc routine vector
107  implicit none
108 
109  real(rp), intent(out) :: valw (ka)
110  real(rp), intent(in) :: mflx (ka)
111  real(rp), intent(in) :: val (ka)
112  real(rp), intent(in) :: gsqrt(ka)
113  real(rp), intent(in) :: cdz (ka)
114 
115  integer :: k
116  !---------------------------------------------------------------------------
117 
118  do k = ks+2, ke-3
119 #ifdef DEBUG
120  call check( __line__, mflx(k) )
121 
122  call check( __line__, val(k) )
123  call check( __line__, val(k+1) )
124 
125  call check( __line__, val(k-1) )
126  call check( __line__, val(k+2) )
127 
128  call check( __line__, val(k-2) )
129  call check( __line__, val(k+3) )
130 
131 #endif
132  valw(k) = f61 * ( val(k+1)+val(k) ) &
133  + f62 * ( val(k+2)+val(k-1) ) &
134  + f63 * ( val(k+3)+val(k-2) )
135  enddo
136 #ifdef DEBUG
137  k = iundef
138 #endif
139 
140 #ifdef DEBUG
141 
142  call check( __line__, mflx(ks) )
143  call check( __line__, val(ks ) )
144  call check( __line__, val(ks+1) )
145  call check( __line__, mflx(ke-1) )
146  call check( __line__, val(ke ) )
147  call check( __line__, val(ke-1) )
148 
149  call check( __line__, mflx(ks+1) )
150  call check( __line__, val(ks+2 ) )
151  call check( __line__, val(ks+3) )
152  call check( __line__, mflx(ke-2) )
153  call check( __line__, val(ke-2 ) )
154  call check( __line__, val(ke-3) )
155 
156 #endif
157 
158  valw(ks) = f2 * ( val(ks+1)+val(ks) )
159  valw(ke-1) = f2 * ( val(ke)+val(ke-1) )
160 
161  valw(ks+1) = f41 * ( val(ks+2)+val(ks+1) ) &
162  + f42 * ( val(ks+3)+val(ks) )
163  valw(ke-2) = f41 * ( val(ke-1)+val(ke-2) ) &
164  + f42 * ( val(ke)+val(ke-3) )
165 
166 
167  return
168  end subroutine atmos_dyn_fvm_flux_valuew_z_cd6
169 
170  !-----------------------------------------------------------------------------
172  subroutine atmos_dyn_fvm_fluxz_xyz_cd6( &
173  flux, &
174  mflx, val, GSQRT, &
175  num_diff, &
176  CDZ, &
177  IIS, IIE, JJS, JJE )
178  use scale_const, only: &
179  eps => const_eps
180  implicit none
181 
182  real(rp), intent(inout) :: flux (ka,ia,ja)
183  real(rp), intent(in) :: mflx (ka,ia,ja)
184  real(rp), intent(in) :: val (ka,ia,ja)
185  real(rp), intent(in) :: gsqrt (ka,ia,ja)
186  real(rp), intent(in) :: num_diff(ka,ia,ja)
187  real(rp), intent(in) :: cdz (ka)
188  integer, intent(in) :: iis, iie, jjs, jje
189 
190  real(rp) :: vel
191  integer :: k, i, j
192  !---------------------------------------------------------------------------
193 
194  !$omp parallel default(none) private(i,j,k, vel) &
195  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mflx,val,flux,GSQRT,num_diff,EPS)
196 
197  !$acc data copy(flux) copyin(mflx, val, GSQRT, num_diff, CDZ)
198 
199  !$omp do OMP_SCHEDULE_ collapse(2)
200  !$acc kernels
201  do j = jjs, jje
202  do i = iis, iie
203  do k = ks+2, ke-3
204 #ifdef DEBUG
205  call check( __line__, mflx(k,i,j) )
206 
207  call check( __line__, val(k,i,j) )
208  call check( __line__, val(k+1,i,j) )
209 
210  call check( __line__, val(k-1,i,j) )
211  call check( __line__, val(k+2,i,j) )
212 
213  call check( __line__, val(k-2,i,j) )
214  call check( __line__, val(k+3,i,j) )
215 
216 #endif
217  vel = mflx(k,i,j)
218  flux(k,i,j) = vel &
219  * ( f61 * ( val(k+1,i,j)+val(k,i,j) ) &
220  + f62 * ( val(k+2,i,j)+val(k-1,i,j) ) &
221  + f63 * ( val(k+3,i,j)+val(k-2,i,j) ) ) &
222  + gsqrt(k,i,j) * num_diff(k,i,j)
223  enddo
224  enddo
225  enddo
226  !$acc end kernels
227  !$omp end do nowait
228 #ifdef DEBUG
229  k = iundef; i = iundef; j = iundef
230 #endif
231 
232  !$omp do OMP_SCHEDULE_ collapse(2)
233  !$acc kernels
234  do j = jjs, jje
235  do i = iis, iie
236 #ifdef DEBUG
237 
238  call check( __line__, mflx(ks,i,j) )
239  call check( __line__, val(ks ,i,j) )
240  call check( __line__, val(ks+1,i,j) )
241  call check( __line__, mflx(ke-1,i,j) )
242  call check( __line__, val(ke ,i,j) )
243  call check( __line__, val(ke-1,i,j) )
244 
245  call check( __line__, mflx(ks+1,i,j) )
246  call check( __line__, val(ks+2 ,i,j) )
247  call check( __line__, val(ks+3,i,j) )
248  call check( __line__, mflx(ke-2,i,j) )
249  call check( __line__, val(ke-2 ,i,j) )
250  call check( __line__, val(ke-3,i,j) )
251 
252 #endif
253  flux(ks-1,i,j) = 0.0_rp
254 
255  vel = mflx(ks,i,j)
256  flux(ks,i,j) = vel &
257  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) ) &
258  + gsqrt(ks,i,j) * num_diff(ks,i,j)
259  vel = mflx(ke-1,i,j)
260  flux(ke-1,i,j) = vel &
261  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) ) &
262  + gsqrt(ke-1,i,j) * num_diff(ke-1,i,j)
263 
264  vel = mflx(ks+1,i,j)
265  flux(ks+1,i,j) = vel &
266  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
267  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) ) &
268  + gsqrt(ks+1,i,j) * num_diff(ks+1,i,j)
269  vel = mflx(ke-2,i,j)
270  flux(ke-2,i,j) = vel &
271  * ( f41 * ( val(ke-1,i,j)+val(ke-2,i,j) ) &
272  + f42 * ( val(ke,i,j)+val(ke-3,i,j) ) ) &
273  + gsqrt(ke-2,i,j) * num_diff(ke-2,i,j)
274 
275  flux(ke ,i,j) = 0.0_rp
276  enddo
277  enddo
278  !$acc end kernels
279  !$omp end do nowait
280 
281  !$acc end data
282 
283  !$omp end parallel
284 #ifdef DEBUG
285  k = iundef; i = iundef; j = iundef
286 #endif
287 
288  return
289  end subroutine atmos_dyn_fvm_fluxz_xyz_cd6
290 
291  !-----------------------------------------------------------------------------
293  subroutine atmos_dyn_fvm_fluxx_xyz_cd6( &
294  flux, &
295  mflx, val, GSQRT, &
296  num_diff, &
297  CDZ, &
298  IIS, IIE, JJS, JJE )
299  implicit none
300 
301  real(rp), intent(inout) :: flux (ka,ia,ja)
302  real(rp), intent(in) :: mflx (ka,ia,ja)
303  real(rp), intent(in) :: val (ka,ia,ja)
304  real(rp), intent(in) :: gsqrt (ka,ia,ja)
305  real(rp), intent(in) :: num_diff(ka,ia,ja)
306  real(rp), intent(in) :: cdz(ka)
307  integer, intent(in) :: iis, iie, jjs, jje
308 
309  real(rp) :: vel
310  integer :: k, i, j
311  !---------------------------------------------------------------------------
312 
313  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
314  !$omp private(vel) &
315  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mflx,val,flux,GSQRT,num_diff)
316  !$acc kernels
317  do j = jjs, jje
318  do i = iis-1, iie
319  do k = ks, ke
320 #ifdef DEBUG
321  call check( __line__, mflx(k,i,j) )
322 
323  call check( __line__, val(k,i,j) )
324  call check( __line__, val(k,i+1,j) )
325 
326  call check( __line__, val(k,i-1,j) )
327  call check( __line__, val(k,i+2,j) )
328 
329  call check( __line__, val(k,i-2,j) )
330  call check( __line__, val(k,i+3,j) )
331 
332 #endif
333  vel = mflx(k,i,j)
334  flux(k,i,j) = vel &
335  * ( f61 * ( val(k,i+1,j)+val(k,i,j) ) &
336  + f62 * ( val(k,i+2,j)+val(k,i-1,j) ) &
337  + f63 * ( val(k,i+3,j)+val(k,i-2,j) ) ) &
338  + gsqrt(k,i,j) * num_diff(k,i,j)
339  enddo
340  enddo
341  enddo
342  !$acc end kernels
343 #ifdef DEBUG
344  k = iundef; i = iundef; j = iundef
345 #endif
346 
347  return
348  end subroutine atmos_dyn_fvm_fluxx_xyz_cd6
349 
350  !-----------------------------------------------------------------------------
352  subroutine atmos_dyn_fvm_fluxy_xyz_cd6( &
353  flux, &
354  mflx, val, GSQRT, &
355  num_diff, &
356  CDZ, &
357  IIS, IIE, JJS, JJE )
358  implicit none
359 
360  real(rp), intent(inout) :: flux (ka,ia,ja)
361  real(rp), intent(in) :: mflx (ka,ia,ja)
362  real(rp), intent(in) :: val (ka,ia,ja)
363  real(rp), intent(in) :: gsqrt (ka,ia,ja)
364  real(rp), intent(in) :: num_diff(ka,ia,ja)
365  real(rp), intent(in) :: cdz(ka)
366  integer, intent(in) :: iis, iie, jjs, jje
367 
368  real(rp) :: vel
369  integer :: k, i, j
370  !---------------------------------------------------------------------------
371 
372  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
373  !$omp private(vel) &
374  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mflx,val,flux,GSQRT,num_diff)
375  !$acc kernels
376  do j = jjs-1, jje
377  do i = iis, iie
378  do k = ks, ke
379 #ifdef DEBUG
380  call check( __line__, mflx(k,i,j) )
381 
382  call check( __line__, val(k,i,j) )
383  call check( __line__, val(k,i,j+1) )
384 
385  call check( __line__, val(k,i,j-1) )
386  call check( __line__, val(k,i,j+2) )
387 
388  call check( __line__, val(k,i,j-2) )
389  call check( __line__, val(k,i,j+3) )
390 
391 #endif
392  vel = mflx(k,i,j)
393  flux(k,i,j) = vel &
394  * ( f61 * ( val(k,i,j+1)+val(k,i,j) ) &
395  + f62 * ( val(k,i,j+2)+val(k,i,j-1) ) &
396  + f63 * ( val(k,i,j+3)+val(k,i,j-2) ) ) &
397  + gsqrt(k,i,j) * num_diff(k,i,j)
398  enddo
399  enddo
400  enddo
401  !$acc end kernels
402 #ifdef DEBUG
403  k = iundef; i = iundef; j = iundef
404 #endif
405 
406  return
407  end subroutine atmos_dyn_fvm_fluxy_xyz_cd6
408 
409 
410  !-----------------------------------------------------------------------------
412  subroutine atmos_dyn_fvm_fluxz_xyw_cd6( &
413  flux, &
414  mom, val, DENS, &
415  GSQRT, J33G, &
416  num_diff, &
417  CDZ, FDZ, &
418  dtrk, &
419  IIS, IIE, JJS, JJE )
420  implicit none
421 
422  real(rp), intent(inout) :: flux (ka,ia,ja)
423  real(rp), intent(in) :: mom (ka,ia,ja)
424  real(rp), intent(in) :: val (ka,ia,ja)
425  real(rp), intent(in) :: dens (ka,ia,ja)
426  real(rp), intent(in) :: gsqrt (ka,ia,ja)
427  real(rp), intent(in) :: j33g
428  real(rp), intent(in) :: num_diff(ka,ia,ja)
429  real(rp), intent(in) :: cdz (ka)
430  real(rp), intent(in) :: fdz (ka-1)
431  real(rp), intent(in) :: dtrk
432  integer, intent(in) :: iis, iie, jjs, jje
433 
434  real(rp) :: vel
435  integer :: k, i, j
436  !---------------------------------------------------------------------------
437 
438  ! note than z-index is added by -1
439 
440  !$omp parallel default(none) private(i,j,k,vel) &
441  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,flux,J33G,GSQRT,num_diff,DENS,FDZ,dtrk)
442 
443  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, num_diff, CDZ, FDZ)
444 
445  !$omp do OMP_SCHEDULE_ collapse(2)
446  !$acc kernels
447  do j = jjs, jje
448  do i = iis, iie
449  do k = ks+3, ke-2
450 #ifdef DEBUG
451  call check( __line__, mom(k-1,i,j) )
452  call check( __line__, mom(k ,i,j) )
453 
454  call check( __line__, val(k-1,i,j) )
455  call check( __line__, val(k,i,j) )
456 
457  call check( __line__, val(k-2,i,j) )
458  call check( __line__, val(k+1,i,j) )
459 
460  call check( __line__, val(k-3,i,j) )
461  call check( __line__, val(k+2,i,j) )
462 
463 #endif
464  vel = ( 0.5_rp * ( mom(k-1,i,j) &
465  + mom(k,i,j) ) ) &
466  / dens(k,i,j)
467  flux(k-1,i,j) = j33g * vel &
468  * ( f61 * ( val(k,i,j)+val(k-1,i,j) ) &
469  + f62 * ( val(k+1,i,j)+val(k-2,i,j) ) &
470  + f63 * ( val(k+2,i,j)+val(k-3,i,j) ) ) &
471  + gsqrt(k,i,j) * num_diff(k,i,j)
472  enddo
473  enddo
474  enddo
475  !$acc end kernels
476  !$omp end do nowait
477 #ifdef DEBUG
478  k = iundef; i = iundef; j = iundef
479 #endif
480 
481  !$omp do OMP_SCHEDULE_ collapse(2)
482  !$acc kernels
483  do j = jjs, jje
484  do i = iis, iie
485 #ifdef DEBUG
486 
487  call check( __line__, val(ks,i,j) )
488  call check( __line__, val(ks+1,i,j) )
489  call check( __line__, val(ks+2,i,j) )
490  call check( __line__, val(ks+3,i,j) )
491 
492 
493  call check( __line__, val(ke-3,i,j) )
494  call check( __line__, val(ke-2,i,j) )
495  call check( __line__, val(ke-1,i,j) )
496 
497 #endif
498  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS.
499  ! The flux at KS can be non-zero.
500  ! To reduce calculations, all the fluxes are set to zero.
501  flux(ks-1,i,j) = 0.0_rp ! k = KS
502 
503  vel = ( 0.5_rp * ( mom(ks,i,j) &
504  + mom(ks+1,i,j) ) ) &
505  / dens(ks+1,i,j)
506  flux(ks,i,j) = j33g * vel &
507  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) ) &
508  + gsqrt(ks+1,i,j) * num_diff(ks+1,i,j) ! k = KS+1
509 
510  vel = ( 0.5_rp * ( mom(ks+1,i,j) &
511  + mom(ks+2,i,j) ) ) &
512  / dens(ks+2,i,j)
513  flux(ks+1,i,j) = j33g * vel &
514  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
515  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) ) &
516  + gsqrt(ks+2,i,j) * num_diff(ks+2,i,j) ! k = KS+2
517 
518 
519 
520  vel = ( 0.5_rp * ( mom(ke-2,i,j) &
521  + mom(ke-1,i,j) ) ) &
522  / dens(ke-1,i,j)
523  flux(ke-2,i,j) = j33g * vel &
524  * ( f41 * ( val(ke-1,i,j)+val(ke-2,i,j) ) &
525  + f42 * ( val(ke,i,j)+val(ke-3,i,j) ) ) &
526  + gsqrt(ke-1,i,j) * num_diff(ke-1,i,j) ! k = KE-1
527 
528  flux(ke-1,i,j) = 0.0_rp ! k = KE
529  flux(ke ,i,j) = 0.0_rp ! k = KE+1
530  enddo
531  enddo
532  !$acc end kernels
533  !$omp end do nowait
534 
535  !$acc end data
536 
537  !$omp end parallel
538 
539  return
540  end subroutine atmos_dyn_fvm_fluxz_xyw_cd6
541 
542 
543  !-----------------------------------------------------------------------------
545  subroutine atmos_dyn_fvm_fluxj13_xyw_cd6( &
546  flux, &
547  mom, val, DENS, &
548  GSQRT, J13G, MAPF, &
549  CDZ, TwoD, &
550  IIS, IIE, JJS, JJE )
551  implicit none
552 
553  real(rp), intent(inout) :: flux (ka,ia,ja)
554  real(rp), intent(in) :: mom (ka,ia,ja)
555  real(rp), intent(in) :: val (ka,ia,ja)
556  real(rp), intent(in) :: dens (ka,ia,ja)
557  real(rp), intent(in) :: gsqrt (ka,ia,ja)
558  real(rp), intent(in) :: j13g (ka,ia,ja)
559  real(rp), intent(in) :: mapf ( ia,ja,2)
560  real(rp), intent(in) :: cdz (ka)
561  logical, intent(in) :: twod
562  integer, intent(in) :: iis, iie, jjs, jje
563 
564  real(rp) :: vel
565  integer :: k, i, j
566  !---------------------------------------------------------------------------
567 
568  !$omp parallel default(none) private(i,j,k,vel) &
569  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J13G,MAPF)
570 
571  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, J13G, MAPF, CDZ)
572 
573  !$omp do OMP_SCHEDULE_ collapse(2)
574  !$acc kernels
575  do j = jjs, jje
576  do i = iis, iie
577  do k = ks+2, ke-2
578  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i-1,j) ) ) &
579  / dens(k,i,j)
580  vel = vel * j13g(k,i,j)
581  flux(k-1,i,j) = vel / mapf(i,j,+2) &
582  * ( f61 * ( val(k,i,j)+val(k-1,i,j) ) &
583  + f62 * ( val(k+1,i,j)+val(k-2,i,j) ) &
584  + f63 * ( val(k+2,i,j)+val(k-3,i,j) ) )
585  enddo
586  enddo
587  enddo
588  !$acc end kernels
589  !$omp end do nowait
590 
591  !$omp do OMP_SCHEDULE_ collapse(2)
592  !$acc kernels
593  do j = jjs, jje
594  do i = iis, iie
595  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS.
596  ! The flux at KS can be non-zero.
597  ! To reduce calculations, all the fluxes are set to zero.
598  flux(ks-1,i,j) = 0.0_rp ! k = KS
599 
600  ! physically incorrect but for numerical stability
601  vel = ( ( 0.5_rp * ( mom(ks+1,i,j)+mom(ks+1,i-1,j) ) ) / dens(ks+1,i,j) &
602  + ( 0.5_rp * ( mom(ks,i,j)+mom(ks,i-1,j) ) ) / dens(ks ,i,j) ) * 0.5_rp
603 ! vel = ( 0.5_RP * ( mom(KS+1,i,j)+mom(KS+1,i-1,j) ) ) &
604 ! / DENS(KS+1,i,j)
605  vel = vel * j13g(ks+1,i,j)
606  flux(ks,i,j) = vel / mapf(i,j,+2) &
607  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) ) ! k = KS+1
608 
609 
610  vel = ( 0.5_rp * ( mom(ke-1,i,j)+mom(ke-1,i-1,j) ) ) &
611  / dens(ke-1,i,j)
612  vel = vel * j13g(ke-1,i,j)
613  flux(ke-2,i,j) = vel / mapf(i,j,+2) &
614  * ( f2 * ( val(ke-1,i,j)+val(ke-2,i,j) ) ) ! k = KE-3
615 
616  flux(ke-1,i,j) = 0.0_rp
617  enddo
618  enddo
619  !$acc end kernels
620  !$omp end do nowait
621 
622  !$acc end data
623 
624  !$omp end parallel
625 
626  return
627  end subroutine atmos_dyn_fvm_fluxj13_xyw_cd6
628 
629  !-----------------------------------------------------------------------------
631  subroutine atmos_dyn_fvm_fluxj23_xyw_cd6( &
632  flux, &
633  mom, val, DENS, &
634  GSQRT, J23G, MAPF, &
635  CDZ, TwoD, &
636  IIS, IIE, JJS, JJE )
637  implicit none
638 
639  real(rp), intent(inout) :: flux (ka,ia,ja)
640  real(rp), intent(in) :: mom (ka,ia,ja)
641  real(rp), intent(in) :: val (ka,ia,ja)
642  real(rp), intent(in) :: dens (ka,ia,ja)
643  real(rp), intent(in) :: gsqrt (ka,ia,ja)
644  real(rp), intent(in) :: j23g (ka,ia,ja)
645  real(rp), intent(in) :: mapf ( ia,ja,2)
646  real(rp), intent(in) :: cdz (ka)
647  logical, intent(in) :: twod
648  integer, intent(in) :: iis, iie, jjs, jje
649 
650  real(rp) :: vel
651  integer :: k, i, j
652  !---------------------------------------------------------------------------
653 
654  !$omp parallel default(none) private(i,j,k,vel) &
655  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J23G,MAPF)
656 
657  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, J23G, MAPF, CDZ)
658 
659  !$omp do OMP_SCHEDULE_ collapse(2)
660  !$acc kernels
661  do j = jjs, jje
662  do i = iis, iie
663  do k = ks+2, ke-2
664  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i,j-1) ) ) &
665  / dens(k,i,j)
666  vel = vel * j23g(k,i,j)
667  flux(k-1,i,j) = vel / mapf(i,j,+1) &
668  * ( f61 * ( val(k,i,j)+val(k-1,i,j) ) &
669  + f62 * ( val(k+1,i,j)+val(k-2,i,j) ) &
670  + f63 * ( val(k+2,i,j)+val(k-3,i,j) ) )
671  enddo
672  enddo
673  enddo
674  !$acc end kernels
675  !$omp end do nowait
676 
677  !$omp do OMP_SCHEDULE_ collapse(2)
678  !$acc kernels
679  do j = jjs, jje
680  do i = iis, iie
681  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS.
682  ! The flux at KS can be non-zero.
683  ! To reduce calculations, all the fluxes are set to zero.
684  flux(ks-1,i,j) = 0.0_rp ! k = KS
685 
686  ! physically incorrect but for numerical stability
687  vel = ( ( 0.5_rp * ( mom(ks+1,i,j)+mom(ks+1,i,j-1) ) ) / dens(ks+1,i,j) &
688  + ( 0.5_rp * ( mom(ks,i,j)+mom(ks,i,j-1) ) ) / dens(ks ,i,j) ) * 0.5_rp
689 ! vel = ( 0.5_RP * ( mom(KS+1,i,j)+mom(KS+1,i,j-1) ) ) &
690 ! / DENS(KS+1,i,j)
691  vel = vel * j23g(ks+1,i,j)
692  flux(ks,i,j) = vel / mapf(i,j,+1) &
693  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) ) ! k = KS+1
694 
695 
696  vel = ( 0.5_rp * ( mom(ke-1,i,j)+mom(ke-1,i,j-1) ) ) &
697  / dens(ke-1,i,j)
698  vel = vel * j23g(ke-1,i,j)
699  flux(ke-2,i,j) = vel / mapf(i,j,+1) &
700  * ( f2 * ( val(ke-1,i,j)+val(ke-2,i,j) ) ) ! k = KE-3
701 
702  flux(ke-1,i,j) = 0.0_rp
703  enddo
704  enddo
705  !$acc end kernels
706  !$omp end do nowait
707 
708  !$acc end data
709 
710  !$omp end parallel
711 
712  return
713  end subroutine atmos_dyn_fvm_fluxj23_xyw_cd6
714 
715 
716  !-----------------------------------------------------------------------------
718  subroutine atmos_dyn_fvm_fluxx_xyw_cd6( &
719  flux, &
720  mom, val, DENS, &
721  GSQRT, MAPF, &
722  num_diff, &
723  CDZ, TwoD, &
724  IIS, IIE, JJS, JJE )
725  implicit none
726 
727  real(rp), intent(inout) :: flux (ka,ia,ja)
728  real(rp), intent(in) :: mom (ka,ia,ja)
729  real(rp), intent(in) :: val (ka,ia,ja)
730  real(rp), intent(in) :: dens (ka,ia,ja)
731  real(rp), intent(in) :: gsqrt (ka,ia,ja)
732  real(rp), intent(in) :: mapf ( ia,ja,2)
733  real(rp), intent(in) :: num_diff(ka,ia,ja)
734  real(rp), intent(in) :: cdz (ka)
735  logical, intent(in) :: twod
736  integer, intent(in) :: iis, iie, jjs, jje
737 
738  real(rp) :: vel
739  integer :: k, i, j
740  !---------------------------------------------------------------------------
741 
742  !$omp parallel default(none) private(i,j,k,vel) &
743  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff) &
744  !$omp shared(CDZ)
745 
746  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, MAPF, num_diff, CDZ)
747 
748  !$omp do OMP_SCHEDULE_ collapse(2)
749  !$acc kernels
750  do j = jjs, jje
751  do i = iis-1, iie
752  do k = ks, ke-1
753 #ifdef DEBUG
754  call check( __line__, mom(k ,i,j) )
755  call check( __line__, mom(k+1,i,j) )
756 
757  call check( __line__, val(k,i,j) )
758  call check( __line__, val(k,i+1,j) )
759 
760  call check( __line__, val(k,i-1,j) )
761  call check( __line__, val(k,i+2,j) )
762 
763  call check( __line__, val(k,i-2,j) )
764  call check( __line__, val(k,i+3,j) )
765 
766 #endif
767  vel = ( f2h(k,1,i_uyz) &
768  * mom(k+1,i,j) &
769  + f2h(k,2,i_uyz) &
770  * mom(k,i,j) ) &
771  / ( f2h(k,1,i_uyz) &
772  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i+1,j) ) &
773  + f2h(k,2,i_uyz) &
774  * 0.5_rp * ( dens(k,i,j)+dens(k,i+1,j) ) )
775  flux(k,i,j) = gsqrt(k,i,j) / mapf(i,j,+2) * vel &
776  * ( f61 * ( val(k,i+1,j)+val(k,i,j) ) &
777  + f62 * ( val(k,i+2,j)+val(k,i-1,j) ) &
778  + f63 * ( val(k,i+3,j)+val(k,i-2,j) ) ) &
779  + gsqrt(k,i,j) * num_diff(k,i,j)
780  enddo
781  enddo
782  enddo
783  !$acc end kernels
784  !$omp end do nowait
785 #ifdef DEBUG
786  k = iundef; i = iundef; j = iundef
787 #endif
788 
789  !$omp do OMP_SCHEDULE_ collapse(2)
790  !$acc kernels
791  do j = jjs, jje
792  do i = iis-1, iie
793  flux(ke,i,j) = 0.0_rp
794  enddo
795  enddo
796  !$acc end kernels
797  !$omp end do nowait
798 
799  !$acc end data
800 
801  !$omp end parallel
802 #ifdef DEBUG
803  k = iundef; i = iundef; j = iundef
804 #endif
805 
806  return
807  end subroutine atmos_dyn_fvm_fluxx_xyw_cd6
808 
809  !-----------------------------------------------------------------------------
811  subroutine atmos_dyn_fvm_fluxy_xyw_cd6( &
812  flux, &
813  mom, val, DENS, &
814  GSQRT, MAPF, &
815  num_diff, &
816  CDZ, TwoD, &
817  IIS, IIE, JJS, JJE )
818  implicit none
819 
820  real(rp), intent(inout) :: flux (ka,ia,ja)
821  real(rp), intent(in) :: mom (ka,ia,ja)
822  real(rp), intent(in) :: val (ka,ia,ja)
823  real(rp), intent(in) :: dens (ka,ia,ja)
824  real(rp), intent(in) :: gsqrt (ka,ia,ja)
825  real(rp), intent(in) :: mapf ( ia,ja,2)
826  real(rp), intent(in) :: num_diff(ka,ia,ja)
827  real(rp), intent(in) :: cdz (ka)
828  logical, intent(in) :: twod
829  integer, intent(in) :: iis, iie, jjs, jje
830 
831  real(rp) :: vel
832  integer :: k, i, j
833  !---------------------------------------------------------------------------
834 
835  !$omp parallel default(none) private(i,j,k,vel) &
836  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff) &
837  !$omp shared(CDZ)
838 
839  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, MAPF, num_diff, CDZ)
840 
841  !$omp do OMP_SCHEDULE_ collapse(2)
842  !$acc kernels
843  do j = jjs-1, jje
844  do i = iis, iie
845  do k = ks, ke-1
846 #ifdef DEBUG
847  call check( __line__, mom(k ,i,j) )
848  call check( __line__, mom(k+1,i,j) )
849 
850  call check( __line__, val(k,i,j) )
851  call check( __line__, val(k,i,j+1) )
852 
853  call check( __line__, val(k,i,j-1) )
854  call check( __line__, val(k,i,j+2) )
855 
856  call check( __line__, val(k,i,j-2) )
857  call check( __line__, val(k,i,j+3) )
858 
859 #endif
860  vel = ( f2h(k,1,i_xvz) &
861  * mom(k+1,i,j) &
862  + f2h(k,2,i_xvz) &
863  * mom(k,i,j) ) &
864  / ( f2h(k,1,i_xvz) &
865  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i,j+1) ) &
866  + f2h(k,2,i_xvz) &
867  * 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) ) )
868  flux(k,i,j) = gsqrt(k,i,j) / mapf(i,j,+1) * vel &
869  * ( f61 * ( val(k,i,j+1)+val(k,i,j) ) &
870  + f62 * ( val(k,i,j+2)+val(k,i,j-1) ) &
871  + f63 * ( val(k,i,j+3)+val(k,i,j-2) ) ) &
872  + gsqrt(k,i,j) * num_diff(k,i,j)
873  enddo
874  enddo
875  enddo
876  !$acc end kernels
877  !$omp end do nowait
878 #ifdef DEBUG
879  k = iundef; i = iundef; j = iundef
880 #endif
881 
882  !$omp do OMP_SCHEDULE_ collapse(2)
883  !$acc kernels
884  do j = jjs-1, jje
885  do i = iis, iie
886  flux(ke,i,j) = 0.0_rp
887  enddo
888  enddo
889  !$acc end kernels
890  !$omp end do nowait
891 
892  !$acc end data
893 
894  !$omp end parallel
895 #ifdef DEBUG
896  k = iundef; i = iundef; j = iundef
897 #endif
898 
899  return
900  end subroutine atmos_dyn_fvm_fluxy_xyw_cd6
901 
902 
903  !-----------------------------------------------------------------------------
905  subroutine atmos_dyn_fvm_fluxz_uyz_cd6( &
906  flux, &
907  mom, val, DENS, &
908  GSQRT, J33G, &
909  num_diff, &
910  CDZ, TwoD, &
911  IIS, IIE, JJS, JJE )
912  implicit none
913 
914  real(rp), intent(inout) :: flux (ka,ia,ja)
915  real(rp), intent(in) :: mom (ka,ia,ja)
916  real(rp), intent(in) :: val (ka,ia,ja)
917  real(rp), intent(in) :: dens (ka,ia,ja)
918  real(rp), intent(in) :: gsqrt (ka,ia,ja)
919  real(rp), intent(in) :: j33g
920  real(rp), intent(in) :: num_diff(ka,ia,ja)
921  real(rp), intent(in) :: cdz (ka)
922  logical, intent(in) :: twod
923  integer, intent(in) :: iis, iie, jjs, jje
924 
925  real(rp) :: vel
926  integer :: k, i, j
927  !---------------------------------------------------------------------------
928 
929  !$omp parallel default(none) private(i,j,k,vel) &
930  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J33G,GSQRT,num_diff) &
931  !$omp shared(CDZ,TwoD)
932 
933  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, num_diff, CDZ)
934 
935 
936  if ( twod ) then
937 
938  !$omp do OMP_SCHEDULE_
939  !$acc kernels
940  do j = jjs, jje
941  do k = ks+2, ke-3
942  i = iis
943 #ifdef DEBUG
944  call check( __line__, mom(k,i,j) )
945 
946  call check( __line__, val(k,i,j) )
947  call check( __line__, val(k+1,i,j) )
948 
949  call check( __line__, val(k-1,i,j) )
950  call check( __line__, val(k+2,i,j) )
951 
952  call check( __line__, val(k-2,i,j) )
953  call check( __line__, val(k+3,i,j) )
954 
955 #endif
956  vel = ( mom(k,i,j) ) &
957  / ( f2h(k,1,i_xyz) &
958  * dens(k+1,i,j) &
959  + f2h(k,2,i_xyz) &
960  * dens(k,i,j) )
961  flux(k,i,j) = j33g * vel &
962  * ( f61 * ( val(k+1,i,j)+val(k,i,j) ) &
963  + f62 * ( val(k+2,i,j)+val(k-1,i,j) ) &
964  + f63 * ( val(k+3,i,j)+val(k-2,i,j) ) ) &
965  + gsqrt(k,i,j) * num_diff(k,i,j)
966  enddo
967  enddo
968  !$acc end kernels
969  !$omp end do nowait
970 #ifdef DEBUG
971  k = iundef; i = iundef; j = iundef
972 #endif
973 
974  !$omp do OMP_SCHEDULE_
975  !$acc kernels
976  do j = jjs, jje
977  i = iis
978 #ifdef DEBUG
979 
980  call check( __line__, mom(ks,i ,j) )
981  call check( __line__, val(ks+1,i,j) )
982  call check( __line__, val(ks,i,j) )
983 
984  call check( __line__, mom(ks+1,i ,j) )
985  call check( __line__, val(ks+3,i,j) )
986  call check( __line__, val(ks+2,i,j) )
987 
988 #endif
989  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
990  ! The flux at KS-1 can be non-zero.
991  ! To reduce calculations, all the fluxes are set to zero.
992  flux(ks-1,i,j) = 0.0_rp
993 
994  vel = ( mom(ks,i,j) ) &
995  / ( f2h(ks,1,i_xyz) &
996  * dens(ks+1,i,j) &
997  + f2h(ks,2,i_xyz) &
998  * dens(ks,i,j) )
999  flux(ks,i,j) = j33g * vel &
1000  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) ) &
1001  + gsqrt(ks,i,j) * num_diff(ks,i,j)
1002  vel = ( mom(ke-1,i,j) ) &
1003  / ( f2h(ke-1,1,i_xyz) &
1004  * dens(ke,i,j) &
1005  + f2h(ke-1,2,i_xyz) &
1006  * dens(ke-1,i,j) )
1007  flux(ke-1,i,j) = j33g * vel &
1008  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) ) &
1009  + gsqrt(ke-1,i,j) * num_diff(ke-1,i,j)
1010 
1011  vel = ( mom(ks+1,i,j) ) &
1012  / ( f2h(ks+1,1,i_xyz) &
1013  * dens(ks+2,i,j) &
1014  + f2h(ks+1,2,i_xyz) &
1015  * dens(ks+1,i,j) )
1016  flux(ks+1,i,j) = j33g * vel &
1017  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
1018  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) ) &
1019  + gsqrt(ks+1,i,j) * num_diff(ks+1,i,j)
1020  vel = ( mom(ke-2,i,j) ) &
1021  / ( f2h(ke-2,1,i_xyz) &
1022  * dens(ke-1,i,j) &
1023  + f2h(ke-2,2,i_xyz) &
1024  * dens(ke-2,i,j) )
1025  flux(ke-2,i,j) = j33g * vel &
1026  * ( f41 * ( val(ke-1,i,j)+val(ke-2,i,j) ) &
1027  + f42 * ( val(ke,i,j)+val(ke-3,i,j) ) ) &
1028  + gsqrt(ke-2,i,j) * num_diff(ke-2,i,j)
1029 
1030  flux(ke,i,j) = 0.0_rp
1031  enddo
1032  !$acc end kernels
1033  !$omp end do nowait
1034 
1035  else
1036 
1037 
1038  !$omp do OMP_SCHEDULE_ collapse(2)
1039  !$acc kernels
1040  do j = jjs, jje
1041  do i = iis, iie
1042  do k = ks+2, ke-3
1043 #ifdef DEBUG
1044  call check( __line__, mom(k,i,j) )
1045  call check( __line__, mom(k,i+1,j) )
1046 
1047  call check( __line__, val(k,i,j) )
1048  call check( __line__, val(k+1,i,j) )
1049 
1050  call check( __line__, val(k-1,i,j) )
1051  call check( __line__, val(k+2,i,j) )
1052 
1053  call check( __line__, val(k-2,i,j) )
1054  call check( __line__, val(k+3,i,j) )
1055 
1056 #endif
1057  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i+1,j) ) ) &
1058  / ( f2h(k,1,i_uyz) &
1059  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i+1,j) ) &
1060  + f2h(k,2,i_uyz) &
1061  * 0.5_rp * ( dens(k,i,j)+dens(k,i+1,j) ) )
1062  flux(k,i,j) = j33g * vel &
1063  * ( f61 * ( val(k+1,i,j)+val(k,i,j) ) &
1064  + f62 * ( val(k+2,i,j)+val(k-1,i,j) ) &
1065  + f63 * ( val(k+3,i,j)+val(k-2,i,j) ) ) &
1066  + gsqrt(k,i,j) * num_diff(k,i,j)
1067  enddo
1068  enddo
1069  enddo
1070  !$acc end kernels
1071  !$omp end do nowait
1072 #ifdef DEBUG
1073  k = iundef; i = iundef; j = iundef
1074 #endif
1075 
1076  !$omp do OMP_SCHEDULE_ collapse(2)
1077  !$acc kernels
1078  do j = jjs, jje
1079  do i = iis, iie
1080 #ifdef DEBUG
1081 
1082  call check( __line__, mom(ks,i ,j) )
1083  call check( __line__, mom(ks,i+1,j) )
1084  call check( __line__, val(ks+1,i,j) )
1085  call check( __line__, val(ks,i,j) )
1086 
1087  call check( __line__, mom(ks+1,i ,j) )
1088  call check( __line__, mom(ks+1,i+1,j) )
1089  call check( __line__, val(ks+3,i,j) )
1090  call check( __line__, val(ks+2,i,j) )
1091 
1092 #endif
1093  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
1094  ! The flux at KS-1 can be non-zero.
1095  ! To reduce calculations, all the fluxes are set to zero.
1096  flux(ks-1,i,j) = 0.0_rp
1097 
1098  vel = ( 0.5_rp * ( mom(ks,i,j)+mom(ks,i+1,j) ) ) &
1099  / ( f2h(ks,1,i_uyz) &
1100  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i+1,j) ) &
1101  + f2h(ks,2,i_uyz) &
1102  * 0.5_rp * ( dens(ks,i,j)+dens(ks,i+1,j) ) )
1103  flux(ks,i,j) = j33g * vel &
1104  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) ) &
1105  + gsqrt(ks,i,j) * num_diff(ks,i,j)
1106  vel = ( 0.5_rp * ( mom(ke-1,i,j)+mom(ke-1,i+1,j) ) ) &
1107  / ( f2h(ke-1,1,i_uyz) &
1108  * 0.5_rp * ( dens(ke,i,j)+dens(ke,i+1,j) ) &
1109  + f2h(ke-1,2,i_uyz) &
1110  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i+1,j) ) )
1111  flux(ke-1,i,j) = j33g * vel &
1112  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) ) &
1113  + gsqrt(ke-1,i,j) * num_diff(ke-1,i,j)
1114 
1115  vel = ( 0.5_rp * ( mom(ks+1,i,j)+mom(ks+1,i+1,j) ) ) &
1116  / ( f2h(ks+1,1,i_uyz) &
1117  * 0.5_rp * ( dens(ks+2,i,j)+dens(ks+2,i+1,j) ) &
1118  + f2h(ks+1,2,i_uyz) &
1119  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i+1,j) ) )
1120  flux(ks+1,i,j) = j33g * vel &
1121  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
1122  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) ) &
1123  + gsqrt(ks+1,i,j) * num_diff(ks+1,i,j)
1124  vel = ( 0.5_rp * ( mom(ke-2,i,j)+mom(ke-2,i+1,j) ) ) &
1125  / ( f2h(ke-2,1,i_uyz) &
1126  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i+1,j) ) &
1127  + f2h(ke-2,2,i_uyz) &
1128  * 0.5_rp * ( dens(ke-2,i,j)+dens(ke-2,i+1,j) ) )
1129  flux(ke-2,i,j) = j33g * vel &
1130  * ( f41 * ( val(ke-1,i,j)+val(ke-2,i,j) ) &
1131  + f42 * ( val(ke,i,j)+val(ke-3,i,j) ) ) &
1132  + gsqrt(ke-2,i,j) * num_diff(ke-2,i,j)
1133 
1134  flux(ke,i,j) = 0.0_rp
1135  enddo
1136  enddo
1137  !$acc end kernels
1138  !$omp end do nowait
1139 
1140  end if
1141 
1142 
1143  !$acc end data
1144 
1145  !$omp end parallel
1146 #ifdef DEBUG
1147  k = iundef; i = iundef; j = iundef
1148 #endif
1149 
1150  return
1151  end subroutine atmos_dyn_fvm_fluxz_uyz_cd6
1152 
1153  !-----------------------------------------------------------------------------
1155  subroutine atmos_dyn_fvm_fluxj13_uyz_cd6( &
1156  flux, &
1157  mom, val, DENS, &
1158  GSQRT, J13G, MAPF, &
1159  CDZ, TwoD, &
1160  IIS, IIE, JJS, JJE )
1161  implicit none
1162 
1163  real(rp), intent(inout) :: flux (ka,ia,ja)
1164  real(rp), intent(in) :: mom (ka,ia,ja)
1165  real(rp), intent(in) :: val (ka,ia,ja)
1166  real(rp), intent(in) :: dens (ka,ia,ja)
1167  real(rp), intent(in) :: gsqrt (ka,ia,ja)
1168  real(rp), intent(in) :: j13g (ka,ia,ja)
1169  real(rp), intent(in) :: mapf ( ia,ja,2)
1170  real(rp), intent(in) :: cdz (ka)
1171  logical, intent(in) :: twod
1172  integer, intent(in) :: iis, iie, jjs, jje
1173 
1174  real(rp) :: vel
1175  integer :: k, i, j
1176  !---------------------------------------------------------------------------
1177 
1178  !$omp parallel default(none) private(i,j,k,vel) &
1179  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J13G,MAPF) &
1180  !$omp shared(GSQRT,CDZ,TwoD)
1181 
1182  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, J13G, MAPF, CDZ)
1183 
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+2, ke-3
1191  vel = ( f2h(k,1,i_uyz) &
1192  * mom(k+1,i,j) &
1193  + f2h(k,2,i_uyz) &
1194  * mom(k,i,j) ) &
1195  / ( f2h(k,1,i_uyz) &
1196  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i+1,j) ) &
1197  + f2h(k,2,i_uyz) &
1198  * 0.5_rp * ( dens(k,i,j)+dens(k,i+1,j) ) )
1199  vel = vel * j13g(k,i,j)
1200  flux(k,i,j) = vel / mapf(i,j,+2) &
1201  * ( f61 * ( val(k+1,i,j)+val(k,i,j) ) &
1202  + f62 * ( val(k+2,i,j)+val(k-1,i,j) ) &
1203  + f63 * ( val(k+3,i,j)+val(k-2,i,j) ) )
1204  enddo
1205  enddo
1206  enddo
1207  !$acc end kernels
1208  !$omp end do nowait
1209 
1210  !$omp do OMP_SCHEDULE_ collapse(2)
1211  !$acc kernels
1212  do j = jjs, jje
1213  do i = iis, iie
1214  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
1215  ! The flux at KS-1 can be non-zero.
1216  ! To reduce calculations, all the fluxes are set to zero.
1217  flux(ks-1,i,j) = 0.0_rp
1218 
1219  vel = ( f2h(ks,1,i_uyz) &
1220  * mom(ks+1,i,j) &
1221  + f2h(ks,2,i_uyz) &
1222  * mom(ks,i,j) ) &
1223  / ( f2h(ks,1,i_uyz) &
1224  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i+1,j) ) &
1225  + f2h(ks,2,i_uyz) &
1226  * 0.5_rp * ( dens(ks,i,j)+dens(ks,i+1,j) ) )
1227  vel = vel * j13g(ks,i,j)
1228  flux(ks,i,j) = vel / mapf(i,j,+2) &
1229  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) )
1230 
1231  vel = ( f2h(ke-1,1,i_uyz) &
1232  * mom(ke,i,j) &
1233  + f2h(ke-1,2,i_uyz) &
1234  * mom(ke-1,i,j) ) &
1235  / ( f2h(ke-1,1,i_uyz) &
1236  * 0.5_rp * ( dens(ke,i,j)+dens(ke,i+1,j) ) &
1237  + f2h(ke-1,2,i_uyz) &
1238  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i+1,j) ) )
1239  vel = vel * j13g(ke-1,i,j)
1240  flux(ke-1,i,j) = vel / mapf(i,j,+2) &
1241  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) )
1242 
1243  vel = ( f2h(ks+1,1,i_uyz) &
1244  * mom(ks+2,i,j) &
1245  + f2h(ks+1,2,i_uyz) &
1246  * mom(ks+1,i,j) ) &
1247  / ( f2h(ks+1,1,i_uyz) &
1248  * 0.5_rp * ( dens(ks+2,i,j)+dens(ks+2,i+1,j) ) &
1249  + f2h(ks+1,2,i_uyz) &
1250  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i+1,j) ) )
1251  vel = vel * j13g(ks+1,i,j)
1252  flux(ks+1,i,j) = vel / mapf(i,j,+2) &
1253  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
1254  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) )
1255 
1256  vel = ( f2h(ke-2,1,i_uyz) &
1257  * mom(ke-1,i,j) &
1258  + f2h(ke-2,2,i_uyz) &
1259  * mom(ke-2,i,j) ) &
1260  / ( f2h(ke-2,1,i_uyz) &
1261  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i+1,j) ) &
1262  + f2h(ke-2,2,i_uyz) &
1263  * 0.5_rp * ( dens(ke-2,i,j)+dens(ke-2,i+1,j) ) )
1264  vel = vel * j13g(ke-2,i,j)
1265  flux(ke-2,i,j) = vel / mapf(i,j,+2) &
1266  * ( f41 * ( val(ke-1,i,j)+val(ke-2,i,j) ) &
1267  + f42 * ( val(ke,i,j)+val(ke-3,i,j) ) )
1268 
1269  flux(ke ,i,j) = 0.0_rp
1270  enddo
1271  enddo
1272  !$acc end kernels
1273  !$omp end do nowait
1274 
1275 
1276 
1277  !$acc end data
1278 
1279  !$omp end parallel
1280  return
1281  end subroutine atmos_dyn_fvm_fluxj13_uyz_cd6
1282 
1283  !-----------------------------------------------------------------------------
1285  subroutine atmos_dyn_fvm_fluxj23_uyz_cd6( &
1286  flux, &
1287  mom, val, DENS, &
1288  GSQRT, J23G, MAPF, &
1289  CDZ, TwoD, &
1290  IIS, IIE, JJS, JJE )
1291  implicit none
1292 
1293  real(rp), intent(inout) :: flux (ka,ia,ja)
1294  real(rp), intent(in) :: mom (ka,ia,ja)
1295  real(rp), intent(in) :: val (ka,ia,ja)
1296  real(rp), intent(in) :: dens (ka,ia,ja)
1297  real(rp), intent(in) :: gsqrt (ka,ia,ja)
1298  real(rp), intent(in) :: j23g (ka,ia,ja)
1299  real(rp), intent(in) :: mapf ( ia,ja,2)
1300  real(rp), intent(in) :: cdz (ka)
1301  logical, intent(in) :: twod
1302  integer, intent(in) :: iis, iie, jjs, jje
1303 
1304  real(rp) :: vel
1305  integer :: k, i, j
1306  !---------------------------------------------------------------------------
1307 
1308  !$omp parallel default(none) private(i,j,k,vel) &
1309  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J23G,MAPF) &
1310  !$omp shared(GSQRT,CDZ,TwoD)
1311 
1312  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, J23G, MAPF, CDZ)
1313 
1314 
1315  if ( twod ) then
1316 
1317  !$omp do OMP_SCHEDULE_
1318  !$acc kernels
1319  do j = jjs, jje
1320  do k = ks+2, ke-3
1321  i = iis
1322  vel = ( f2h(k,1,i_xyz) &
1323  * 0.5_rp * ( mom(k+1,i,j)+mom(k+1,i,j-1) ) &
1324  + f2h(k,2,i_xyz) &
1325  * 0.5_rp * ( mom(k,i,j)+mom(k,i,j-1) ) ) &
1326  / ( f2h(k,1,i_xyz) &
1327  * dens(k+1,i,j) &
1328  + f2h(k,2,i_xyz) &
1329  * dens(k,i,j) )
1330  vel = vel * j23g(k,i,j)
1331  flux(k,i,j) = vel * ( f61 * ( val(k+1,i,j)+val(k,i,j) ) &
1332  + f62 * ( val(k+2,i,j)+val(k-1,i,j) ) &
1333  + f63 * ( val(k+3,i,j)+val(k-2,i,j) ) )
1334  enddo
1335  enddo
1336  !$acc end kernels
1337  !$omp end do nowait
1338 
1339  !$omp do OMP_SCHEDULE_
1340  !$acc kernels
1341  do j = jjs, jje
1342  i = iis
1343  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
1344  ! The flux at KS-1 can be non-zero.
1345  ! To reduce calculations, all the fluxes are set to zero.
1346  flux(ks-1,i,j) = 0.0_rp
1347 
1348  vel = ( f2h(ks,1,i_xyz) &
1349  * 0.5_rp * ( mom(ks+1,i,j)+mom(ks+1,i,j-1) ) &
1350  + f2h(ks,2,i_xyz) &
1351  * 0.5_rp * ( mom(ks,i,j)+mom(ks,i,j-1) ) ) &
1352  / ( f2h(ks,1,i_xyz) &
1353  * dens(ks+1,i,j) &
1354  + f2h(ks,2,i_xyz) &
1355  * dens(ks,i,j) )
1356  vel = vel * j23g(ks,i,j)
1357  flux(ks,i,j) = vel / mapf(i,j,+1) &
1358  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) )
1359 
1360  vel = ( f2h(ke-1,1,i_xyz) &
1361  * 0.5_rp * ( mom(ke,i,j)+mom(ke,i,j-1) ) &
1362  + f2h(ke-1,2,i_xyz) &
1363  * 0.5_rp * ( mom(ke-1,i,j)+mom(ke-1,i,j-1) ) ) &
1364  / ( f2h(ke-1,1,i_xyz) &
1365  * dens(ke,i,j) &
1366  + f2h(ke-1,2,i_xyz) &
1367  * dens(ke-1,i,j) )
1368  vel = vel * j23g(ke-1,i,j)
1369  flux(ke-1,i,j) = vel / mapf(i,j,+1) &
1370  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) )
1371 
1372  vel = ( f2h(ks+1,1,i_xyz) &
1373  * 0.5_rp * ( mom(ks+2,i,j)+mom(ks+2,i,j-1) ) &
1374  + f2h(ks+1,2,i_xyz) &
1375  * 0.5_rp * ( mom(ks+1,i,j)+mom(ks+1,i,j-1) ) ) &
1376  / ( f2h(ks+1,1,i_xyz) &
1377  * dens(ks+2,i,j) &
1378  + f2h(ks+1,2,i_xyz) &
1379  * dens(ks+1,i,j) )
1380  vel = vel * j23g(ks+1,i,j)
1381  flux(ks+1,i,j) = vel / mapf(i,j,+1) &
1382  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
1383  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) )
1384 
1385  vel = ( f2h(ke-2,1,i_xyz) &
1386  * 0.5_rp * ( mom(ke-1,i,j)+mom(ke-1,i,j-1) ) &
1387  + f2h(ke-2,2,i_xyz) &
1388  * 0.5_rp * ( mom(ke-2,i,j)+mom(ke-2,i,j-1) ) ) &
1389  / ( f2h(ke-2,1,i_xyz) &
1390  * dens(ke-1,i,j) &
1391  + f2h(ke-2,2,i_xyz) &
1392  * dens(ke-2,i,j) )
1393  vel = vel * j23g(ke-2,i,j)
1394  flux(ke-2,i,j) = vel / mapf(i,j,+1) &
1395  * ( f41 * ( val(ke-1,i,j)+val(ke-2,i,j) ) &
1396  + f42 * ( val(ke,i,j)+val(ke-3,i,j) ) )
1397 
1398  flux(ke ,i,j) = 0.0_rp
1399  enddo
1400  !$acc end kernels
1401  !$omp end do nowait
1402 
1403  else
1404 
1405 
1406  !$omp do OMP_SCHEDULE_ collapse(2)
1407  !$acc kernels
1408  do j = jjs, jje
1409  do i = iis, iie
1410  do k = ks+2, ke-3
1411  vel = ( f2h(k,1,i_uyz) &
1412  * 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) ) &
1413  + f2h(k,2,i_uyz) &
1414  * 0.25_rp * ( mom(k,i,j)+mom(k,i+1,j)+mom(k,i,j-1)+mom(k,i+1,j-1) ) ) &
1415  / ( f2h(k,1,i_uyz) &
1416  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i+1,j) ) &
1417  + f2h(k,2,i_uyz) &
1418  * 0.5_rp * ( dens(k,i,j)+dens(k,i+1,j) ) )
1419  vel = vel * j23g(k,i,j)
1420  flux(k,i,j) = vel / mapf(i,j,+1) &
1421  * ( f61 * ( val(k+1,i,j)+val(k,i,j) ) &
1422  + f62 * ( val(k+2,i,j)+val(k-1,i,j) ) &
1423  + f63 * ( val(k+3,i,j)+val(k-2,i,j) ) )
1424  enddo
1425  enddo
1426  enddo
1427  !$acc end kernels
1428  !$omp end do nowait
1429 
1430  !$omp do OMP_SCHEDULE_ collapse(2)
1431  !$acc kernels
1432  do j = jjs, jje
1433  do i = iis, iie
1434  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
1435  ! The flux at KS-1 can be non-zero.
1436  ! To reduce calculations, all the fluxes are set to zero.
1437  flux(ks-1,i,j) = 0.0_rp
1438 
1439  vel = ( f2h(ks,1,i_uyz) &
1440  * 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) ) &
1441  + f2h(ks,2,i_uyz) &
1442  * 0.25_rp * ( mom(ks,i,j)+mom(ks,i+1,j)+mom(ks,i,j-1)+mom(ks,i+1,j-1) ) ) &
1443  / ( f2h(ks,1,i_uyz) &
1444  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i+1,j) ) &
1445  + f2h(ks,2,i_uyz) &
1446  * 0.5_rp * ( dens(ks,i,j)+dens(ks,i+1,j) ) )
1447  vel = vel * j23g(ks,i,j)
1448  flux(ks,i,j) = vel / mapf(i,j,+1) &
1449  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) )
1450 
1451  vel = ( f2h(ke-1,1,i_uyz) &
1452  * 0.25_rp * ( mom(ke,i,j)+mom(ke,i+1,j)+mom(ke,i,j-1)+mom(ke,i+1,j-1) ) &
1453  + f2h(ke-1,2,i_uyz) &
1454  * 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) ) ) &
1455  / ( f2h(ke-1,1,i_uyz) &
1456  * 0.5_rp * ( dens(ke,i,j)+dens(ke,i+1,j) ) &
1457  + f2h(ke-1,2,i_uyz) &
1458  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i+1,j) ) )
1459  vel = vel * j23g(ke-1,i,j)
1460  flux(ke-1,i,j) = vel / mapf(i,j,+1) &
1461  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) )
1462 
1463  vel = ( f2h(ks+1,1,i_uyz) &
1464  * 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) ) &
1465  + f2h(ks+1,2,i_uyz) &
1466  * 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) ) ) &
1467  / ( f2h(ks+1,1,i_uyz) &
1468  * 0.5_rp * ( dens(ks+2,i,j)+dens(ks+2,i+1,j) ) &
1469  + f2h(ks+1,2,i_uyz) &
1470  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i+1,j) ) )
1471  vel = vel * j23g(ks+1,i,j)
1472  flux(ks+1,i,j) = vel / mapf(i,j,+1) &
1473  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
1474  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) )
1475 
1476  vel = ( f2h(ke-2,1,i_uyz) &
1477  * 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) ) &
1478  + f2h(ke-2,2,i_uyz) &
1479  * 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) ) ) &
1480  / ( f2h(ke-2,1,i_uyz) &
1481  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i+1,j) ) &
1482  + f2h(ke-2,2,i_uyz) &
1483  * 0.5_rp * ( dens(ke-2,i,j)+dens(ke-2,i+1,j) ) )
1484  vel = vel * j23g(ke-2,i,j)
1485  flux(ke-2,i,j) = vel / mapf(i,j,+1) &
1486  * ( f41 * ( val(ke-1,i,j)+val(ke-2,i,j) ) &
1487  + f42 * ( val(ke,i,j)+val(ke-3,i,j) ) )
1488 
1489  flux(ke ,i,j) = 0.0_rp
1490  enddo
1491  enddo
1492  !$acc end kernels
1493  !$omp end do nowait
1494 
1495 
1496  end if
1497 
1498 
1499  !$acc end data
1500 
1501  !$omp end parallel
1502  return
1503  end subroutine atmos_dyn_fvm_fluxj23_uyz_cd6
1504 
1505  !-----------------------------------------------------------------------------
1507  subroutine atmos_dyn_fvm_fluxx_uyz_cd6( &
1508  flux, &
1509  mom, val, DENS, &
1510  GSQRT, MAPF, &
1511  num_diff, &
1512  CDZ, TwoD, &
1513  IIS, IIE, JJS, JJE )
1514  implicit none
1515 
1516  real(rp), intent(inout) :: flux (ka,ia,ja)
1517  real(rp), intent(in) :: mom (ka,ia,ja)
1518  real(rp), intent(in) :: val (ka,ia,ja)
1519  real(rp), intent(in) :: dens (ka,ia,ja)
1520  real(rp), intent(in) :: gsqrt (ka,ia,ja)
1521  real(rp), intent(in) :: mapf ( ia,ja,2)
1522  real(rp), intent(in) :: num_diff(ka,ia,ja)
1523  real(rp), intent(in) :: cdz (ka)
1524  logical, intent(in) :: twod
1525  integer, intent(in) :: iis, iie, jjs, jje
1526 
1527  real(rp) :: vel
1528  integer :: k, i, j
1529  !---------------------------------------------------------------------------
1530 
1531  ! note that x-index is added by -1
1532 
1533 
1534 
1535  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
1536  !$omp private(vel) &
1537  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff)
1538  !$acc kernels
1539  do j = jjs, jje
1540  do i = iis, iie+1
1541  do k = ks, ke
1542 #ifdef DEBUG
1543  call check( __line__, mom(k,i ,j) )
1544  call check( __line__, mom(k,i-1,j) )
1545 
1546  call check( __line__, val(k,i-1,j) )
1547  call check( __line__, val(k,i,j) )
1548 
1549  call check( __line__, val(k,i-2,j) )
1550  call check( __line__, val(k,i+1,j) )
1551 
1552  call check( __line__, val(k,i-3,j) )
1553  call check( __line__, val(k,i+2,j) )
1554 
1555 #endif
1556  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i-1,j) ) ) &
1557  / ( dens(k,i,j) )
1558  flux(k,i-1,j) = gsqrt(k,i,j) / mapf(i,j,+2) * vel &
1559  * ( f61 * ( val(k,i,j)+val(k,i-1,j) ) &
1560  + f62 * ( val(k,i+1,j)+val(k,i-2,j) ) &
1561  + f63 * ( val(k,i+2,j)+val(k,i-3,j) ) ) &
1562  + gsqrt(k,i,j) * num_diff(k,i,j)
1563  enddo
1564  enddo
1565  enddo
1566  !$acc end kernels
1567 #ifdef DEBUG
1568  k = iundef; i = iundef; j = iundef
1569 #endif
1570 
1571 
1572 
1573  return
1574  end subroutine atmos_dyn_fvm_fluxx_uyz_cd6
1575 
1576  !-----------------------------------------------------------------------------
1578  subroutine atmos_dyn_fvm_fluxy_uyz_cd6( &
1579  flux, &
1580  mom, val, DENS, &
1581  GSQRT, MAPF, &
1582  num_diff, &
1583  CDZ, TwoD, &
1584  IIS, IIE, JJS, JJE )
1585  implicit none
1586 
1587  real(rp), intent(inout) :: flux (ka,ia,ja)
1588  real(rp), intent(in) :: mom (ka,ia,ja)
1589  real(rp), intent(in) :: val (ka,ia,ja)
1590  real(rp), intent(in) :: dens (ka,ia,ja)
1591  real(rp), intent(in) :: gsqrt (ka,ia,ja)
1592  real(rp), intent(in) :: mapf ( ia,ja,2)
1593  real(rp), intent(in) :: num_diff(ka,ia,ja)
1594  real(rp), intent(in) :: cdz (ka)
1595  logical, intent(in) :: twod
1596  integer, intent(in) :: iis, iie, jjs, jje
1597 
1598  real(rp) :: vel
1599  integer :: k, i, j
1600  !---------------------------------------------------------------------------
1601 
1602 
1603 
1604  if ( twod ) then
1605 
1606  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ &
1607  !$omp private(vel) &
1608  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff,TwoD)
1609  !$acc kernels
1610  do j = jjs-1, jje
1611  do k = ks, ke
1612  i = iis
1613 #ifdef DEBUG
1614  call check( __line__, mom(k,i ,j) )
1615 
1616  call check( __line__, val(k,i,j) )
1617  call check( __line__, val(k,i,j+1) )
1618 
1619  call check( __line__, val(k,i,j-1) )
1620  call check( __line__, val(k,i,j+2) )
1621 
1622  call check( __line__, val(k,i,j-2) )
1623  call check( __line__, val(k,i,j+3) )
1624 
1625 #endif
1626  vel = ( mom(k,i,j) ) &
1627  / ( 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) ) )
1628  flux(k,i,j) = gsqrt(k,i,j) / mapf(i,j,+1) * vel &
1629  * ( f61 * ( val(k,i,j+1)+val(k,i,j) ) &
1630  + f62 * ( val(k,i,j+2)+val(k,i,j-1) ) &
1631  + f63 * ( val(k,i,j+3)+val(k,i,j-2) ) ) &
1632  + gsqrt(k,i,j) * num_diff(k,i,j)
1633  enddo
1634  enddo
1635  !$acc end kernels
1636 #ifdef DEBUG
1637  k = iundef; i = iundef; j = iundef
1638 #endif
1639 
1640  else
1641 
1642 
1643  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
1644  !$omp private(vel) &
1645  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff)
1646  !$acc kernels
1647  do j = jjs-1, jje
1648  do i = iis, iie
1649  do k = ks, ke
1650 #ifdef DEBUG
1651  call check( __line__, mom(k,i ,j) )
1652  call check( __line__, mom(k,i-1,j) )
1653 
1654  call check( __line__, val(k,i,j) )
1655  call check( __line__, val(k,i,j+1) )
1656 
1657  call check( __line__, val(k,i,j-1) )
1658  call check( __line__, val(k,i,j+2) )
1659 
1660  call check( __line__, val(k,i,j-2) )
1661  call check( __line__, val(k,i,j+3) )
1662 
1663 #endif
1664  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i+1,j) ) ) &
1665  / ( 0.25_rp * ( dens(k,i,j)+dens(k,i+1,j)+dens(k,i,j+1)+dens(k,i+1,j+1) ) )
1666  flux(k,i,j) = gsqrt(k,i,j) / mapf(i,j,+1) * vel &
1667  * ( f61 * ( val(k,i,j+1)+val(k,i,j) ) &
1668  + f62 * ( val(k,i,j+2)+val(k,i,j-1) ) &
1669  + f63 * ( val(k,i,j+3)+val(k,i,j-2) ) ) &
1670  + gsqrt(k,i,j) * num_diff(k,i,j)
1671  enddo
1672  enddo
1673  enddo
1674  !$acc end kernels
1675 #ifdef DEBUG
1676  k = iundef; i = iundef; j = iundef
1677 #endif
1678 
1679 
1680  end if
1681 
1682 
1683  return
1684  end subroutine atmos_dyn_fvm_fluxy_uyz_cd6
1685 
1686 
1687 
1688  !-----------------------------------------------------------------------------
1690  subroutine atmos_dyn_fvm_fluxz_xvz_cd6( &
1691  flux, &
1692  mom, val, DENS, &
1693  GSQRT, J33G, &
1694  num_diff, &
1695  CDZ, TwoD, &
1696  IIS, IIE, JJS, JJE )
1697  implicit none
1698 
1699  real(rp), intent(inout) :: flux (ka,ia,ja)
1700  real(rp), intent(in) :: mom (ka,ia,ja)
1701  real(rp), intent(in) :: val (ka,ia,ja)
1702  real(rp), intent(in) :: dens (ka,ia,ja)
1703  real(rp), intent(in) :: gsqrt (ka,ia,ja)
1704  real(rp), intent(in) :: j33g
1705  real(rp), intent(in) :: num_diff(ka,ia,ja)
1706  real(rp), intent(in) :: cdz (ka)
1707  logical, intent(in) :: twod
1708  integer, intent(in) :: iis, iie, jjs, jje
1709 
1710  real(rp) :: vel
1711  integer :: k, i, j
1712  !---------------------------------------------------------------------------
1713 
1714  !$omp parallel default(none) private(i,j,k,vel) &
1715  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J33G,GSQRT,num_diff) &
1716  !$omp shared(CDZ,TwoD)
1717 
1718  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, num_diff, CDZ)
1719 
1720 
1721  !$omp do OMP_SCHEDULE_ collapse(2)
1722  !$acc kernels
1723  do j = jjs, jje
1724  do i = iis, iie
1725  do k = ks+2, ke-3
1726 #ifdef DEBUG
1727  call check( __line__, mom(k,i,j) )
1728  call check( __line__, mom(k,i,j+1) )
1729 
1730  call check( __line__, val(k,i,j) )
1731  call check( __line__, val(k+1,i,j) )
1732 
1733  call check( __line__, val(k-1,i,j) )
1734  call check( __line__, val(k+2,i,j) )
1735 
1736  call check( __line__, val(k-2,i,j) )
1737  call check( __line__, val(k+3,i,j) )
1738 
1739 #endif
1740  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i,j+1) ) ) &
1741  / ( f2h(k,1,i_xvz) &
1742  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i,j+1) ) &
1743  + f2h(k,2,i_xvz) &
1744  * 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) ) )
1745  flux(k,i,j) = j33g * vel &
1746  * ( f61 * ( val(k+1,i,j)+val(k,i,j) ) &
1747  + f62 * ( val(k+2,i,j)+val(k-1,i,j) ) &
1748  + f63 * ( val(k+3,i,j)+val(k-2,i,j) ) ) &
1749  + gsqrt(k,i,j) * num_diff(k,i,j)
1750  enddo
1751  enddo
1752  enddo
1753  !$acc end kernels
1754  !$omp end do nowait
1755 #ifdef DEBUG
1756  k = iundef; i = iundef; j = iundef
1757 #endif
1758 
1759  !$omp do OMP_SCHEDULE_ collapse(2)
1760  !$acc kernels
1761  do j = jjs, jje
1762  do i = iis, iie
1763 #ifdef DEBUG
1764 
1765  call check( __line__, mom(ks,i ,j) )
1766  call check( __line__, mom(ks,i,j+1) )
1767  call check( __line__, val(ks+1,i,j) )
1768  call check( __line__, val(ks,i,j) )
1769 
1770  call check( __line__, mom(ks+1,i ,j) )
1771  call check( __line__, mom(ks+1,i,j+1) )
1772  call check( __line__, val(ks+3,i,j) )
1773  call check( __line__, val(ks+2,i,j) )
1774 
1775 #endif
1776  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
1777  ! The flux at KS-1 can be non-zero.
1778  ! To reduce calculations, all the fluxes are set to zero.
1779  flux(ks-1,i,j) = 0.0_rp
1780 
1781  vel = ( 0.5_rp * ( mom(ks,i,j)+mom(ks,i,j+1) ) ) &
1782  / ( f2h(ks,1,i_xvz) &
1783  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i,j+1) ) &
1784  + f2h(ks,2,i_xvz) &
1785  * 0.5_rp * ( dens(ks,i,j)+dens(ks,i,j+1) ) )
1786  flux(ks,i,j) = j33g * vel &
1787  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) ) &
1788  + gsqrt(ks,i,j) * num_diff(ks,i,j)
1789  vel = ( 0.5_rp * ( mom(ke-1,i,j)+mom(ke-1,i,j+1) ) ) &
1790  / ( f2h(ke-1,1,i_xvz) &
1791  * 0.5_rp * ( dens(ke,i,j)+dens(ke,i,j+1) ) &
1792  + f2h(ke-1,2,i_xvz) &
1793  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i,j+1) ) )
1794  flux(ke-1,i,j) = j33g * vel &
1795  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) ) &
1796  + gsqrt(ke-1,i,j) * num_diff(ke-1,i,j)
1797 
1798  vel = ( 0.5_rp * ( mom(ks+1,i,j)+mom(ks+1,i,j+1) ) ) &
1799  / ( f2h(ks+1,1,i_xvz) &
1800  * 0.5_rp * ( dens(ks+2,i,j)+dens(ks+2,i,j+1) ) &
1801  + f2h(ks+1,2,i_xvz) &
1802  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i,j+1) ) )
1803  flux(ks+1,i,j) = j33g * vel &
1804  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
1805  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) ) &
1806  + gsqrt(ks+1,i,j) * num_diff(ks+1,i,j)
1807  vel = ( 0.5_rp * ( mom(ke-2,i,j)+mom(ke-2,i,j+1) ) ) &
1808  / ( f2h(ke-2,1,i_xvz) &
1809  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i,j+1) ) &
1810  + f2h(ke-2,2,i_xvz) &
1811  * 0.5_rp * ( dens(ke-2,i,j)+dens(ke-2,i,j+1) ) )
1812  flux(ke-2,i,j) = j33g * vel &
1813  * ( f41 * ( val(ke-1,i,j)+val(ke-2,i,j) ) &
1814  + f42 * ( val(ke,i,j)+val(ke-3,i,j) ) ) &
1815  + gsqrt(ke-2,i,j) * num_diff(ke-2,i,j)
1816 
1817  flux(ke,i,j) = 0.0_rp
1818  enddo
1819  enddo
1820  !$acc end kernels
1821  !$omp end do nowait
1822 
1823 
1824  !$acc end data
1825 
1826  !$omp end parallel
1827 #ifdef DEBUG
1828  k = iundef; i = iundef; j = iundef
1829 #endif
1830 
1831  return
1832  end subroutine atmos_dyn_fvm_fluxz_xvz_cd6
1833 
1834  !-----------------------------------------------------------------------------
1836  subroutine atmos_dyn_fvm_fluxj13_xvz_cd6( &
1837  flux, &
1838  mom, val, DENS, &
1839  GSQRT, J13G, MAPF, &
1840  CDZ, TwoD, &
1841  IIS, IIE, JJS, JJE )
1842  implicit none
1843 
1844  real(rp), intent(inout) :: flux (ka,ia,ja)
1845  real(rp), intent(in) :: mom (ka,ia,ja)
1846  real(rp), intent(in) :: val (ka,ia,ja)
1847  real(rp), intent(in) :: dens (ka,ia,ja)
1848  real(rp), intent(in) :: gsqrt (ka,ia,ja)
1849  real(rp), intent(in) :: j13g (ka,ia,ja)
1850  real(rp), intent(in) :: mapf ( ia,ja,2)
1851  real(rp), intent(in) :: cdz (ka)
1852  logical, intent(in) :: twod
1853  integer, intent(in) :: iis, iie, jjs, jje
1854 
1855  real(rp) :: vel
1856  integer :: k, i, j
1857  !---------------------------------------------------------------------------
1858 
1859  !$omp parallel default(none) private(i,j,k,vel) &
1860  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J13G,MAPF) &
1861  !$omp shared(GSQRT,CDZ,TwoD)
1862 
1863  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, J13G, MAPF, CDZ)
1864 
1865 
1866 
1867  !$omp do OMP_SCHEDULE_ collapse(2)
1868  !$acc kernels
1869  do j = jjs, jje
1870  do i = iis, iie
1871  do k = ks+2, ke-3
1872  vel = ( f2h(k,1,i_xvz) &
1873  * 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) ) &
1874  + f2h(k,2,i_xvz) &
1875  * 0.25_rp * ( mom(k,i,j)+mom(k,i-1,j)+mom(k,i,j+1)+mom(k,i-1,j+1) ) ) &
1876  / ( f2h(k,1,i_xvz) &
1877  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i,j+1) ) &
1878  + f2h(k,2,i_xvz) &
1879  * 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) ) )
1880  vel = vel * j13g(k,i,j)
1881  flux(k,i,j) = vel / mapf(i,j,+2) &
1882  * ( f61 * ( val(k+1,i,j)+val(k,i,j) ) &
1883  + f62 * ( val(k+2,i,j)+val(k-1,i,j) ) &
1884  + f63 * ( val(k+3,i,j)+val(k-2,i,j) ) )
1885  enddo
1886  enddo
1887  enddo
1888  !$acc end kernels
1889  !$omp end do nowait
1890 
1891  !$omp do OMP_SCHEDULE_ collapse(2)
1892  !$acc kernels
1893  do j = jjs, jje
1894  do i = iis, iie
1895  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
1896  ! The flux at KS-1 can be non-zero.
1897  ! To reduce calculations, all the fluxes are set to zero.
1898  flux(ks-1,i,j) = 0.0_rp
1899 
1900  vel = ( f2h(ks,1,i_xvz) &
1901  * 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) ) &
1902  + f2h(ks,2,i_xvz) &
1903  * 0.25_rp * ( mom(ks,i,j)+mom(ks,i-1,j)+mom(ks,i,j+1)+mom(ks,i-1,j+1) ) ) &
1904  / ( f2h(ks,1,i_xvz) &
1905  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i,j+1) ) &
1906  + f2h(ks,2,i_xvz) &
1907  * 0.5_rp * ( dens(ks,i,j)+dens(ks,i,j+1) ) )
1908  vel = vel * j13g(ks,i,j)
1909  flux(ks,i,j) = vel / mapf(i,j,+2) &
1910  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) )
1911 
1912  vel = ( f2h(ke-1,1,i_xvz) &
1913  * 0.25_rp * ( mom(ke,i,j)+mom(ke,i-1,j)+mom(ke,i,j+1)+mom(ke,i-1,j+1) ) &
1914  + f2h(ke-1,2,i_xvz) &
1915  * 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) ) ) &
1916  / ( f2h(ke-1,1,i_xvz) &
1917  * 0.5_rp * ( dens(ke,i,j)+dens(ke,i,j+1) ) &
1918  + f2h(ke-1,2,i_xvz) &
1919  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i,j+1) ) )
1920  vel = vel * j13g(ke-1,i,j)
1921  flux(ke-1,i,j) = vel / mapf(i,j,+2) &
1922  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) )
1923 
1924  vel = ( f2h(ks+1,1,i_xvz) &
1925  * 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) ) &
1926  + f2h(ks+1,2,i_xvz) &
1927  * 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) ) ) &
1928  / ( f2h(ks+1,1,i_xvz) &
1929  * 0.5_rp * ( dens(ks+2,i,j)+dens(ks+2,i,j+1) ) &
1930  + f2h(ks+1,2,i_xvz) &
1931  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i,j+1) ) )
1932  vel = vel * j13g(ks+1,i,j)
1933  flux(ks+1,i,j) = vel / mapf(i,j,+2) &
1934  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
1935  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) )
1936 
1937  vel = ( f2h(ke-2,1,i_xvz) &
1938  * 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) ) &
1939  + f2h(ke-2,2,i_xvz) &
1940  * 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) ) ) &
1941  / ( f2h(ke-2,1,i_xvz) &
1942  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i,j+1) ) &
1943  + f2h(ke-2,2,i_xvz) &
1944  * 0.5_rp * ( dens(ke-2,i,j)+dens(ke-2,i,j+1) ) )
1945  vel = vel * j13g(ke-2,i,j)
1946  flux(ke-2,i,j) = vel / mapf(i,j,+2) &
1947  * ( f41 * ( val(ke-1,i,j)+val(ke-2,i,j) ) &
1948  + f42 * ( val(ke,i,j)+val(ke-3,i,j) ) )
1949 
1950  flux(ke ,i,j) = 0.0_rp
1951  enddo
1952  enddo
1953  !$acc end kernels
1954  !$omp end do nowait
1955 
1956 
1957 
1958  !$acc end data
1959 
1960  !$omp end parallel
1961  return
1962  end subroutine atmos_dyn_fvm_fluxj13_xvz_cd6
1963 
1964  !-----------------------------------------------------------------------------
1966  subroutine atmos_dyn_fvm_fluxj23_xvz_cd6( &
1967  flux, &
1968  mom, val, DENS, &
1969  GSQRT, J23G, MAPF, &
1970  CDZ, TwoD, &
1971  IIS, IIE, JJS, JJE )
1972  implicit none
1973 
1974  real(rp), intent(inout) :: flux (ka,ia,ja)
1975  real(rp), intent(in) :: mom (ka,ia,ja)
1976  real(rp), intent(in) :: val (ka,ia,ja)
1977  real(rp), intent(in) :: dens (ka,ia,ja)
1978  real(rp), intent(in) :: gsqrt (ka,ia,ja)
1979  real(rp), intent(in) :: j23g (ka,ia,ja)
1980  real(rp), intent(in) :: mapf ( ia,ja,2)
1981  real(rp), intent(in) :: cdz (ka)
1982  logical, intent(in) :: twod
1983  integer, intent(in) :: iis, iie, jjs, jje
1984 
1985  real(rp) :: vel
1986  integer :: k, i, j
1987  !---------------------------------------------------------------------------
1988 
1989  !$omp parallel default(none) private(i,j,k,vel) &
1990  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J23G,MAPF) &
1991  !$omp shared(GSQRT,CDZ,TwoD)
1992 
1993  !$acc data copy(flux) copyin(mom, val, DENS, GSQRT, J23G, MAPF, CDZ)
1994 
1995 
1996 
1997  !$omp do OMP_SCHEDULE_ collapse(2)
1998  !$acc kernels
1999  do j = jjs, jje
2000  do i = iis, iie
2001  do k = ks+2, ke-3
2002  vel = ( f2h(k,1,i_xvz) &
2003  * mom(k+1,i,j) &
2004  + f2h(k,2,i_xvz) &
2005  * mom(k,i,j) ) &
2006  / ( f2h(k,1,i_xvz) &
2007  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i,j+1) ) &
2008  + f2h(k,2,i_xvz) &
2009  * 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) ) )
2010  vel = vel * j23g(k,i,j)
2011  flux(k,i,j) = vel / mapf(i,j,+1) &
2012  * ( f61 * ( val(k+1,i,j)+val(k,i,j) ) &
2013  + f62 * ( val(k+2,i,j)+val(k-1,i,j) ) &
2014  + f63 * ( val(k+3,i,j)+val(k-2,i,j) ) )
2015  enddo
2016  enddo
2017  enddo
2018  !$acc end kernels
2019  !$omp end do nowait
2020 
2021  !$omp do OMP_SCHEDULE_ collapse(2)
2022  !$acc kernels
2023  do j = jjs, jje
2024  do i = iis, iie
2025  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
2026  ! The flux at KS-1 can be non-zero.
2027  ! To reduce calculations, all the fluxes are set to zero.
2028  flux(ks-1,i,j) = 0.0_rp
2029 
2030  vel = ( f2h(ks,1,i_xvz) &
2031  * mom(ks+1,i,j) &
2032  + f2h(ks,2,i_xvz) &
2033  * mom(ks,i,j) ) &
2034  / ( f2h(ks,1,i_xvz) &
2035  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i,j+1) ) &
2036  + f2h(ks,2,i_xvz) &
2037  * 0.5_rp * ( dens(ks,i,j)+dens(ks,i,j+1) ) )
2038  vel = vel * j23g(ks,i,j)
2039  flux(ks,i,j) = vel / mapf(i,j,+1) &
2040  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) )
2041 
2042  vel = ( f2h(ke-1,1,i_xvz) &
2043  * mom(ke,i,j) &
2044  + f2h(ke-1,2,i_xvz) &
2045  * mom(ke-1,i,j) ) &
2046  / ( f2h(ke-1,1,i_xvz) &
2047  * 0.5_rp * ( dens(ke,i,j)+dens(ke,i,j+1) ) &
2048  + f2h(ke-1,2,i_xvz) &
2049  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i,j+1) ) )
2050  vel = vel * j23g(ke-1,i,j)
2051  flux(ke-1,i,j) = vel / mapf(i,j,+1) &
2052  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) )
2053 
2054  vel = ( f2h(ks+1,1,i_xvz) &
2055  * mom(ks+2,i,j) &
2056  + f2h(ks+1,2,i_xvz) &
2057  * mom(ks+1,i,j) ) &
2058  / ( f2h(ks+1,1,i_xvz) &
2059  * 0.5_rp * ( dens(ks+2,i,j)+dens(ks+2,i,j+1) ) &
2060  + f2h(ks+1,2,i_xvz) &
2061  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i,j+1) ) )
2062  vel = vel * j23g(ks+1,i,j)
2063  flux(ks+1,i,j) = vel / mapf(i,j,+1) &
2064  * ( f41 * ( val(ks+2,i,j)+val(ks+1,i,j) ) &
2065  + f42 * ( val(ks+3,i,j)+val(ks,i,j) ) )
2066 
2067  vel = ( f2h(ke-2,1,i_xvz) &
2068  * mom(ke-1,i,j) &
2069  + f2h(ke-2,2,i_xvz) &
2070  * mom(ke-2,i,j) ) &
2071  / ( f2h(ke-2,1,i_xvz) &
2072  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i,j+1) ) &
2073  + f2h(ke-2,2,i_xvz) &
2074  * 0.5_rp * ( dens(ke-2,i,j)+dens(ke-2,i,j+1) ) )
2075  vel = vel * j23g(ke-2,i,j)
2076  flux(ke-2,i,j) = vel / mapf(i,j,+1) &
2077  * ( f41 * ( val(ke-1,i,j)+val(ke-2,i,j) ) &
2078  + f42 * ( val(ke,i,j)+val(ke-3,i,j) ) )
2079 
2080  flux(ke ,i,j) = 0.0_rp
2081  enddo
2082  enddo
2083  !$acc end kernels
2084  !$omp end do nowait
2085 
2086 
2087 
2088  !$acc end data
2089 
2090  !$omp end parallel
2091  return
2092  end subroutine atmos_dyn_fvm_fluxj23_xvz_cd6
2093 
2094  !-----------------------------------------------------------------------------
2096  subroutine atmos_dyn_fvm_fluxx_xvz_cd6( &
2097  flux, &
2098  mom, val, DENS, &
2099  GSQRT, MAPF, &
2100  num_diff, &
2101  CDZ, TwoD, &
2102  IIS, IIE, JJS, JJE )
2103  implicit none
2104 
2105  real(rp), intent(inout) :: flux (ka,ia,ja)
2106  real(rp), intent(in) :: mom (ka,ia,ja)
2107  real(rp), intent(in) :: val (ka,ia,ja)
2108  real(rp), intent(in) :: dens (ka,ia,ja)
2109  real(rp), intent(in) :: gsqrt (ka,ia,ja)
2110  real(rp), intent(in) :: mapf ( ia,ja,2)
2111  real(rp), intent(in) :: num_diff(ka,ia,ja)
2112  real(rp), intent(in) :: cdz (ka)
2113  logical, intent(in) :: twod
2114  integer, intent(in) :: iis, iie, jjs, jje
2115 
2116  real(rp) :: vel
2117  integer :: k, i, j
2118  !---------------------------------------------------------------------------
2119 
2120 
2121 
2122  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
2123  !$omp private(vel) &
2124  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff)
2125  !$acc kernels
2126  do j = jjs, jje
2127  do i = iis-1, iie
2128  do k = ks, ke
2129 #ifdef DEBUG
2130  call check( __line__, mom(k,i ,j) )
2131  call check( __line__, mom(k,i,j-1) )
2132 
2133  call check( __line__, val(k,i,j) )
2134  call check( __line__, val(k,i+1,j) )
2135 
2136  call check( __line__, val(k,i-1,j) )
2137  call check( __line__, val(k,i+2,j) )
2138 
2139  call check( __line__, val(k,i-2,j) )
2140  call check( __line__, val(k,i+3,j) )
2141 
2142 #endif
2143  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i,j+1) ) ) &
2144  / ( 0.25_rp * ( dens(k,i,j)+dens(k,i+1,j)+dens(k,i,j+1)+dens(k,i+1,j+1) ) )
2145  flux(k,i,j) = gsqrt(k,i,j) / mapf(i,j,+2) * vel &
2146  * ( f61 * ( val(k,i+1,j)+val(k,i,j) ) &
2147  + f62 * ( val(k,i+2,j)+val(k,i-1,j) ) &
2148  + f63 * ( val(k,i+3,j)+val(k,i-2,j) ) ) &
2149  + gsqrt(k,i,j) * num_diff(k,i,j)
2150  enddo
2151  enddo
2152  enddo
2153  !$acc end kernels
2154 #ifdef DEBUG
2155  k = iundef; i = iundef; j = iundef
2156 #endif
2157 
2158 
2159 
2160  return
2161  end subroutine atmos_dyn_fvm_fluxx_xvz_cd6
2162 
2163  !-----------------------------------------------------------------------------
2165  subroutine atmos_dyn_fvm_fluxy_xvz_cd6( &
2166  flux, &
2167  mom, val, DENS, &
2168  GSQRT, MAPF, &
2169  num_diff, &
2170  CDZ, TwoD, &
2171  IIS, IIE, JJS, JJE )
2172  implicit none
2173 
2174  real(rp), intent(inout) :: flux (ka,ia,ja)
2175  real(rp), intent(in) :: mom (ka,ia,ja)
2176  real(rp), intent(in) :: val (ka,ia,ja)
2177  real(rp), intent(in) :: dens (ka,ia,ja)
2178  real(rp), intent(in) :: gsqrt (ka,ia,ja)
2179  real(rp), intent(in) :: mapf ( ia,ja,2)
2180  real(rp), intent(in) :: num_diff(ka,ia,ja)
2181  real(rp), intent(in) :: cdz (ka)
2182  logical, intent(in) :: twod
2183  integer, intent(in) :: iis, iie, jjs, jje
2184 
2185  real(rp) :: vel
2186  integer :: k, i, j
2187  !---------------------------------------------------------------------------
2188 
2189  ! note that y-index is added by -1
2190 
2191 
2192 
2193  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
2194  !$omp private(vel) &
2195  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff)
2196  !$acc kernels
2197  do j = jjs, jje+1
2198  do i = iis, iie
2199  do k = ks, ke
2200 #ifdef DEBUG
2201  call check( __line__, mom(k,i ,j) )
2202  call check( __line__, mom(k,i,j-1) )
2203 
2204  call check( __line__, val(k,i,j-1) )
2205  call check( __line__, val(k,i,j) )
2206 
2207  call check( __line__, val(k,i,j-2) )
2208  call check( __line__, val(k,i,j+1) )
2209 
2210  call check( __line__, val(k,i,j-3) )
2211  call check( __line__, val(k,i,j+2) )
2212 
2213 #endif
2214  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i,j-1) ) ) &
2215  / ( dens(k,i,j) )
2216  flux(k,i,j-1) = gsqrt(k,i,j) / mapf(i,j,+1) * vel &
2217  * ( f61 * ( val(k,i,j)+val(k,i,j-1) ) &
2218  + f62 * ( val(k,i,j+1)+val(k,i,j-2) ) &
2219  + f63 * ( val(k,i,j+2)+val(k,i,j-3) ) ) &
2220  + gsqrt(k,i,j) * num_diff(k,i,j)
2221  enddo
2222  enddo
2223  enddo
2224  !$acc end kernels
2225 #ifdef DEBUG
2226  k = iundef; i = iundef; j = iundef
2227 #endif
2228 
2229 
2230 
2231  return
2232  end subroutine atmos_dyn_fvm_fluxy_xvz_cd6
2233 
2234 
2235 
2236 
2237 
2238 
2239 
2241 
2242 !--
2243 ! vi:set readonly sw=4 ts=8
2244 !
2245 !Local Variables:
2246 !mode: f90
2247 !buffer-read-only: t
2248 !End:
2249 !
2250 !++
scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxj13_xvz_cd6
subroutine, public atmos_dyn_fvm_fluxj13_xvz_cd6(flux, mom, val, DENS, GSQRT, J13G, MAPF, CDZ, TwoD, IIS, IIE, JJS, JJE)
calculation J13-flux at XVZ
Definition: scale_atmos_dyn_fvm_flux_cd6.F90:1842
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_cd6::atmos_dyn_fvm_fluxj23_xyw_cd6
subroutine, public atmos_dyn_fvm_fluxj23_xyw_cd6(flux, mom, val, DENS, GSQRT, J23G, MAPF, CDZ, TwoD, IIS, IIE, JJS, JJE)
calculation J23-flux at XYW
Definition: scale_atmos_dyn_fvm_flux_cd6.F90:637
scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxy_xyz_cd6
subroutine, public atmos_dyn_fvm_fluxy_xyz_cd6(flux, mflx, val, GSQRT, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation Y-flux at XYZ
Definition: scale_atmos_dyn_fvm_flux_cd6.F90:358
scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxx_xyz_cd6
subroutine, public atmos_dyn_fvm_fluxx_xyz_cd6(flux, mflx, val, GSQRT, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation X-flux at XYZ
Definition: scale_atmos_dyn_fvm_flux_cd6.F90:299
scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxx_uyz_cd6
subroutine, public atmos_dyn_fvm_fluxx_uyz_cd6(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_cd6.F90:1514
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_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxj13_uyz_cd6
subroutine, public atmos_dyn_fvm_fluxj13_uyz_cd6(flux, mom, val, DENS, GSQRT, J13G, MAPF, CDZ, TwoD, IIS, IIE, JJS, JJE)
calculation J13-flux at UYZ
Definition: scale_atmos_dyn_fvm_flux_cd6.F90:1161
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_cd6::atmos_dyn_fvm_fluxy_xvz_cd6
subroutine, public atmos_dyn_fvm_fluxy_xvz_cd6(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_cd6.F90:2172
scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxx_xyw_cd6
subroutine, public atmos_dyn_fvm_fluxx_xyw_cd6(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_cd6.F90:725
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:35
scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxz_xvz_cd6
subroutine, public atmos_dyn_fvm_fluxz_xvz_cd6(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_cd6.F90:1697
scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxy_xyw_cd6
subroutine, public atmos_dyn_fvm_fluxy_xyw_cd6(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_cd6.F90:818
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_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_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxy_uyz_cd6
subroutine, public atmos_dyn_fvm_fluxy_uyz_cd6(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_cd6.F90:1585
scale_atmos_grid_cartesc_index::i_xvz
integer, public i_xvz
Definition: scale_atmos_grid_cartesC_index.F90:96
scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxj13_xyw_cd6
subroutine, public atmos_dyn_fvm_fluxj13_xyw_cd6(flux, mom, val, DENS, GSQRT, J13G, MAPF, CDZ, TwoD, IIS, IIE, JJS, JJE)
calculation J13-flux at XYW
Definition: scale_atmos_dyn_fvm_flux_cd6.F90:551
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
scale_atmos_dyn_fvm_flux_cd6
module scale_atmos_dyn_fvm_flux_cd6
Definition: scale_atmos_dyn_fvm_flux_cd6.F90:16
scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxj23_xvz_cd6
subroutine, public atmos_dyn_fvm_fluxj23_xvz_cd6(flux, mom, val, DENS, GSQRT, J23G, MAPF, CDZ, TwoD, IIS, IIE, JJS, JJE)
calculation J23-flux at XVZ
Definition: scale_atmos_dyn_fvm_flux_cd6.F90:1972
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxx_xvz_cd6
subroutine, public atmos_dyn_fvm_fluxx_xvz_cd6(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_cd6.F90:2103
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_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_flux_valuew_z_cd6
subroutine, public atmos_dyn_fvm_flux_valuew_z_cd6(valW, mflx, val, GSQRT, CDZ)
value at XYW
Definition: scale_atmos_dyn_fvm_flux_cd6.F90:106
scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxj23_uyz_cd6
subroutine, public atmos_dyn_fvm_fluxj23_uyz_cd6(flux, mom, val, DENS, GSQRT, J23G, MAPF, CDZ, TwoD, IIS, IIE, JJS, JJE)
calculation J23-flux at UYZ
Definition: scale_atmos_dyn_fvm_flux_cd6.F90:1291
scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxz_uyz_cd6
subroutine, public atmos_dyn_fvm_fluxz_uyz_cd6(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_cd6.F90:912
scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxz_xyz_cd6
subroutine, public atmos_dyn_fvm_fluxz_xyz_cd6(flux, mflx, val, GSQRT, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation z-flux at XYZ
Definition: scale_atmos_dyn_fvm_flux_cd6.F90:178
scale_atmos_dyn_fvm_flux_cd6::atmos_dyn_fvm_fluxz_xyw_cd6
subroutine, public atmos_dyn_fvm_fluxz_xyw_cd6(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_cd6.F90:420