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