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