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