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