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