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