SCALE-RM
scale_atmos_dyn_fvm_flux_ud3.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 :: F31 = -1.0_rp/12.0_rp
91  real(RP), parameter :: F32 = 7.0_rp/12.0_rp
92  real(RP), parameter :: F33 = 3.0_rp/12.0_rp
93 
94 
95 
96 
97 contains
98  !-----------------------------------------------------------------------------
100 !OCL SERIAL
101  subroutine atmos_dyn_fvm_flux_valuew_z_ud3( &
102  valW, &
103  mflx, val, GSQRT, &
104  CDZ )
105  implicit none
106 
107  real(RP), intent(out) :: valw (ka)
108  real(RP), intent(in) :: mflx (ka)
109  real(RP), intent(in) :: val (ka)
110  real(RP), intent(in) :: gsqrt(ka)
111  real(RP), intent(in) :: cdz (ka)
112 
113  integer :: k
114  !---------------------------------------------------------------------------
115 
116  do k = ks+1, ke-2
117 #ifdef DEBUG
118  call check( __line__, mflx(k) )
119 
120  call check( __line__, val(k) )
121  call check( __line__, val(k+1) )
122 
123  call check( __line__, val(k-1) )
124  call check( __line__, val(k+2) )
125 
126 #endif
127  valw(k) = ( f31 * ( val(k+2)+val(k-1) ) + f32 * ( val(k+1)+val(k) ) ) &
128  - ( f31 * ( val(k+2)-val(k-1) ) + f33 * ( val(k+1)-val(k) ) ) * sign(1.0_rp,mflx(k))
129  enddo
130 #ifdef DEBUG
131  k = iundef
132 #endif
133 
134 #ifdef DEBUG
135 
136  call check( __line__, mflx(ks) )
137  call check( __line__, val(ks ) )
138  call check( __line__, val(ks+1) )
139  call check( __line__, mflx(ke-1) )
140  call check( __line__, val(ke ) )
141  call check( __line__, val(ke-1) )
142 
143 #endif
144 
145 
146 
147  valw(ks) = f2 * ( val(ks+1)+val(ks) )
148  valw(ke-1) = f2 * ( val(ke)+val(ke-1) )
149 
150 
151  return
152  end subroutine atmos_dyn_fvm_flux_valuew_z_ud3
153 
154  !-----------------------------------------------------------------------------
156  subroutine atmos_dyn_fvm_fluxz_xyz_ud3( &
157  flux, &
158  mflx, val, GSQRT, &
159  num_diff, &
160  CDZ, &
161  IIS, IIE, JJS, JJE )
162  implicit none
163 
164  real(RP), intent(inout) :: flux (ka,ia,ja)
165  real(RP), intent(in) :: mflx (ka,ia,ja)
166  real(RP), intent(in) :: val (ka,ia,ja)
167  real(RP), intent(in) :: gsqrt (ka,ia,ja)
168  real(RP), intent(in) :: num_diff(ka,ia,ja)
169  real(RP), intent(in) :: cdz (ka)
170  integer, intent(in) :: iis, iie, jjs, jje
171 
172  real(RP) :: vel
173  integer :: k, i, j
174  !---------------------------------------------------------------------------
175 
176  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
177  !$omp private(vel) &
178  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mflx,val,flux,GSQRT,num_diff)
179  do j = jjs, jje
180  do i = iis, iie
181  do k = ks+1, ke-2
182 #ifdef DEBUG
183  call check( __line__, mflx(k,i,j) )
184 
185  call check( __line__, val(k,i,j) )
186  call check( __line__, val(k+1,i,j) )
187 
188  call check( __line__, val(k-1,i,j) )
189  call check( __line__, val(k+2,i,j) )
190 
191 #endif
192  vel = mflx(k,i,j)
193  flux(k,i,j) = vel &
194  * ( ( f31 * ( val(k+2,i,j)+val(k-1,i,j) ) + f32 * ( val(k+1,i,j)+val(k,i,j) ) ) &
195  - ( f31 * ( val(k+2,i,j)-val(k-1,i,j) ) + f33 * ( val(k+1,i,j)-val(k,i,j) ) ) * sign(1.0_rp,vel) ) &
196  + gsqrt(k,i,j) * num_diff(k,i,j)
197  enddo
198  enddo
199  enddo
200 #ifdef DEBUG
201  k = iundef; i = iundef; j = iundef
202 #endif
203 
204  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
205  !$omp private(vel) &
206  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mflx,val,flux,GSQRT,num_diff)
207  do j = jjs, jje
208  do i = iis, iie
209 #ifdef DEBUG
210 
211  call check( __line__, mflx(ks,i,j) )
212  call check( __line__, val(ks ,i,j) )
213  call check( __line__, val(ks+1,i,j) )
214  call check( __line__, mflx(ke-1,i,j) )
215  call check( __line__, val(ke ,i,j) )
216  call check( __line__, val(ke-1,i,j) )
217 
218 #endif
219  flux(ks-1,i,j) = 0.0_rp
220 
221  vel = mflx(ks,i,j)
222  flux(ks,i,j) = vel &
223  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) ) &
224  + gsqrt(ks,i,j) * num_diff(ks,i,j)
225  vel = mflx(ke-1,i,j)
226  flux(ke-1,i,j) = vel &
227  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) ) &
228  + gsqrt(ke-1,i,j) * num_diff(ke-1,i,j)
229 
230  flux(ke ,i,j) = 0.0_rp
231  enddo
232  enddo
233 #ifdef DEBUG
234  k = iundef; i = iundef; j = iundef
235 #endif
236 
237  return
238  end subroutine atmos_dyn_fvm_fluxz_xyz_ud3
239 
240  !-----------------------------------------------------------------------------
242  subroutine atmos_dyn_fvm_fluxx_xyz_ud3( &
243  flux, &
244  mflx, val, GSQRT, &
245  num_diff, &
246  CDZ, &
247  IIS, IIE, JJS, JJE )
248  implicit none
249 
250  real(RP), intent(inout) :: flux (ka,ia,ja)
251  real(RP), intent(in) :: mflx (ka,ia,ja)
252  real(RP), intent(in) :: val (ka,ia,ja)
253  real(RP), intent(in) :: gsqrt (ka,ia,ja)
254  real(RP), intent(in) :: num_diff(ka,ia,ja)
255  real(RP), intent(in) :: cdz(ka)
256  integer, intent(in) :: iis, iie, jjs, jje
257 
258  real(RP) :: vel
259  integer :: k, i, j
260  !---------------------------------------------------------------------------
261 
262  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
263  !$omp private(vel) &
264  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mflx,val,flux,GSQRT,num_diff)
265  do j = jjs, jje
266  do i = iis-1, iie
267  do k = ks, ke
268 #ifdef DEBUG
269  call check( __line__, mflx(k,i,j) )
270 
271  call check( __line__, val(k,i,j) )
272  call check( __line__, val(k,i+1,j) )
273 
274  call check( __line__, val(k,i-1,j) )
275  call check( __line__, val(k,i+2,j) )
276 
277 #endif
278  vel = mflx(k,i,j)
279  flux(k,i,j) = vel &
280  * ( ( f31 * ( val(k,i+2,j)+val(k,i-1,j) ) + f32 * ( val(k,i+1,j)+val(k,i,j) ) ) &
281  - ( f31 * ( val(k,i+2,j)-val(k,i-1,j) ) + f33 * ( val(k,i+1,j)-val(k,i,j) ) ) * sign(1.0_rp,vel) ) &
282  + gsqrt(k,i,j) * num_diff(k,i,j)
283  enddo
284  enddo
285  enddo
286 #ifdef DEBUG
287  k = iundef; i = iundef; j = iundef
288 #endif
289 
290  return
291  end subroutine atmos_dyn_fvm_fluxx_xyz_ud3
292 
293  !-----------------------------------------------------------------------------
295  subroutine atmos_dyn_fvm_fluxy_xyz_ud3( &
296  flux, &
297  mflx, val, GSQRT, &
298  num_diff, &
299  CDZ, &
300  IIS, IIE, JJS, JJE )
301  implicit none
302 
303  real(RP), intent(inout) :: flux (ka,ia,ja)
304  real(RP), intent(in) :: mflx (ka,ia,ja)
305  real(RP), intent(in) :: val (ka,ia,ja)
306  real(RP), intent(in) :: gsqrt (ka,ia,ja)
307  real(RP), intent(in) :: num_diff(ka,ia,ja)
308  real(RP), intent(in) :: cdz(ka)
309  integer, intent(in) :: iis, iie, jjs, jje
310 
311  real(RP) :: vel
312  integer :: k, i, j
313  !---------------------------------------------------------------------------
314 
315  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
316  !$omp private(vel) &
317  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mflx,val,flux,GSQRT,num_diff)
318  do j = jjs-1, jje
319  do i = iis, iie
320  do k = ks, ke
321 #ifdef DEBUG
322  call check( __line__, mflx(k,i,j) )
323 
324  call check( __line__, val(k,i,j) )
325  call check( __line__, val(k,i,j+1) )
326 
327  call check( __line__, val(k,i,j-1) )
328  call check( __line__, val(k,i,j+2) )
329 
330 #endif
331  vel = mflx(k,i,j)
332  flux(k,i,j) = vel &
333  * ( ( f31 * ( val(k,i,j+2)+val(k,i,j-1) ) + f32 * ( val(k,i,j+1)+val(k,i,j) ) ) &
334  - ( f31 * ( val(k,i,j+2)-val(k,i,j-1) ) + f33 * ( val(k,i,j+1)-val(k,i,j) ) ) * sign(1.0_rp,vel) ) &
335  + gsqrt(k,i,j) * num_diff(k,i,j)
336  enddo
337  enddo
338  enddo
339 #ifdef DEBUG
340  k = iundef; i = iundef; j = iundef
341 #endif
342 
343  return
344  end subroutine atmos_dyn_fvm_fluxy_xyz_ud3
345 
346 
347  !-----------------------------------------------------------------------------
349  subroutine atmos_dyn_fvm_fluxz_xyw_ud3( &
350  flux, &
351  mom, val, DENS, &
352  GSQRT, J33G, &
353  num_diff, &
354  CDZ, FDZ, &
355  dtrk, &
356  IIS, IIE, JJS, JJE )
357  implicit none
358 
359  real(RP), intent(inout) :: flux (ka,ia,ja)
360  real(RP), intent(in) :: mom (ka,ia,ja)
361  real(RP), intent(in) :: val (ka,ia,ja)
362  real(RP), intent(in) :: dens (ka,ia,ja)
363  real(RP), intent(in) :: gsqrt (ka,ia,ja)
364  real(RP), intent(in) :: j33g
365  real(RP), intent(in) :: num_diff(ka,ia,ja)
366  real(RP), intent(in) :: cdz (ka)
367  real(RP), intent(in) :: fdz (ka-1)
368  real(RP), intent(in) :: dtrk
369  integer, intent(in) :: iis, iie, jjs, jje
370 
371  real(RP) :: vel
372  integer :: k, i, j
373  !---------------------------------------------------------------------------
374 
375  ! note than z-index is added by -1
376 
377  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
378  !$omp private(vel) &
379  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,flux,J33G,GSQRT,num_diff,DENS)
380  do j = jjs, jje
381  do i = iis, iie
382  do k = ks+2, ke-1
383 #ifdef DEBUG
384  call check( __line__, mom(k-1,i,j) )
385  call check( __line__, mom(k ,i,j) )
386 
387  call check( __line__, val(k-1,i,j) )
388  call check( __line__, val(k,i,j) )
389 
390  call check( __line__, val(k-2,i,j) )
391  call check( __line__, val(k+1,i,j) )
392 
393 #endif
394  vel = ( 0.5_rp * ( mom(k-1,i,j) &
395  + mom(k,i,j) ) ) &
396  / dens(k,i,j)
397  flux(k-1,i,j) = j33g * vel &
398  * ( ( f31 * ( val(k+1,i,j)+val(k-2,i,j) ) + f32 * ( val(k,i,j)+val(k-1,i,j) ) ) &
399  - ( f31 * ( val(k+1,i,j)-val(k-2,i,j) ) + f33 * ( val(k,i,j)-val(k-1,i,j) ) ) * sign(1.0_rp,vel) ) &
400  + gsqrt(k,i,j) * num_diff(k,i,j)
401  enddo
402  enddo
403  enddo
404 #ifdef DEBUG
405  k = iundef; i = iundef; j = iundef
406 #endif
407 
408  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
409  !$omp private(vel) &
410  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J33G,GSQRT,num_diff,FDZ,dtrk)
411  do j = jjs, jje
412  do i = iis, iie
413 #ifdef DEBUG
414 
415  call check( __line__, val(ks ,i,j) )
416  call check( __line__, val(ks+1,i,j) )
417 
418 
419 #endif
420  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS.
421  ! The flux at KS can be non-zero.
422  ! To reduce calculations, all the fluxes are set to zero.
423  flux(ks-1,i,j) = 0.0_rp ! k = KS
424 
425  vel = ( 0.5_rp * ( mom(ks,i,j) &
426  + mom(ks+1,i,j) ) ) &
427  / dens(ks+1,i,j)
428  flux(ks,i,j) = j33g * vel &
429  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) ) &
430  + gsqrt(ks+1,i,j) * num_diff(ks+1,i,j) ! k = KS+1
431 
432 
433 
434  flux(ke-1,i,j) = 0.0_rp ! k = KE
435  flux(ke ,i,j) = 0.0_rp ! k = KE+1
436  enddo
437  enddo
438 
439  return
440  end subroutine atmos_dyn_fvm_fluxz_xyw_ud3
441 
442 
443  !-----------------------------------------------------------------------------
445  subroutine atmos_dyn_fvm_fluxj13_xyw_ud3( &
446  flux, &
447  mom, val, DENS, &
448  GSQRT, J13G, MAPF, &
449  CDZ, &
450  IIS, IIE, JJS, JJE )
451  implicit none
452 
453  real(RP), intent(inout) :: flux (ka,ia,ja)
454  real(RP), intent(in) :: mom (ka,ia,ja)
455  real(RP), intent(in) :: val (ka,ia,ja)
456  real(RP), intent(in) :: dens (ka,ia,ja)
457  real(RP), intent(in) :: gsqrt (ka,ia,ja)
458  real(RP), intent(in) :: j13g (ka,ia,ja)
459  real(RP), intent(in) :: mapf ( ia,ja,2)
460  real(RP), intent(in) :: cdz (ka)
461  integer, intent(in) :: iis, iie, jjs, jje
462 
463  real(RP) :: vel
464  integer :: k, i, j
465  !---------------------------------------------------------------------------
466 
467  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
468  !$omp private(vel) &
469  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J13G,MAPF)
470  do j = jjs, jje
471  do i = iis, iie
472  do k = ks+2, ke-1
473  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i-1,j) ) ) &
474  / dens(k,i,j)
475  flux(k-1,i,j) = j13g(k,i,j) / mapf(i,j,+2) * vel &
476  * ( ( f31 * ( val(k+1,i,j)+val(k-2,i,j) ) + f32 * ( val(k,i,j)+val(k-1,i,j) ) ) &
477  - ( f31 * ( val(k+1,i,j)-val(k-2,i,j) ) + f33 * ( val(k,i,j)-val(k-1,i,j) ) ) * sign(1.0_rp,vel) )
478  enddo
479  enddo
480  enddo
481 
482  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
483  !$omp private(vel) &
484  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J13G,MAPF)
485  do j = jjs, jje
486  do i = iis, iie
487  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS.
488  ! The flux at KS can be non-zero.
489  ! To reduce calculations, all the fluxes are set to zero.
490  flux(ks-1,i,j) = 0.0_rp ! k = KS
491 
492  ! physically incorrect but for numerical stability
493  vel = ( ( 0.5_rp * ( mom(ks+1,i,j)+mom(ks+1,i-1,j) ) ) / dens(ks+1,i,j) &
494  + ( 0.5_rp * ( mom(ks,i,j)+mom(ks,i-1,j) ) ) / dens(ks ,i,j) ) * 0.5_rp
495 ! vel = ( 0.5_RP * ( mom(KS+1,i,j)+mom(KS+1,i-1,j) ) ) &
496 ! / DENS(KS+1,i,j)
497  flux(ks,i,j) = j13g(ks+1,i,j) / mapf(i,j,+2) * vel &
498  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) ) ! k = KS+1
499 
500 
501  flux(ke-1,i,j) = 0.0_rp
502  enddo
503  enddo
504 
505  return
506  end subroutine atmos_dyn_fvm_fluxj13_xyw_ud3
507 
508  !-----------------------------------------------------------------------------
510  subroutine atmos_dyn_fvm_fluxj23_xyw_ud3( &
511  flux, &
512  mom, val, DENS, &
513  GSQRT, J23G, MAPF, &
514  CDZ, &
515  IIS, IIE, JJS, JJE )
516  implicit none
517 
518  real(RP), intent(inout) :: flux (ka,ia,ja)
519  real(RP), intent(in) :: mom (ka,ia,ja)
520  real(RP), intent(in) :: val (ka,ia,ja)
521  real(RP), intent(in) :: dens (ka,ia,ja)
522  real(RP), intent(in) :: gsqrt (ka,ia,ja)
523  real(RP), intent(in) :: j23g (ka,ia,ja)
524  real(RP), intent(in) :: mapf ( ia,ja,2)
525  real(RP), intent(in) :: cdz (ka)
526  integer, intent(in) :: iis, iie, jjs, jje
527 
528  real(RP) :: vel
529  integer :: k, i, j
530  !---------------------------------------------------------------------------
531 
532  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
533  !$omp private(vel) &
534  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J23G,MAPF)
535  do j = jjs, jje
536  do i = iis, iie
537  do k = ks+2, ke-1
538  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i,j-1) ) ) &
539  / dens(k,i,j)
540  flux(k-1,i,j) = j23g(k,i,j) / mapf(i,j,+1) * vel &
541  * ( ( f31 * ( val(k+1,i,j)+val(k-2,i,j) ) + f32 * ( val(k,i,j)+val(k-1,i,j) ) ) &
542  - ( f31 * ( val(k+1,i,j)-val(k-2,i,j) ) + f33 * ( val(k,i,j)-val(k-1,i,j) ) ) * sign(1.0_rp,vel) )
543  enddo
544  enddo
545  enddo
546 
547  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
548  !$omp private(vel) &
549  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J23G,MAPF)
550  do j = jjs, jje
551  do i = iis, iie
552  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS.
553  ! The flux at KS can be non-zero.
554  ! To reduce calculations, all the fluxes are set to zero.
555  flux(ks-1,i,j) = 0.0_rp ! k = KS
556 
557  ! physically incorrect but for numerical stability
558  vel = ( ( 0.5_rp * ( mom(ks+1,i,j)+mom(ks+1,i,j-1) ) ) / dens(ks+1,i,j) &
559  + ( 0.5_rp * ( mom(ks,i,j)+mom(ks,i,j-1) ) ) / dens(ks ,i,j) ) * 0.5_rp
560 ! vel = ( 0.5_RP * ( mom(KS+1,i,j)+mom(KS+1,i,j-1) ) ) &
561 ! / DENS(KS+1,i,j)
562  flux(ks,i,j) = j23g(ks+1,i,j) / mapf(i,j,+1) * vel &
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_ud3
572 
573 
574  !-----------------------------------------------------------------------------
576  subroutine atmos_dyn_fvm_fluxx_xyw_ud3( &
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  * ( ( f31 * ( val(k,i+2,j)+val(k,i-1,j) ) + f32 * ( val(k,i+1,j)+val(k,i,j) ) ) &
627  - ( f31 * ( val(k,i+2,j)-val(k,i-1,j) ) + f33 * ( val(k,i+1,j)-val(k,i,j) ) ) * sign(1.0_rp,vel) ) &
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_ud3
650 
651  !-----------------------------------------------------------------------------
653  subroutine atmos_dyn_fvm_fluxy_xyw_ud3( &
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  * ( ( f31 * ( val(k,i,j+2)+val(k,i,j-1) ) + f32 * ( val(k,i,j+1)+val(k,i,j) ) ) &
704  - ( f31 * ( val(k,i,j+2)-val(k,i,j-1) ) + f33 * ( val(k,i,j+1)-val(k,i,j) ) ) * sign(1.0_rp,vel) ) &
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_ud3
727 
728 
729  !-----------------------------------------------------------------------------
731  subroutine atmos_dyn_fvm_fluxz_uyz_ud3( &
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  * ( ( f31 * ( val(k+2,i,j)+val(k-1,i,j) ) + f32 * ( val(k+1,i,j)+val(k,i,j) ) ) &
779  - ( f31 * ( val(k+2,i,j)-val(k-1,i,j) ) + f33 * ( val(k+1,i,j)-val(k,i,j) ) ) * sign(1.0_rp,vel) ) &
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_ud3
832 
833  !-----------------------------------------------------------------------------
835  subroutine atmos_dyn_fvm_fluxj13_uyz_ud3( &
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  flux(k,i,j) = j13g(k,i,j) / mapf(i,j,+2) * vel &
873  * ( ( f31 * ( val(k+2,i,j)+val(k-1,i,j) ) + f32 * ( val(k+1,i,j)+val(k,i,j) ) ) &
874  - ( f31 * ( val(k+2,i,j)-val(k-1,i,j) ) + f33 * ( val(k+1,i,j)-val(k,i,j) ) ) * sign(1.0_rp,vel) )
875  enddo
876  enddo
877  enddo
878 
879  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
880  !$omp private(vel) &
881  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J13G,MAPF) &
882  !$omp shared(GSQRT,CDZ)
883  do j = jjs, jje
884  do i = iis, iie
885  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
886  ! The flux at KS-1 can be non-zero.
887  ! To reduce calculations, all the fluxes are set to zero.
888  flux(ks-1,i,j) = 0.0_rp
889 
890  vel = ( f2h(ks,1,i_uyz) &
891  * mom(ks+1,i,j) &
892  + f2h(ks,2,i_uyz) &
893  * mom(ks,i,j) ) &
894  / ( f2h(ks,1,i_uyz) &
895  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i+1,j) ) &
896  + f2h(ks,2,i_uyz) &
897  * 0.5_rp * ( dens(ks,i,j)+dens(ks,i+1,j) ) )
898  flux(ks,i,j) = j13g(ks,i,j) / mapf(i,j,+2) * vel &
899  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) )
900 
901  vel = ( f2h(ke-1,1,i_uyz) &
902  * mom(ke,i,j) &
903  + f2h(ke-1,2,i_uyz) &
904  * mom(ke-1,i,j) ) &
905  / ( f2h(ke-1,1,i_uyz) &
906  * 0.5_rp * ( dens(ke,i,j)+dens(ke,i+1,j) ) &
907  + f2h(ke-1,2,i_uyz) &
908  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i+1,j) ) )
909  flux(ke-1,i,j) = j13g(ke-1,i,j) / mapf(i,j,+2) * vel &
910  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) )
911 
912  flux(ke ,i,j) = 0.0_rp
913  enddo
914  enddo
915 
916  return
917  end subroutine atmos_dyn_fvm_fluxj13_uyz_ud3
918 
919  !-----------------------------------------------------------------------------
921  subroutine atmos_dyn_fvm_fluxj23_uyz_ud3( &
922  flux, &
923  mom, val, DENS, &
924  GSQRT, J23G, MAPF, &
925  CDZ, &
926  IIS, IIE, JJS, JJE )
927  implicit none
928 
929  real(RP), intent(inout) :: flux (ka,ia,ja)
930  real(RP), intent(in) :: mom (ka,ia,ja)
931  real(RP), intent(in) :: val (ka,ia,ja)
932  real(RP), intent(in) :: dens (ka,ia,ja)
933  real(RP), intent(in) :: gsqrt (ka,ia,ja)
934  real(RP), intent(in) :: j23g (ka,ia,ja)
935  real(RP), intent(in) :: mapf ( ia,ja,2)
936  real(RP), intent(in) :: cdz (ka)
937  integer, intent(in) :: iis, iie, jjs, jje
938 
939  real(RP) :: vel
940  integer :: k, i, j
941  !---------------------------------------------------------------------------
942 
943  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
944  !$omp private(vel) &
945  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J23G,MAPF) &
946  !$omp shared(GSQRT,CDZ)
947  do j = jjs, jje
948  do i = iis, iie
949  do k = ks+1, ke-2
950  vel = ( f2h(k,1,i_uyz) &
951  * 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) ) &
952  + f2h(k,2,i_uyz) &
953  * 0.25_rp * ( mom(k,i,j)+mom(k,i+1,j)+mom(k,i,j-1)+mom(k,i+1,j-1) ) ) &
954  / ( f2h(k,1,i_uyz) &
955  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i+1,j) ) &
956  + f2h(k,2,i_uyz) &
957  * 0.5_rp * ( dens(k,i,j)+dens(k,i+1,j) ) )
958  flux(k,i,j) = j23g(k,i,j) / mapf(i,j,+1) * vel &
959  * ( ( f31 * ( val(k+2,i,j)+val(k-1,i,j) ) + f32 * ( val(k+1,i,j)+val(k,i,j) ) ) &
960  - ( f31 * ( val(k+2,i,j)-val(k-1,i,j) ) + f33 * ( val(k+1,i,j)-val(k,i,j) ) ) * sign(1.0_rp,vel) )
961  enddo
962  enddo
963  enddo
964 
965  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
966  !$omp private(vel) &
967  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J23G,MAPF) &
968  !$omp shared(GSQRT,CDZ)
969  do j = jjs, jje
970  do i = iis, iie
971  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
972  ! The flux at KS-1 can be non-zero.
973  ! To reduce calculations, all the fluxes are set to zero.
974  flux(ks-1,i,j) = 0.0_rp
975 
976  vel = ( f2h(ks,1,i_uyz) &
977  * 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) ) &
978  + f2h(ks,2,i_uyz) &
979  * 0.25_rp * ( mom(ks,i,j)+mom(ks,i+1,j)+mom(ks,i,j-1)+mom(ks,i+1,j-1) ) ) &
980  / ( f2h(ks,1,i_uyz) &
981  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i+1,j) ) &
982  + f2h(ks,2,i_uyz) &
983  * 0.5_rp * ( dens(ks,i,j)+dens(ks,i+1,j) ) )
984  flux(ks,i,j) = j23g(ks,i,j) / mapf(i,j,+1) * vel &
985  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) )
986 
987  vel = ( f2h(ke-1,1,i_uyz) &
988  * 0.25_rp * ( mom(ke,i,j)+mom(ke,i+1,j)+mom(ke,i,j-1)+mom(ke,i+1,j-1) ) &
989  + f2h(ke-1,2,i_uyz) &
990  * 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) ) ) &
991  / ( f2h(ke-1,1,i_uyz) &
992  * 0.5_rp * ( dens(ke,i,j)+dens(ke,i+1,j) ) &
993  + f2h(ke-1,2,i_uyz) &
994  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i+1,j) ) )
995  flux(ke-1,i,j) = j23g(ke-1,i,j) / mapf(i,j,+1) * vel &
996  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) )
997 
998  flux(ke ,i,j) = 0.0_rp
999  enddo
1000  enddo
1001 
1002  return
1003  end subroutine atmos_dyn_fvm_fluxj23_uyz_ud3
1004 
1005  !-----------------------------------------------------------------------------
1007  subroutine atmos_dyn_fvm_fluxx_uyz_ud3( &
1008  flux, &
1009  mom, val, DENS, &
1010  GSQRT, MAPF, &
1011  num_diff, &
1012  CDZ, &
1013  IIS, IIE, JJS, JJE )
1014  implicit none
1015 
1016  real(RP), intent(inout) :: flux (ka,ia,ja)
1017  real(RP), intent(in) :: mom (ka,ia,ja)
1018  real(RP), intent(in) :: val (ka,ia,ja)
1019  real(RP), intent(in) :: dens (ka,ia,ja)
1020  real(RP), intent(in) :: gsqrt (ka,ia,ja)
1021  real(RP), intent(in) :: mapf ( ia,ja,2)
1022  real(RP), intent(in) :: num_diff(ka,ia,ja)
1023  real(RP), intent(in) :: cdz (ka)
1024  integer, intent(in) :: iis, iie, jjs, jje
1025 
1026  real(RP) :: vel
1027  integer :: k, i, j
1028  !---------------------------------------------------------------------------
1029 
1030  ! note that x-index is added by -1
1031 
1032  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
1033  !$omp private(vel) &
1034  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff)
1035  do j = jjs, jje
1036  do i = iis, iie+1
1037  do k = ks, ke
1038 #ifdef DEBUG
1039  call check( __line__, mom(k,i ,j) )
1040  call check( __line__, mom(k,i-1,j) )
1041 
1042  call check( __line__, val(k,i-1,j) )
1043  call check( __line__, val(k,i,j) )
1044 
1045  call check( __line__, val(k,i-2,j) )
1046  call check( __line__, val(k,i+1,j) )
1047 
1048 #endif
1049  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i-1,j) ) ) &
1050  / ( dens(k,i,j) )
1051  flux(k,i-1,j) = gsqrt(k,i,j) / mapf(i,j,+2) * vel &
1052  * ( ( f31 * ( val(k,i+1,j)+val(k,i-2,j) ) + f32 * ( val(k,i,j)+val(k,i-1,j) ) ) &
1053  - ( f31 * ( val(k,i+1,j)-val(k,i-2,j) ) + f33 * ( val(k,i,j)-val(k,i-1,j) ) ) * sign(1.0_rp,vel) ) &
1054  + gsqrt(k,i,j) * num_diff(k,i,j)
1055  enddo
1056  enddo
1057  enddo
1058 #ifdef DEBUG
1059  k = iundef; i = iundef; j = iundef
1060 #endif
1061 
1062  return
1063  end subroutine atmos_dyn_fvm_fluxx_uyz_ud3
1064 
1065  !-----------------------------------------------------------------------------
1067  subroutine atmos_dyn_fvm_fluxy_uyz_ud3( &
1068  flux, &
1069  mom, val, DENS, &
1070  GSQRT, MAPF, &
1071  num_diff, &
1072  CDZ, &
1073  IIS, IIE, JJS, JJE )
1074  implicit none
1075 
1076  real(RP), intent(inout) :: flux (ka,ia,ja)
1077  real(RP), intent(in) :: mom (ka,ia,ja)
1078  real(RP), intent(in) :: val (ka,ia,ja)
1079  real(RP), intent(in) :: dens (ka,ia,ja)
1080  real(RP), intent(in) :: gsqrt (ka,ia,ja)
1081  real(RP), intent(in) :: mapf ( ia,ja,2)
1082  real(RP), intent(in) :: num_diff(ka,ia,ja)
1083  real(RP), intent(in) :: cdz (ka)
1084  integer, intent(in) :: iis, iie, jjs, jje
1085 
1086  real(RP) :: vel
1087  integer :: k, i, j
1088  !---------------------------------------------------------------------------
1089 
1090  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
1091  !$omp private(vel) &
1092  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff)
1093  do j = jjs-1, jje
1094  do i = iis, iie
1095  do k = ks, ke
1096 #ifdef DEBUG
1097  call check( __line__, mom(k,i ,j) )
1098  call check( __line__, mom(k,i-1,j) )
1099 
1100  call check( __line__, val(k,i,j) )
1101  call check( __line__, val(k,i,j+1) )
1102 
1103  call check( __line__, val(k,i,j-1) )
1104  call check( __line__, val(k,i,j+2) )
1105 
1106 #endif
1107  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i+1,j) ) ) &
1108  / ( 0.25_rp * ( dens(k,i,j)+dens(k,i+1,j)+dens(k,i,j+1)+dens(k,i+1,j+1) ) )
1109  flux(k,i,j) = gsqrt(k,i,j) / mapf(i,j,+1) * vel &
1110  * ( ( f31 * ( val(k,i,j+2)+val(k,i,j-1) ) + f32 * ( val(k,i,j+1)+val(k,i,j) ) ) &
1111  - ( f31 * ( val(k,i,j+2)-val(k,i,j-1) ) + f33 * ( val(k,i,j+1)-val(k,i,j) ) ) * sign(1.0_rp,vel) ) &
1112  + gsqrt(k,i,j) * num_diff(k,i,j)
1113  enddo
1114  enddo
1115  enddo
1116 #ifdef DEBUG
1117  k = iundef; i = iundef; j = iundef
1118 #endif
1119 
1120  return
1121  end subroutine atmos_dyn_fvm_fluxy_uyz_ud3
1122 
1123 
1124 
1125  !-----------------------------------------------------------------------------
1127  subroutine atmos_dyn_fvm_fluxz_xvz_ud3( &
1128  flux, &
1129  mom, val, DENS, &
1130  GSQRT, J33G, &
1131  num_diff, &
1132  CDZ, &
1133  IIS, IIE, JJS, JJE )
1134  implicit none
1135 
1136  real(RP), intent(inout) :: flux (ka,ia,ja)
1137  real(RP), intent(in) :: mom (ka,ia,ja)
1138  real(RP), intent(in) :: val (ka,ia,ja)
1139  real(RP), intent(in) :: dens (ka,ia,ja)
1140  real(RP), intent(in) :: gsqrt (ka,ia,ja)
1141  real(RP), intent(in) :: j33g
1142  real(RP), intent(in) :: num_diff(ka,ia,ja)
1143  real(RP), intent(in) :: cdz (ka)
1144  integer, intent(in) :: iis, iie, jjs, jje
1145 
1146  real(RP) :: vel
1147  integer :: k, i, j
1148  !---------------------------------------------------------------------------
1149 
1150  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
1151  !$omp private(vel) &
1152  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J33G,GSQRT,num_diff) &
1153  !$omp shared(CDZ)
1154  do j = jjs, jje
1155  do i = iis, iie
1156  do k = ks+1, ke-2
1157 #ifdef DEBUG
1158  call check( __line__, mom(k,i,j) )
1159  call check( __line__, mom(k,i,j+1) )
1160 
1161  call check( __line__, val(k,i,j) )
1162  call check( __line__, val(k+1,i,j) )
1163 
1164  call check( __line__, val(k-1,i,j) )
1165  call check( __line__, val(k+2,i,j) )
1166 
1167 #endif
1168  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i,j+1) ) ) &
1169  / ( f2h(k,1,i_xvz) &
1170  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i,j+1) ) &
1171  + f2h(k,2,i_xvz) &
1172  * 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) ) )
1173  flux(k,i,j) = j33g * vel &
1174  * ( ( f31 * ( val(k+2,i,j)+val(k-1,i,j) ) + f32 * ( val(k+1,i,j)+val(k,i,j) ) ) &
1175  - ( f31 * ( val(k+2,i,j)-val(k-1,i,j) ) + f33 * ( val(k+1,i,j)-val(k,i,j) ) ) * sign(1.0_rp,vel) ) &
1176  + gsqrt(k,i,j) * num_diff(k,i,j)
1177  enddo
1178  enddo
1179  enddo
1180 #ifdef DEBUG
1181  k = iundef; i = iundef; j = iundef
1182 #endif
1183 
1184  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
1185  !$omp private(vel) &
1186  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,flux,J33G,GSQRT,num_diff,DENS,CDZ)
1187  do j = jjs, jje
1188  do i = iis, iie
1189 #ifdef DEBUG
1190 
1191  call check( __line__, mom(ks,i ,j) )
1192  call check( __line__, mom(ks,i,j+1) )
1193  call check( __line__, val(ks+1,i,j) )
1194  call check( __line__, val(ks,i,j) )
1195 
1196 #endif
1197  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
1198  ! The flux at KS-1 can be non-zero.
1199  ! To reduce calculations, all the fluxes are set to zero.
1200  flux(ks-1,i,j) = 0.0_rp
1201 
1202  vel = ( 0.5_rp * ( mom(ks,i,j)+mom(ks,i,j+1) ) ) &
1203  / ( f2h(ks,1,i_xvz) &
1204  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i,j+1) ) &
1205  + f2h(ks,2,i_xvz) &
1206  * 0.5_rp * ( dens(ks,i,j)+dens(ks,i,j+1) ) )
1207  flux(ks,i,j) = j33g * vel &
1208  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) ) &
1209  + gsqrt(ks,i,j) * num_diff(ks,i,j)
1210  vel = ( 0.5_rp * ( mom(ke-1,i,j)+mom(ke-1,i,j+1) ) ) &
1211  / ( f2h(ke-1,1,i_xvz) &
1212  * 0.5_rp * ( dens(ke,i,j)+dens(ke,i,j+1) ) &
1213  + f2h(ke-1,2,i_xvz) &
1214  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i,j+1) ) )
1215  flux(ke-1,i,j) = j33g * vel &
1216  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) ) &
1217  + gsqrt(ke-1,i,j) * num_diff(ke-1,i,j)
1218 
1219  flux(ke,i,j) = 0.0_rp
1220  enddo
1221  enddo
1222 #ifdef DEBUG
1223  k = iundef; i = iundef; j = iundef
1224 #endif
1225 
1226  return
1227  end subroutine atmos_dyn_fvm_fluxz_xvz_ud3
1228 
1229  !-----------------------------------------------------------------------------
1231  subroutine atmos_dyn_fvm_fluxj13_xvz_ud3( &
1232  flux, &
1233  mom, val, DENS, &
1234  GSQRT, J13G, MAPF, &
1235  CDZ, &
1236  IIS, IIE, JJS, JJE )
1237  implicit none
1238 
1239  real(RP), intent(inout) :: flux (ka,ia,ja)
1240  real(RP), intent(in) :: mom (ka,ia,ja)
1241  real(RP), intent(in) :: val (ka,ia,ja)
1242  real(RP), intent(in) :: dens (ka,ia,ja)
1243  real(RP), intent(in) :: gsqrt (ka,ia,ja)
1244  real(RP), intent(in) :: j13g (ka,ia,ja)
1245  real(RP), intent(in) :: mapf ( ia,ja,2)
1246  real(RP), intent(in) :: cdz (ka)
1247  integer, intent(in) :: iis, iie, jjs, jje
1248 
1249  real(RP) :: vel
1250  integer :: k, i, j
1251  !---------------------------------------------------------------------------
1252 
1253  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
1254  !$omp private(vel) &
1255  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J13G,MAPF) &
1256  !$omp shared(GSQRT,CDZ)
1257  do j = jjs, jje
1258  do i = iis, iie
1259  do k = ks+1, ke-2
1260  vel = ( f2h(k,1,i_xvz) &
1261  * 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) ) &
1262  + f2h(k,2,i_xvz) &
1263  * 0.25_rp * ( mom(k,i,j)+mom(k,i-1,j)+mom(k,i,j+1)+mom(k,i-1,j+1) ) ) &
1264  / ( f2h(k,1,i_xvz) &
1265  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i,j+1) ) &
1266  + f2h(k,2,i_xvz) &
1267  * 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) ) )
1268  flux(k,i,j) = j13g(k,i,j) / mapf(i,j,+2) * vel &
1269  * ( ( f31 * ( val(k+2,i,j)+val(k-1,i,j) ) + f32 * ( val(k+1,i,j)+val(k,i,j) ) ) &
1270  - ( f31 * ( val(k+2,i,j)-val(k-1,i,j) ) + f33 * ( val(k+1,i,j)-val(k,i,j) ) ) * sign(1.0_rp,vel) )
1271  enddo
1272  enddo
1273  enddo
1274 
1275  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
1276  !$omp private(vel) &
1277  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J13G,MAPF) &
1278  !$omp shared(GSQRT,CDZ)
1279  do j = jjs, jje
1280  do i = iis, iie
1281  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
1282  ! The flux at KS-1 can be non-zero.
1283  ! To reduce calculations, all the fluxes are set to zero.
1284  flux(ks-1,i,j) = 0.0_rp
1285 
1286  vel = ( f2h(ks,1,i_xvz) &
1287  * 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) ) &
1288  + f2h(ks,2,i_xvz) &
1289  * 0.25_rp * ( mom(ks,i,j)+mom(ks,i-1,j)+mom(ks,i,j+1)+mom(ks,i-1,j+1) ) ) &
1290  / ( f2h(ks,1,i_xvz) &
1291  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i,j+1) ) &
1292  + f2h(ks,2,i_xvz) &
1293  * 0.5_rp * ( dens(ks,i,j)+dens(ks,i,j+1) ) )
1294  flux(ks,i,j) = j13g(ks,i,j) / mapf(i,j,+2) * vel &
1295  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) )
1296 
1297  vel = ( f2h(ke-1,1,i_xvz) &
1298  * 0.25_rp * ( mom(ke,i,j)+mom(ke,i-1,j)+mom(ke,i,j+1)+mom(ke,i-1,j+1) ) &
1299  + f2h(ke-1,2,i_xvz) &
1300  * 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) ) ) &
1301  / ( f2h(ke-1,1,i_xvz) &
1302  * 0.5_rp * ( dens(ke,i,j)+dens(ke,i,j+1) ) &
1303  + f2h(ke-1,2,i_xvz) &
1304  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i,j+1) ) )
1305  flux(ke-1,i,j) = j13g(ke-1,i,j) / mapf(i,j,+2) * vel &
1306  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) )
1307 
1308  flux(ke ,i,j) = 0.0_rp
1309  enddo
1310  enddo
1311 
1312  return
1313  end subroutine atmos_dyn_fvm_fluxj13_xvz_ud3
1314 
1315  !-----------------------------------------------------------------------------
1317  subroutine atmos_dyn_fvm_fluxj23_xvz_ud3( &
1318  flux, &
1319  mom, val, DENS, &
1320  GSQRT, J23G, MAPF, &
1321  CDZ, &
1322  IIS, IIE, JJS, JJE )
1323  implicit none
1324 
1325  real(RP), intent(inout) :: flux (ka,ia,ja)
1326  real(RP), intent(in) :: mom (ka,ia,ja)
1327  real(RP), intent(in) :: val (ka,ia,ja)
1328  real(RP), intent(in) :: dens (ka,ia,ja)
1329  real(RP), intent(in) :: gsqrt (ka,ia,ja)
1330  real(RP), intent(in) :: j23g (ka,ia,ja)
1331  real(RP), intent(in) :: mapf ( ia,ja,2)
1332  real(RP), intent(in) :: cdz (ka)
1333  integer, intent(in) :: iis, iie, jjs, jje
1334 
1335  real(RP) :: vel
1336  integer :: k, i, j
1337  !---------------------------------------------------------------------------
1338 
1339  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
1340  !$omp private(vel) &
1341  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J23G,MAPF) &
1342  !$omp shared(GSQRT,CDZ)
1343  do j = jjs, jje
1344  do i = iis, iie
1345  do k = ks+1, ke-2
1346  vel = ( f2h(k,1,i_xvz) &
1347  * mom(k+1,i,j) &
1348  + f2h(k,2,i_xvz) &
1349  * mom(k,i,j) ) &
1350  / ( f2h(k,1,i_xvz) &
1351  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i,j+1) ) &
1352  + f2h(k,2,i_xvz) &
1353  * 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) ) )
1354  flux(k,i,j) = j23g(k,i,j) / mapf(i,j,+1) * vel &
1355  * ( ( f31 * ( val(k+2,i,j)+val(k-1,i,j) ) + f32 * ( val(k+1,i,j)+val(k,i,j) ) ) &
1356  - ( f31 * ( val(k+2,i,j)-val(k-1,i,j) ) + f33 * ( val(k+1,i,j)-val(k,i,j) ) ) * sign(1.0_rp,vel) )
1357  enddo
1358  enddo
1359  enddo
1360 
1361  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
1362  !$omp private(vel) &
1363  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J23G,MAPF) &
1364  !$omp shared(GSQRT,CDZ)
1365  do j = jjs, jje
1366  do i = iis, iie
1367  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
1368  ! The flux at KS-1 can be non-zero.
1369  ! To reduce calculations, all the fluxes are set to zero.
1370  flux(ks-1,i,j) = 0.0_rp
1371 
1372  vel = ( f2h(ks,1,i_xvz) &
1373  * mom(ks+1,i,j) &
1374  + f2h(ks,2,i_xvz) &
1375  * mom(ks,i,j) ) &
1376  / ( f2h(ks,1,i_xvz) &
1377  * 0.5_rp * ( dens(ks+1,i,j)+dens(ks+1,i,j+1) ) &
1378  + f2h(ks,2,i_xvz) &
1379  * 0.5_rp * ( dens(ks,i,j)+dens(ks,i,j+1) ) )
1380  flux(ks,i,j) = j23g(ks,i,j) / mapf(i,j,+1) * vel &
1381  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) )
1382 
1383  vel = ( f2h(ke-1,1,i_xvz) &
1384  * mom(ke,i,j) &
1385  + f2h(ke-1,2,i_xvz) &
1386  * mom(ke-1,i,j) ) &
1387  / ( f2h(ke-1,1,i_xvz) &
1388  * 0.5_rp * ( dens(ke,i,j)+dens(ke,i,j+1) ) &
1389  + f2h(ke-1,2,i_xvz) &
1390  * 0.5_rp * ( dens(ke-1,i,j)+dens(ke-1,i,j+1) ) )
1391  flux(ke-1,i,j) = j23g(ke-1,i,j) / mapf(i,j,+1) * vel &
1392  * ( f2 * ( val(ke,i,j)+val(ke-1,i,j) ) )
1393 
1394  flux(ke ,i,j) = 0.0_rp
1395  enddo
1396  enddo
1397 
1398  return
1399  end subroutine atmos_dyn_fvm_fluxj23_xvz_ud3
1400 
1401  !-----------------------------------------------------------------------------
1403  subroutine atmos_dyn_fvm_fluxx_xvz_ud3( &
1404  flux, &
1405  mom, val, DENS, &
1406  GSQRT, MAPF, &
1407  num_diff, &
1408  CDZ, &
1409  IIS, IIE, JJS, JJE )
1410  implicit none
1411 
1412  real(RP), intent(inout) :: flux (ka,ia,ja)
1413  real(RP), intent(in) :: mom (ka,ia,ja)
1414  real(RP), intent(in) :: val (ka,ia,ja)
1415  real(RP), intent(in) :: dens (ka,ia,ja)
1416  real(RP), intent(in) :: gsqrt (ka,ia,ja)
1417  real(RP), intent(in) :: mapf ( ia,ja,2)
1418  real(RP), intent(in) :: num_diff(ka,ia,ja)
1419  real(RP), intent(in) :: cdz (ka)
1420  integer, intent(in) :: iis, iie, jjs, jje
1421 
1422  real(RP) :: vel
1423  integer :: k, i, j
1424  !---------------------------------------------------------------------------
1425 
1426  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
1427  !$omp private(vel) &
1428  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff)
1429  do j = jjs, jje
1430  do i = iis-1, iie
1431  do k = ks, ke
1432 #ifdef DEBUG
1433  call check( __line__, mom(k,i ,j) )
1434  call check( __line__, mom(k,i,j-1) )
1435 
1436  call check( __line__, val(k,i,j) )
1437  call check( __line__, val(k,i+1,j) )
1438 
1439  call check( __line__, val(k,i-1,j) )
1440  call check( __line__, val(k,i+2,j) )
1441 
1442 #endif
1443  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i,j+1) ) ) &
1444  / ( 0.25_rp * ( dens(k,i,j)+dens(k,i+1,j)+dens(k,i,j+1)+dens(k,i+1,j+1) ) )
1445  flux(k,i,j) = gsqrt(k,i,j) / mapf(i,j,+2) * vel &
1446  * ( ( f31 * ( val(k,i+2,j)+val(k,i-1,j) ) + f32 * ( val(k,i+1,j)+val(k,i,j) ) ) &
1447  - ( f31 * ( val(k,i+2,j)-val(k,i-1,j) ) + f33 * ( val(k,i+1,j)-val(k,i,j) ) ) * sign(1.0_rp,vel) ) &
1448  + gsqrt(k,i,j) * num_diff(k,i,j)
1449  enddo
1450  enddo
1451  enddo
1452 #ifdef DEBUG
1453  k = iundef; i = iundef; j = iundef
1454 #endif
1455 
1456  return
1457  end subroutine atmos_dyn_fvm_fluxx_xvz_ud3
1458 
1459  !-----------------------------------------------------------------------------
1461  subroutine atmos_dyn_fvm_fluxy_xvz_ud3( &
1462  flux, &
1463  mom, val, DENS, &
1464  GSQRT, MAPF, &
1465  num_diff, &
1466  CDZ, &
1467  IIS, IIE, JJS, JJE )
1468  implicit none
1469 
1470  real(RP), intent(inout) :: flux (ka,ia,ja)
1471  real(RP), intent(in) :: mom (ka,ia,ja)
1472  real(RP), intent(in) :: val (ka,ia,ja)
1473  real(RP), intent(in) :: dens (ka,ia,ja)
1474  real(RP), intent(in) :: gsqrt (ka,ia,ja)
1475  real(RP), intent(in) :: mapf ( ia,ja,2)
1476  real(RP), intent(in) :: num_diff(ka,ia,ja)
1477  real(RP), intent(in) :: cdz (ka)
1478  integer, intent(in) :: iis, iie, jjs, jje
1479 
1480  real(RP) :: vel
1481  integer :: k, i, j
1482  !---------------------------------------------------------------------------
1483 
1484  ! note that y-index is added by -1
1485 
1486  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
1487  !$omp private(vel) &
1488  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff)
1489  do j = jjs, jje+1
1490  do i = iis, iie
1491  do k = ks, ke
1492 #ifdef DEBUG
1493  call check( __line__, mom(k,i ,j) )
1494  call check( __line__, mom(k,i,j-1) )
1495 
1496  call check( __line__, val(k,i,j-1) )
1497  call check( __line__, val(k,i,j) )
1498 
1499  call check( __line__, val(k,i,j-2) )
1500  call check( __line__, val(k,i,j+1) )
1501 
1502 #endif
1503  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i,j-1) ) ) &
1504  / ( dens(k,i,j) )
1505  flux(k,i,j-1) = gsqrt(k,i,j) / mapf(i,j,+1) * vel &
1506  * ( ( f31 * ( val(k,i,j+1)+val(k,i,j-2) ) + f32 * ( val(k,i,j)+val(k,i,j-1) ) ) &
1507  - ( f31 * ( val(k,i,j+1)-val(k,i,j-2) ) + f33 * ( val(k,i,j)-val(k,i,j-1) ) ) * sign(1.0_rp,vel) ) &
1508  + gsqrt(k,i,j) * num_diff(k,i,j)
1509  enddo
1510  enddo
1511  enddo
1512 #ifdef DEBUG
1513  k = iundef; i = iundef; j = iundef
1514 #endif
1515 
1516  return
1517  end subroutine atmos_dyn_fvm_fluxy_xvz_ud3
1518 
1519 
1520 
1521 
1522 
1523 
1524 
1526 
1527 !--
1528 ! vi:set readonly sw=4 ts=8
1529 !
1530 !Local Variables:
1531 !mode: f90
1532 !buffer-read-only: t
1533 !End:
1534 !
1535 !++
module DEBUG
Definition: scale_debug.F90:13
subroutine, public atmos_dyn_fvm_fluxj13_xvz_ud3(flux, mom, val, DENS, GSQRT, J13G, MAPF, CDZ, IIS, IIE, JJS, JJE)
calculation J13-flux at XVZ
subroutine, public atmos_dyn_fvm_fluxj23_uyz_ud3(flux, mom, val, DENS, GSQRT, J23G, MAPF, CDZ, IIS, IIE, JJS, JJE)
calculation J23-flux at UYZ
subroutine, public atmos_dyn_fvm_fluxz_xyz_ud3(flux, mflx, val, GSQRT, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation z-flux at XYZ
subroutine, public atmos_dyn_fvm_flux_valuew_z_ud3(valW, mflx, val, GSQRT, CDZ)
value at XYW
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_fluxx_xyw_ud3(flux, mom, val, DENS, GSQRT, MAPF, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation X-flux at XYW
real(rp), public const_undef
Definition: scale_const.F90:43
module grid index
module TRACER
module Index
Definition: scale_index.F90:14
integer, public ia
of whole cells: x, local, with HALO
subroutine, public atmos_dyn_fvm_fluxj23_xvz_ud3(flux, mom, val, DENS, GSQRT, J23G, MAPF, CDZ, IIS, IIE, JJS, JJE)
calculation J23-flux at XVZ
subroutine, public atmos_dyn_fvm_fluxx_xvz_ud3(flux, mom, val, DENS, GSQRT, MAPF, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation X-flux at XV
integer, public ka
of whole cells: z, local, with HALO
subroutine, public atmos_dyn_fvm_fluxz_xvz_ud3(flux, mom, val, DENS, GSQRT, J33G, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation z-flux at XV
subroutine, public atmos_dyn_fvm_fluxx_xyz_ud3(flux, mflx, val, GSQRT, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation X-flux at XYZ
subroutine, public atmos_dyn_fvm_fluxy_xyw_ud3(flux, mom, val, DENS, GSQRT, MAPF, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation Y-flux at XYW
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:40
module PROCESS
subroutine, public atmos_dyn_fvm_fluxy_xyz_ud3(flux, mflx, val, GSQRT, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation Y-flux at XYZ
subroutine, public atmos_dyn_fvm_fluxj23_xyw_ud3(flux, mom, val, DENS, GSQRT, J23G, MAPF, CDZ, IIS, IIE, JJS, JJE)
calculation J23-flux at XYW
module CONSTANT
Definition: scale_const.F90:14
subroutine, public atmos_dyn_fvm_fluxy_xvz_ud3(flux, mom, val, DENS, GSQRT, MAPF, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation Y-flux at XV
integer, public ks
start point of inner domain: z, local
subroutine, public atmos_dyn_fvm_fluxy_uyz_ud3(flux, mom, val, DENS, GSQRT, MAPF, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation Y-flux at UY
module profiler
Definition: scale_prof.F90:10
subroutine, public atmos_dyn_fvm_fluxz_uyz_ud3(flux, mom, val, DENS, GSQRT, J33G, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation z-flux at UY
subroutine, public atmos_dyn_fvm_fluxj13_xyw_ud3(flux, mom, val, DENS, GSQRT, J13G, MAPF, CDZ, IIS, IIE, JJS, JJE)
calculation J13-flux at XYW
subroutine, public atmos_dyn_fvm_fluxx_uyz_ud3(flux, mom, val, DENS, GSQRT, MAPF, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation X-flux at UY
module PRECISION
subroutine, public atmos_dyn_fvm_fluxj13_uyz_ud3(flux, mom, val, DENS, GSQRT, J13G, MAPF, CDZ, IIS, IIE, JJS, JJE)
calculation J13-flux at UYZ
subroutine, public atmos_dyn_fvm_fluxz_xyw_ud3(flux, mom, val, DENS, GSQRT, J33G, num_diff, CDZ, FDZ, dtrk, IIS, IIE, JJS, JJE)
calculation z-flux at XYW
module scale_atmos_dyn_fvm_flux_ud3
integer, public ja
of whole cells: y, local, with HALO