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