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