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