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