SCALE-RM
scale_atmos_dyn_fvm_flux_cd2.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 
88 
89 
90 contains
91  !-----------------------------------------------------------------------------
93 !OCL SERIAL
95  valW, &
96  mflx, val, GSQRT, &
97  CDZ )
98  implicit none
99 
100  real(RP), intent(out) :: valW (ka)
101  real(RP), intent(in) :: mflx (ka)
102  real(RP), intent(in) :: val (ka)
103  real(RP), intent(in) :: GSQRT(ka)
104  real(RP), intent(in) :: CDZ (ka)
105 
106  integer :: k
107  !---------------------------------------------------------------------------
108 
109  do k = ks, ke-1
110 #ifdef DEBUG
111  call check( __line__, mflx(k) )
112 
113  call check( __line__, val(k) )
114  call check( __line__, val(k+1) )
115 
116 #endif
117  valw(k) = f2 * ( val(k+1)+val(k) )
118  enddo
119 #ifdef DEBUG
120  k = iundef
121 #endif
122 
123 #ifdef DEBUG
124 
125 #endif
126 
127 
128  return
129  end subroutine atmos_dyn_fvm_flux_valuew_z_cd2
130 
131  !-----------------------------------------------------------------------------
133  subroutine atmos_dyn_fvm_fluxz_xyz_cd2( &
134  flux, &
135  mflx, val, GSQRT, &
136  num_diff, &
137  CDZ, &
138  IIS, IIE, JJS, JJE )
139  use scale_const, only: &
140  eps => const_eps
141  implicit none
142 
143  real(RP), intent(inout) :: flux (ka,ia,ja)
144  real(RP), intent(in) :: mflx (ka,ia,ja)
145  real(RP), intent(in) :: val (ka,ia,ja)
146  real(RP), intent(in) :: GSQRT (ka,ia,ja)
147  real(RP), intent(in) :: num_diff(ka,ia,ja)
148  real(RP), intent(in) :: CDZ (ka)
149  integer, intent(in) :: IIS, IIE, JJS, JJE
150 
151  real(RP) :: vel
152  integer :: k, i, j
153  !---------------------------------------------------------------------------
154 
155  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
156  !$omp private(vel) &
157  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mflx,val,flux,GSQRT,num_diff)
158  do j = jjs, jje
159  do i = iis, iie
160  do k = ks, ke-1
161 #ifdef DEBUG
162  call check( __line__, mflx(k,i,j) )
163 
164  call check( __line__, val(k,i,j) )
165  call check( __line__, val(k+1,i,j) )
166 
167 #endif
168  vel = mflx(k,i,j)
169  flux(k,i,j) = vel &
170  * ( f2 * ( val(k+1,i,j)+val(k,i,j) ) ) &
171  + gsqrt(k,i,j) * num_diff(k,i,j)
172  enddo
173  enddo
174  enddo
175 #ifdef DEBUG
176  k = iundef; i = iundef; j = iundef
177 #endif
178 
179  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
180  !$omp private(vel) &
181  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mflx,val,flux,GSQRT,num_diff,EPS)
182  do j = jjs, jje
183  do i = iis, iie
184 #ifdef DEBUG
185 
186 #endif
187  flux(ks-1,i,j) = 0.0_rp
188 
189  flux(ke ,i,j) = 0.0_rp
190  enddo
191  enddo
192 #ifdef DEBUG
193  k = iundef; i = iundef; j = iundef
194 #endif
195 
196  return
197  end subroutine atmos_dyn_fvm_fluxz_xyz_cd2
198 
199  !-----------------------------------------------------------------------------
201  subroutine atmos_dyn_fvm_fluxx_xyz_cd2( &
202  flux, &
203  mflx, val, GSQRT, &
204  num_diff, &
205  CDZ, &
206  IIS, IIE, JJS, JJE )
207  implicit none
208 
209  real(RP), intent(inout) :: flux (ka,ia,ja)
210  real(RP), intent(in) :: mflx (ka,ia,ja)
211  real(RP), intent(in) :: val (ka,ia,ja)
212  real(RP), intent(in) :: GSQRT (ka,ia,ja)
213  real(RP), intent(in) :: num_diff(ka,ia,ja)
214  real(RP), intent(in) :: CDZ(ka)
215  integer, intent(in) :: IIS, IIE, JJS, JJE
216 
217  real(RP) :: vel
218  integer :: k, i, j
219  !---------------------------------------------------------------------------
220 
221  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
222  !$omp private(vel) &
223  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mflx,val,flux,GSQRT,num_diff)
224  do j = jjs, jje
225  do i = iis-1, iie
226  do k = ks, ke
227 #ifdef DEBUG
228  call check( __line__, mflx(k,i,j) )
229 
230  call check( __line__, val(k,i,j) )
231  call check( __line__, val(k,i+1,j) )
232 
233 #endif
234  vel = mflx(k,i,j)
235  flux(k,i,j) = vel &
236  * ( f2 * ( val(k,i+1,j)+val(k,i,j) ) ) &
237  + gsqrt(k,i,j) * num_diff(k,i,j)
238  enddo
239  enddo
240  enddo
241 #ifdef DEBUG
242  k = iundef; i = iundef; j = iundef
243 #endif
244 
245  return
246  end subroutine atmos_dyn_fvm_fluxx_xyz_cd2
247 
248  !-----------------------------------------------------------------------------
250  subroutine atmos_dyn_fvm_fluxy_xyz_cd2( &
251  flux, &
252  mflx, val, GSQRT, &
253  num_diff, &
254  CDZ, &
255  IIS, IIE, JJS, JJE )
256  implicit none
257 
258  real(RP), intent(inout) :: flux (ka,ia,ja)
259  real(RP), intent(in) :: mflx (ka,ia,ja)
260  real(RP), intent(in) :: val (ka,ia,ja)
261  real(RP), intent(in) :: GSQRT (ka,ia,ja)
262  real(RP), intent(in) :: num_diff(ka,ia,ja)
263  real(RP), intent(in) :: CDZ(ka)
264  integer, intent(in) :: IIS, IIE, JJS, JJE
265 
266  real(RP) :: vel
267  integer :: k, i, j
268  !---------------------------------------------------------------------------
269 
270  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
271  !$omp private(vel) &
272  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mflx,val,flux,GSQRT,num_diff)
273  do j = jjs-1, jje
274  do i = iis, iie
275  do k = ks, ke
276 #ifdef DEBUG
277  call check( __line__, mflx(k,i,j) )
278 
279  call check( __line__, val(k,i,j) )
280  call check( __line__, val(k,i,j+1) )
281 
282 #endif
283  vel = mflx(k,i,j)
284  flux(k,i,j) = vel &
285  * ( f2 * ( val(k,i,j+1)+val(k,i,j) ) ) &
286  + gsqrt(k,i,j) * num_diff(k,i,j)
287  enddo
288  enddo
289  enddo
290 #ifdef DEBUG
291  k = iundef; i = iundef; j = iundef
292 #endif
293 
294  return
295  end subroutine atmos_dyn_fvm_fluxy_xyz_cd2
296 
297 
298  !-----------------------------------------------------------------------------
300  subroutine atmos_dyn_fvm_fluxz_xyw_cd2( &
301  flux, &
302  mom, val, DENS, &
303  GSQRT, J33G, &
304  num_diff, &
305  CDZ, FDZ, &
306  dtrk, &
307  IIS, IIE, JJS, JJE )
308  implicit none
309 
310  real(RP), intent(inout) :: flux (ka,ia,ja)
311  real(RP), intent(in) :: mom (ka,ia,ja)
312  real(RP), intent(in) :: val (ka,ia,ja)
313  real(RP), intent(in) :: DENS (ka,ia,ja)
314  real(RP), intent(in) :: GSQRT (ka,ia,ja)
315  real(RP), intent(in) :: J33G
316  real(RP), intent(in) :: num_diff(ka,ia,ja)
317  real(RP), intent(in) :: CDZ (ka)
318  real(RP), intent(in) :: FDZ (ka-1)
319  real(RP), intent(in) :: dtrk
320  integer, intent(in) :: IIS, IIE, JJS, JJE
321 
322  real(RP) :: vel
323  integer :: k, i, j
324  !---------------------------------------------------------------------------
325 
326  ! note than z-index is added by -1
327 
328  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
329  !$omp private(vel) &
330  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,flux,J33G,GSQRT,num_diff,DENS)
331  do j = jjs, jje
332  do i = iis, iie
333  do k = ks+1, ke-1
334 #ifdef DEBUG
335  call check( __line__, mom(k-1,i,j) )
336  call check( __line__, mom(k ,i,j) )
337 
338  call check( __line__, val(k-1,i,j) )
339  call check( __line__, val(k,i,j) )
340 
341 #endif
342  vel = ( 0.5_rp * ( mom(k-1,i,j) &
343  + mom(k,i,j) ) ) &
344  / dens(k,i,j)
345  flux(k-1,i,j) = j33g * vel &
346  * ( f2 * ( val(k,i,j)+val(k-1,i,j) ) ) &
347  + gsqrt(k,i,j) * num_diff(k,i,j)
348  enddo
349  enddo
350  enddo
351 #ifdef DEBUG
352  k = iundef; i = iundef; j = iundef
353 #endif
354 
355  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
356  !$omp private(vel) &
357  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J33G,GSQRT,num_diff,FDZ,dtrk)
358  do j = jjs, jje
359  do i = iis, iie
360 #ifdef DEBUG
361 
362 
363 #endif
364  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS.
365  ! The flux at KS can be non-zero.
366  ! To reduce calculations, all the fluxes are set to zero.
367  flux(ks-1,i,j) = 0.0_rp ! k = KS
368 
369 
370 
371  flux(ke-1,i,j) = 0.0_rp ! k = KE
372  flux(ke ,i,j) = 0.0_rp ! k = KE+1
373  enddo
374  enddo
375 
376  return
377  end subroutine atmos_dyn_fvm_fluxz_xyw_cd2
378 
379 
380  !-----------------------------------------------------------------------------
382  subroutine atmos_dyn_fvm_fluxj13_xyw_cd2( &
383  flux, &
384  mom, val, DENS, &
385  GSQRT, J13G, MAPF, &
386  CDZ, &
387  IIS, IIE, JJS, JJE )
388  implicit none
389 
390  real(RP), intent(inout) :: flux (ka,ia,ja)
391  real(RP), intent(in) :: mom (ka,ia,ja)
392  real(RP), intent(in) :: val (ka,ia,ja)
393  real(RP), intent(in) :: DENS (ka,ia,ja)
394  real(RP), intent(in) :: GSQRT (ka,ia,ja)
395  real(RP), intent(in) :: J13G (ka,ia,ja)
396  real(RP), intent(in) :: MAPF ( ia,ja,2)
397  real(RP), intent(in) :: CDZ (ka)
398  integer, intent(in) :: IIS, IIE, JJS, JJE
399 
400  real(RP) :: vel
401  integer :: k, i, j
402  !---------------------------------------------------------------------------
403 
404  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
405  !$omp private(vel) &
406  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J13G,MAPF)
407  do j = jjs, jje
408  do i = iis, iie
409  do k = ks+2, ke-1
410  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i-1,j) ) ) &
411  / dens(k,i,j)
412  vel = vel * j13g(k,i,j)
413  flux(k-1,i,j) = vel / mapf(i,j,+2) &
414  * ( f2 * ( val(k,i,j)+val(k-1,i,j) ) )
415  enddo
416  enddo
417  enddo
418 
419  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
420  !$omp private(vel) &
421  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J13G,MAPF)
422  do j = jjs, jje
423  do i = iis, iie
424  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS.
425  ! The flux at KS can be non-zero.
426  ! To reduce calculations, all the fluxes are set to zero.
427  flux(ks-1,i,j) = 0.0_rp ! k = KS
428 
429  ! physically incorrect but for numerical stability
430  vel = ( ( 0.5_rp * ( mom(ks+1,i,j)+mom(ks+1,i-1,j) ) ) / dens(ks+1,i,j) &
431  + ( 0.5_rp * ( mom(ks,i,j)+mom(ks,i-1,j) ) ) / dens(ks ,i,j) ) * 0.5_rp
432 ! vel = ( 0.5_RP * ( mom(KS+1,i,j)+mom(KS+1,i-1,j) ) ) &
433 ! / DENS(KS+1,i,j)
434  vel = vel * j13g(ks+1,i,j)
435  flux(ks,i,j) = vel / mapf(i,j,+2) &
436  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) ) ! k = KS+1
437 
438 
439  flux(ke-1,i,j) = 0.0_rp
440  enddo
441  enddo
442 
443  return
444  end subroutine atmos_dyn_fvm_fluxj13_xyw_cd2
445 
446  !-----------------------------------------------------------------------------
448  subroutine atmos_dyn_fvm_fluxj23_xyw_cd2( &
449  flux, &
450  mom, val, DENS, &
451  GSQRT, J23G, MAPF, &
452  CDZ, &
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) :: J23G (ka,ia,ja)
462  real(RP), intent(in) :: MAPF ( ia,ja,2)
463  real(RP), intent(in) :: CDZ (ka)
464  integer, intent(in) :: IIS, IIE, JJS, JJE
465 
466  real(RP) :: vel
467  integer :: k, i, j
468  !---------------------------------------------------------------------------
469 
470  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
471  !$omp private(vel) &
472  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J23G,MAPF)
473  do j = jjs, jje
474  do i = iis, iie
475  do k = ks+2, ke-1
476  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i,j-1) ) ) &
477  / dens(k,i,j)
478  vel = vel * j23g(k,i,j)
479  flux(k-1,i,j) = vel / mapf(i,j,+1) &
480  * ( f2 * ( val(k,i,j)+val(k-1,i,j) ) )
481  enddo
482  enddo
483  enddo
484 
485  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
486  !$omp private(vel) &
487  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J23G,MAPF)
488  do j = jjs, jje
489  do i = iis, iie
490  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS.
491  ! The flux at KS can be non-zero.
492  ! To reduce calculations, all the fluxes are set to zero.
493  flux(ks-1,i,j) = 0.0_rp ! k = KS
494 
495  ! physically incorrect but for numerical stability
496  vel = ( ( 0.5_rp * ( mom(ks+1,i,j)+mom(ks+1,i,j-1) ) ) / dens(ks+1,i,j) &
497  + ( 0.5_rp * ( mom(ks,i,j)+mom(ks,i,j-1) ) ) / dens(ks ,i,j) ) * 0.5_rp
498 ! vel = ( 0.5_RP * ( mom(KS+1,i,j)+mom(KS+1,i,j-1) ) ) &
499 ! / DENS(KS+1,i,j)
500  vel = vel * j23g(ks+1,i,j)
501  flux(ks,i,j) = vel / mapf(i,j,+1) &
502  * ( f2 * ( val(ks+1,i,j)+val(ks,i,j) ) ) ! k = KS+1
503 
504 
505  flux(ke-1,i,j) = 0.0_rp
506  enddo
507  enddo
508 
509  return
510  end subroutine atmos_dyn_fvm_fluxj23_xyw_cd2
511 
512 
513  !-----------------------------------------------------------------------------
515  subroutine atmos_dyn_fvm_fluxx_xyw_cd2( &
516  flux, &
517  mom, val, DENS, &
518  GSQRT, MAPF, &
519  num_diff, &
520  CDZ, &
521  IIS, IIE, JJS, JJE )
522  implicit none
523 
524  real(RP), intent(inout) :: flux (ka,ia,ja)
525  real(RP), intent(in) :: mom (ka,ia,ja)
526  real(RP), intent(in) :: val (ka,ia,ja)
527  real(RP), intent(in) :: DENS (ka,ia,ja)
528  real(RP), intent(in) :: GSQRT (ka,ia,ja)
529  real(RP), intent(in) :: MAPF ( ia,ja,2)
530  real(RP), intent(in) :: num_diff(ka,ia,ja)
531  real(RP), intent(in) :: CDZ (ka)
532  integer, intent(in) :: IIS, IIE, JJS, JJE
533 
534  real(RP) :: vel
535  integer :: k, i, j
536  !---------------------------------------------------------------------------
537 
538  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
539  !$omp private(vel) &
540  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff) &
541  !$omp shared(CDZ)
542  do j = jjs, jje
543  do i = iis-1, iie
544  do k = ks, ke-1
545 #ifdef DEBUG
546  call check( __line__, mom(k ,i,j) )
547  call check( __line__, mom(k+1,i,j) )
548 
549  call check( __line__, val(k,i,j) )
550  call check( __line__, val(k,i+1,j) )
551 
552 #endif
553  vel = ( f2h(k,1,i_uyz) &
554  * mom(k+1,i,j) &
555  + f2h(k,2,i_uyz) &
556  * mom(k,i,j) ) &
557  / ( f2h(k,1,i_uyz) &
558  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i+1,j) ) &
559  + f2h(k,2,i_uyz) &
560  * 0.5_rp * ( dens(k,i,j)+dens(k,i+1,j) ) )
561  flux(k,i,j) = gsqrt(k,i,j) / mapf(i,j,+2) * vel &
562  * ( f2 * ( val(k,i+1,j)+val(k,i,j) ) ) &
563  + gsqrt(k,i,j) * num_diff(k,i,j)
564  enddo
565  enddo
566  enddo
567 #ifdef DEBUG
568  k = iundef; i = iundef; j = iundef
569 #endif
570 
571  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
572  !$omp private(vel) &
573  !$omp shared(JJS,JJE,IIS,IIE,KE,flux)
574  do j = jjs, jje
575  do i = iis-1, iie
576  flux(ke,i,j) = 0.0_rp
577  enddo
578  enddo
579 #ifdef DEBUG
580  k = iundef; i = iundef; j = iundef
581 #endif
582 
583  return
584  end subroutine atmos_dyn_fvm_fluxx_xyw_cd2
585 
586  !-----------------------------------------------------------------------------
588  subroutine atmos_dyn_fvm_fluxy_xyw_cd2( &
589  flux, &
590  mom, val, DENS, &
591  GSQRT, MAPF, &
592  num_diff, &
593  CDZ, &
594  IIS, IIE, JJS, JJE )
595  implicit none
596 
597  real(RP), intent(inout) :: flux (ka,ia,ja)
598  real(RP), intent(in) :: mom (ka,ia,ja)
599  real(RP), intent(in) :: val (ka,ia,ja)
600  real(RP), intent(in) :: DENS (ka,ia,ja)
601  real(RP), intent(in) :: GSQRT (ka,ia,ja)
602  real(RP), intent(in) :: MAPF ( ia,ja,2)
603  real(RP), intent(in) :: num_diff(ka,ia,ja)
604  real(RP), intent(in) :: CDZ (ka)
605  integer, intent(in) :: IIS, IIE, JJS, JJE
606 
607  real(RP) :: vel
608  integer :: k, i, j
609  !---------------------------------------------------------------------------
610 
611  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
612  !$omp private(vel) &
613  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff) &
614  !$omp shared(CDZ)
615  do j = jjs-1, jje
616  do i = iis, iie
617  do k = ks, ke-1
618 #ifdef DEBUG
619  call check( __line__, mom(k ,i,j) )
620  call check( __line__, mom(k+1,i,j) )
621 
622  call check( __line__, val(k,i,j) )
623  call check( __line__, val(k,i,j+1) )
624 
625 #endif
626  vel = ( f2h(k,1,i_xvz) &
627  * mom(k+1,i,j) &
628  + f2h(k,2,i_xvz) &
629  * mom(k,i,j) ) &
630  / ( f2h(k,1,i_xvz) &
631  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i,j+1) ) &
632  + f2h(k,2,i_xvz) &
633  * 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) ) )
634  flux(k,i,j) = gsqrt(k,i,j) / mapf(i,j,+1) * vel &
635  * ( f2 * ( val(k,i,j+1)+val(k,i,j) ) ) &
636  + gsqrt(k,i,j) * num_diff(k,i,j)
637  enddo
638  enddo
639  enddo
640 #ifdef DEBUG
641  k = iundef; i = iundef; j = iundef
642 #endif
643 
644  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
645  !$omp private(vel) &
646  !$omp shared(JJS,JJE,IIS,IIE,KE,flux)
647  do j = jjs-1, jje
648  do i = iis, iie
649  flux(ke,i,j) = 0.0_rp
650  enddo
651  enddo
652 #ifdef DEBUG
653  k = iundef; i = iundef; j = iundef
654 #endif
655 
656  return
657  end subroutine atmos_dyn_fvm_fluxy_xyw_cd2
658 
659 
660  !-----------------------------------------------------------------------------
662  subroutine atmos_dyn_fvm_fluxz_uyz_cd2( &
663  flux, &
664  mom, val, DENS, &
665  GSQRT, J33G, &
666  num_diff, &
667  CDZ, &
668  IIS, IIE, JJS, JJE )
669  implicit none
670 
671  real(RP), intent(inout) :: flux (ka,ia,ja)
672  real(RP), intent(in) :: mom (ka,ia,ja)
673  real(RP), intent(in) :: val (ka,ia,ja)
674  real(RP), intent(in) :: DENS (ka,ia,ja)
675  real(RP), intent(in) :: GSQRT (ka,ia,ja)
676  real(RP), intent(in) :: J33G
677  real(RP), intent(in) :: num_diff(ka,ia,ja)
678  real(RP), intent(in) :: CDZ (ka)
679  integer, intent(in) :: IIS, IIE, JJS, JJE
680 
681  real(RP) :: vel
682  integer :: k, i, j
683  !---------------------------------------------------------------------------
684 
685  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
686  !$omp private(vel) &
687  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J33G,GSQRT,num_diff) &
688  !$omp shared(CDZ)
689  do j = jjs, jje
690  do i = iis, iie
691  do k = ks, ke-1
692 #ifdef DEBUG
693  call check( __line__, mom(k,i,j) )
694  call check( __line__, mom(k,i+1,j) )
695 
696  call check( __line__, val(k,i,j) )
697  call check( __line__, val(k+1,i,j) )
698 
699 #endif
700  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i+1,j) ) ) &
701  / ( f2h(k,1,i_uyz) &
702  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i+1,j) ) &
703  + f2h(k,2,i_uyz) &
704  * 0.5_rp * ( dens(k,i,j)+dens(k,i+1,j) ) )
705  flux(k,i,j) = j33g * vel &
706  * ( f2 * ( val(k+1,i,j)+val(k,i,j) ) ) &
707  + gsqrt(k,i,j) * num_diff(k,i,j)
708  enddo
709  enddo
710  enddo
711 #ifdef DEBUG
712  k = iundef; i = iundef; j = iundef
713 #endif
714 
715  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
716  !$omp private(vel) &
717  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,flux,J33G,GSQRT,num_diff,DENS,CDZ)
718  do j = jjs, jje
719  do i = iis, iie
720 #ifdef DEBUG
721 
722 #endif
723  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
724  ! The flux at KS-1 can be non-zero.
725  ! To reduce calculations, all the fluxes are set to zero.
726  flux(ks-1,i,j) = 0.0_rp
727 
728  flux(ke,i,j) = 0.0_rp
729  enddo
730  enddo
731 #ifdef DEBUG
732  k = iundef; i = iundef; j = iundef
733 #endif
734 
735  return
736  end subroutine atmos_dyn_fvm_fluxz_uyz_cd2
737 
738  !-----------------------------------------------------------------------------
740  subroutine atmos_dyn_fvm_fluxj13_uyz_cd2( &
741  flux, &
742  mom, val, DENS, &
743  GSQRT, J13G, MAPF, &
744  CDZ, &
745  IIS, IIE, JJS, JJE )
746  implicit none
747 
748  real(RP), intent(inout) :: flux (ka,ia,ja)
749  real(RP), intent(in) :: mom (ka,ia,ja)
750  real(RP), intent(in) :: val (ka,ia,ja)
751  real(RP), intent(in) :: DENS (ka,ia,ja)
752  real(RP), intent(in) :: GSQRT (ka,ia,ja)
753  real(RP), intent(in) :: J13G (ka,ia,ja)
754  real(RP), intent(in) :: MAPF ( ia,ja,2)
755  real(RP), intent(in) :: CDZ (ka)
756  integer, intent(in) :: IIS, IIE, JJS, JJE
757 
758  real(RP) :: vel
759  integer :: k, i, j
760  !---------------------------------------------------------------------------
761 
762  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
763  !$omp private(vel) &
764  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J13G,MAPF) &
765  !$omp shared(GSQRT,CDZ)
766  do j = jjs, jje
767  do i = iis, iie
768  do k = ks, ke-1
769  vel = ( f2h(k,1,i_uyz) &
770  * mom(k+1,i,j) &
771  + f2h(k,2,i_uyz) &
772  * mom(k,i,j) ) &
773  / ( f2h(k,1,i_uyz) &
774  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i+1,j) ) &
775  + f2h(k,2,i_uyz) &
776  * 0.5_rp * ( dens(k,i,j)+dens(k,i+1,j) ) )
777  vel = vel * j13g(k,i,j)
778  flux(k,i,j) = vel / mapf(i,j,+2) &
779  * ( f2 * ( val(k+1,i,j)+val(k,i,j) ) )
780  enddo
781  enddo
782  enddo
783 
784  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
785  !$omp private(vel) &
786  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J13G,MAPF) &
787  !$omp shared(GSQRT,CDZ)
788  do j = jjs, jje
789  do i = iis, iie
790  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
791  ! The flux at KS-1 can be non-zero.
792  ! To reduce calculations, all the fluxes are set to zero.
793  flux(ks-1,i,j) = 0.0_rp
794 
795  flux(ke ,i,j) = 0.0_rp
796  enddo
797  enddo
798 
799  return
800  end subroutine atmos_dyn_fvm_fluxj13_uyz_cd2
801 
802  !-----------------------------------------------------------------------------
804  subroutine atmos_dyn_fvm_fluxj23_uyz_cd2( &
805  flux, &
806  mom, val, DENS, &
807  GSQRT, J23G, MAPF, &
808  CDZ, &
809  IIS, IIE, JJS, JJE )
810  implicit none
811 
812  real(RP), intent(inout) :: flux (ka,ia,ja)
813  real(RP), intent(in) :: mom (ka,ia,ja)
814  real(RP), intent(in) :: val (ka,ia,ja)
815  real(RP), intent(in) :: DENS (ka,ia,ja)
816  real(RP), intent(in) :: GSQRT (ka,ia,ja)
817  real(RP), intent(in) :: J23G (ka,ia,ja)
818  real(RP), intent(in) :: MAPF ( ia,ja,2)
819  real(RP), intent(in) :: CDZ (ka)
820  integer, intent(in) :: IIS, IIE, JJS, JJE
821 
822  real(RP) :: vel
823  integer :: k, i, j
824  !---------------------------------------------------------------------------
825 
826  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
827  !$omp private(vel) &
828  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J23G,MAPF) &
829  !$omp shared(GSQRT,CDZ)
830  do j = jjs, jje
831  do i = iis, iie
832  do k = ks, ke-1
833  vel = ( f2h(k,1,i_uyz) &
834  * 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) ) &
835  + f2h(k,2,i_uyz) &
836  * 0.25_rp * ( mom(k,i,j)+mom(k,i+1,j)+mom(k,i,j-1)+mom(k,i+1,j-1) ) ) &
837  / ( f2h(k,1,i_uyz) &
838  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i+1,j) ) &
839  + f2h(k,2,i_uyz) &
840  * 0.5_rp * ( dens(k,i,j)+dens(k,i+1,j) ) )
841  vel = vel * j23g(k,i,j)
842  flux(k,i,j) = vel / mapf(i,j,+1) &
843  * ( f2 * ( val(k+1,i,j)+val(k,i,j) ) )
844  enddo
845  enddo
846  enddo
847 
848  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
849  !$omp private(vel) &
850  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J23G,MAPF) &
851  !$omp shared(GSQRT,CDZ)
852  do j = jjs, jje
853  do i = iis, iie
854  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
855  ! The flux at KS-1 can be non-zero.
856  ! To reduce calculations, all the fluxes are set to zero.
857  flux(ks-1,i,j) = 0.0_rp
858 
859  flux(ke ,i,j) = 0.0_rp
860  enddo
861  enddo
862 
863  return
864  end subroutine atmos_dyn_fvm_fluxj23_uyz_cd2
865 
866  !-----------------------------------------------------------------------------
868  subroutine atmos_dyn_fvm_fluxx_uyz_cd2( &
869  flux, &
870  mom, val, DENS, &
871  GSQRT, MAPF, &
872  num_diff, &
873  CDZ, &
874  IIS, IIE, JJS, JJE )
875  implicit none
876 
877  real(RP), intent(inout) :: flux (ka,ia,ja)
878  real(RP), intent(in) :: mom (ka,ia,ja)
879  real(RP), intent(in) :: val (ka,ia,ja)
880  real(RP), intent(in) :: DENS (ka,ia,ja)
881  real(RP), intent(in) :: GSQRT (ka,ia,ja)
882  real(RP), intent(in) :: MAPF ( ia,ja,2)
883  real(RP), intent(in) :: num_diff(ka,ia,ja)
884  real(RP), intent(in) :: CDZ (ka)
885  integer, intent(in) :: IIS, IIE, JJS, JJE
886 
887  real(RP) :: vel
888  integer :: k, i, j
889  !---------------------------------------------------------------------------
890 
891  ! note that x-index is added by -1
892 
893  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
894  !$omp private(vel) &
895  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff)
896  do j = jjs, jje
897  do i = iis, iie+1
898  do k = ks, ke
899 #ifdef DEBUG
900  call check( __line__, mom(k,i ,j) )
901  call check( __line__, mom(k,i-1,j) )
902 
903  call check( __line__, val(k,i-1,j) )
904  call check( __line__, val(k,i,j) )
905 
906 #endif
907  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i-1,j) ) ) &
908  / ( dens(k,i,j) )
909  flux(k,i-1,j) = gsqrt(k,i,j) / mapf(i,j,+2) * vel &
910  * ( f2 * ( val(k,i,j)+val(k,i-1,j) ) ) &
911  + gsqrt(k,i,j) * num_diff(k,i,j)
912  enddo
913  enddo
914  enddo
915 #ifdef DEBUG
916  k = iundef; i = iundef; j = iundef
917 #endif
918 
919  return
920  end subroutine atmos_dyn_fvm_fluxx_uyz_cd2
921 
922  !-----------------------------------------------------------------------------
924  subroutine atmos_dyn_fvm_fluxy_uyz_cd2( &
925  flux, &
926  mom, val, DENS, &
927  GSQRT, MAPF, &
928  num_diff, &
929  CDZ, &
930  IIS, IIE, JJS, JJE )
931  implicit none
932 
933  real(RP), intent(inout) :: flux (ka,ia,ja)
934  real(RP), intent(in) :: mom (ka,ia,ja)
935  real(RP), intent(in) :: val (ka,ia,ja)
936  real(RP), intent(in) :: DENS (ka,ia,ja)
937  real(RP), intent(in) :: GSQRT (ka,ia,ja)
938  real(RP), intent(in) :: MAPF ( ia,ja,2)
939  real(RP), intent(in) :: num_diff(ka,ia,ja)
940  real(RP), intent(in) :: CDZ (ka)
941  integer, intent(in) :: IIS, IIE, JJS, JJE
942 
943  real(RP) :: vel
944  integer :: k, i, j
945  !---------------------------------------------------------------------------
946 
947  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
948  !$omp private(vel) &
949  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff)
950  do j = jjs-1, jje
951  do i = iis, iie
952  do k = ks, ke
953 #ifdef DEBUG
954  call check( __line__, mom(k,i ,j) )
955  call check( __line__, mom(k,i-1,j) )
956 
957  call check( __line__, val(k,i,j) )
958  call check( __line__, val(k,i,j+1) )
959 
960 #endif
961  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i+1,j) ) ) &
962  / ( 0.25_rp * ( dens(k,i,j)+dens(k,i+1,j)+dens(k,i,j+1)+dens(k,i+1,j+1) ) )
963  flux(k,i,j) = gsqrt(k,i,j) / mapf(i,j,+1) * vel &
964  * ( f2 * ( val(k,i,j+1)+val(k,i,j) ) ) &
965  + gsqrt(k,i,j) * num_diff(k,i,j)
966  enddo
967  enddo
968  enddo
969 #ifdef DEBUG
970  k = iundef; i = iundef; j = iundef
971 #endif
972 
973  return
974  end subroutine atmos_dyn_fvm_fluxy_uyz_cd2
975 
976 
977 
978  !-----------------------------------------------------------------------------
980  subroutine atmos_dyn_fvm_fluxz_xvz_cd2( &
981  flux, &
982  mom, val, DENS, &
983  GSQRT, J33G, &
984  num_diff, &
985  CDZ, &
986  IIS, IIE, JJS, JJE )
987  implicit none
988 
989  real(RP), intent(inout) :: flux (ka,ia,ja)
990  real(RP), intent(in) :: mom (ka,ia,ja)
991  real(RP), intent(in) :: val (ka,ia,ja)
992  real(RP), intent(in) :: DENS (ka,ia,ja)
993  real(RP), intent(in) :: GSQRT (ka,ia,ja)
994  real(RP), intent(in) :: J33G
995  real(RP), intent(in) :: num_diff(ka,ia,ja)
996  real(RP), intent(in) :: CDZ (ka)
997  integer, intent(in) :: IIS, IIE, JJS, JJE
998 
999  real(RP) :: vel
1000  integer :: k, i, j
1001  !---------------------------------------------------------------------------
1002 
1003  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
1004  !$omp private(vel) &
1005  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J33G,GSQRT,num_diff) &
1006  !$omp shared(CDZ)
1007  do j = jjs, jje
1008  do i = iis, iie
1009  do k = ks, ke-1
1010 #ifdef DEBUG
1011  call check( __line__, mom(k,i,j) )
1012  call check( __line__, mom(k,i,j+1) )
1013 
1014  call check( __line__, val(k,i,j) )
1015  call check( __line__, val(k+1,i,j) )
1016 
1017 #endif
1018  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i,j+1) ) ) &
1019  / ( f2h(k,1,i_xvz) &
1020  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i,j+1) ) &
1021  + f2h(k,2,i_xvz) &
1022  * 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) ) )
1023  flux(k,i,j) = j33g * vel &
1024  * ( f2 * ( val(k+1,i,j)+val(k,i,j) ) ) &
1025  + gsqrt(k,i,j) * num_diff(k,i,j)
1026  enddo
1027  enddo
1028  enddo
1029 #ifdef DEBUG
1030  k = iundef; i = iundef; j = iundef
1031 #endif
1032 
1033  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
1034  !$omp private(vel) &
1035  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,flux,J33G,GSQRT,num_diff,DENS,CDZ)
1036  do j = jjs, jje
1037  do i = iis, iie
1038 #ifdef DEBUG
1039 
1040 #endif
1041  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
1042  ! The flux at KS-1 can be non-zero.
1043  ! To reduce calculations, all the fluxes are set to zero.
1044  flux(ks-1,i,j) = 0.0_rp
1045 
1046  flux(ke,i,j) = 0.0_rp
1047  enddo
1048  enddo
1049 #ifdef DEBUG
1050  k = iundef; i = iundef; j = iundef
1051 #endif
1052 
1053  return
1054  end subroutine atmos_dyn_fvm_fluxz_xvz_cd2
1055 
1056  !-----------------------------------------------------------------------------
1058  subroutine atmos_dyn_fvm_fluxj13_xvz_cd2( &
1059  flux, &
1060  mom, val, DENS, &
1061  GSQRT, J13G, MAPF, &
1062  CDZ, &
1063  IIS, IIE, JJS, JJE )
1064  implicit none
1065 
1066  real(RP), intent(inout) :: flux (ka,ia,ja)
1067  real(RP), intent(in) :: mom (ka,ia,ja)
1068  real(RP), intent(in) :: val (ka,ia,ja)
1069  real(RP), intent(in) :: DENS (ka,ia,ja)
1070  real(RP), intent(in) :: GSQRT (ka,ia,ja)
1071  real(RP), intent(in) :: J13G (ka,ia,ja)
1072  real(RP), intent(in) :: MAPF ( ia,ja,2)
1073  real(RP), intent(in) :: CDZ (ka)
1074  integer, intent(in) :: IIS, IIE, JJS, JJE
1075 
1076  real(RP) :: vel
1077  integer :: k, i, j
1078  !---------------------------------------------------------------------------
1079 
1080  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
1081  !$omp private(vel) &
1082  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J13G,MAPF) &
1083  !$omp shared(GSQRT,CDZ)
1084  do j = jjs, jje
1085  do i = iis, iie
1086  do k = ks, ke-1
1087  vel = ( f2h(k,1,i_xvz) &
1088  * 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) ) &
1089  + f2h(k,2,i_xvz) &
1090  * 0.25_rp * ( mom(k,i,j)+mom(k,i-1,j)+mom(k,i,j+1)+mom(k,i-1,j+1) ) ) &
1091  / ( f2h(k,1,i_xvz) &
1092  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i,j+1) ) &
1093  + f2h(k,2,i_xvz) &
1094  * 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) ) )
1095  vel = vel * j13g(k,i,j)
1096  flux(k,i,j) = vel / mapf(i,j,+2) &
1097  * ( f2 * ( val(k+1,i,j)+val(k,i,j) ) )
1098  enddo
1099  enddo
1100  enddo
1101 
1102  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
1103  !$omp private(vel) &
1104  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J13G,MAPF) &
1105  !$omp shared(GSQRT,CDZ)
1106  do j = jjs, jje
1107  do i = iis, iie
1108  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
1109  ! The flux at KS-1 can be non-zero.
1110  ! To reduce calculations, all the fluxes are set to zero.
1111  flux(ks-1,i,j) = 0.0_rp
1112 
1113  flux(ke ,i,j) = 0.0_rp
1114  enddo
1115  enddo
1116 
1117  return
1118  end subroutine atmos_dyn_fvm_fluxj13_xvz_cd2
1119 
1120  !-----------------------------------------------------------------------------
1122  subroutine atmos_dyn_fvm_fluxj23_xvz_cd2( &
1123  flux, &
1124  mom, val, DENS, &
1125  GSQRT, J23G, MAPF, &
1126  CDZ, &
1127  IIS, IIE, JJS, JJE )
1128  implicit none
1129 
1130  real(RP), intent(inout) :: 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) :: J23G (ka,ia,ja)
1136  real(RP), intent(in) :: MAPF ( ia,ja,2)
1137  real(RP), intent(in) :: CDZ (ka)
1138  integer, intent(in) :: IIS, IIE, JJS, JJE
1139 
1140  real(RP) :: vel
1141  integer :: k, i, j
1142  !---------------------------------------------------------------------------
1143 
1144  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
1145  !$omp private(vel) &
1146  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J23G,MAPF) &
1147  !$omp shared(GSQRT,CDZ)
1148  do j = jjs, jje
1149  do i = iis, iie
1150  do k = ks, ke-1
1151  vel = ( f2h(k,1,i_xvz) &
1152  * mom(k+1,i,j) &
1153  + f2h(k,2,i_xvz) &
1154  * mom(k,i,j) ) &
1155  / ( f2h(k,1,i_xvz) &
1156  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i,j+1) ) &
1157  + f2h(k,2,i_xvz) &
1158  * 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) ) )
1159  vel = vel * j23g(k,i,j)
1160  flux(k,i,j) = vel / mapf(i,j,+1) &
1161  * ( f2 * ( val(k+1,i,j)+val(k,i,j) ) )
1162  enddo
1163  enddo
1164  enddo
1165 
1166  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
1167  !$omp private(vel) &
1168  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,J23G,MAPF) &
1169  !$omp shared(GSQRT,CDZ)
1170  do j = jjs, jje
1171  do i = iis, iie
1172  ! The boundary condition is qflx_hi + qflxJ13 + qfluxJ23 = 0 at KS-1.
1173  ! The flux at KS-1 can be non-zero.
1174  ! To reduce calculations, all the fluxes are set to zero.
1175  flux(ks-1,i,j) = 0.0_rp
1176 
1177  flux(ke ,i,j) = 0.0_rp
1178  enddo
1179  enddo
1180 
1181  return
1182  end subroutine atmos_dyn_fvm_fluxj23_xvz_cd2
1183 
1184  !-----------------------------------------------------------------------------
1186  subroutine atmos_dyn_fvm_fluxx_xvz_cd2( &
1187  flux, &
1188  mom, val, DENS, &
1189  GSQRT, MAPF, &
1190  num_diff, &
1191  CDZ, &
1192  IIS, IIE, JJS, JJE )
1193  implicit none
1194 
1195  real(RP), intent(inout) :: flux (ka,ia,ja)
1196  real(RP), intent(in) :: mom (ka,ia,ja)
1197  real(RP), intent(in) :: val (ka,ia,ja)
1198  real(RP), intent(in) :: DENS (ka,ia,ja)
1199  real(RP), intent(in) :: GSQRT (ka,ia,ja)
1200  real(RP), intent(in) :: MAPF ( ia,ja,2)
1201  real(RP), intent(in) :: num_diff(ka,ia,ja)
1202  real(RP), intent(in) :: CDZ (ka)
1203  integer, intent(in) :: IIS, IIE, JJS, JJE
1204 
1205  real(RP) :: vel
1206  integer :: k, i, j
1207  !---------------------------------------------------------------------------
1208 
1209  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
1210  !$omp private(vel) &
1211  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff)
1212  do j = jjs, jje
1213  do i = iis-1, iie
1214  do k = ks, ke
1215 #ifdef DEBUG
1216  call check( __line__, mom(k,i ,j) )
1217  call check( __line__, mom(k,i,j-1) )
1218 
1219  call check( __line__, val(k,i,j) )
1220  call check( __line__, val(k,i+1,j) )
1221 
1222 #endif
1223  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i,j+1) ) ) &
1224  / ( 0.25_rp * ( dens(k,i,j)+dens(k,i+1,j)+dens(k,i,j+1)+dens(k,i+1,j+1) ) )
1225  flux(k,i,j) = gsqrt(k,i,j) / mapf(i,j,+2) * vel &
1226  * ( f2 * ( val(k,i+1,j)+val(k,i,j) ) ) &
1227  + gsqrt(k,i,j) * num_diff(k,i,j)
1228  enddo
1229  enddo
1230  enddo
1231 #ifdef DEBUG
1232  k = iundef; i = iundef; j = iundef
1233 #endif
1234 
1235  return
1236  end subroutine atmos_dyn_fvm_fluxx_xvz_cd2
1237 
1238  !-----------------------------------------------------------------------------
1240  subroutine atmos_dyn_fvm_fluxy_xvz_cd2( &
1241  flux, &
1242  mom, val, DENS, &
1243  GSQRT, MAPF, &
1244  num_diff, &
1245  CDZ, &
1246  IIS, IIE, JJS, JJE )
1247  implicit none
1248 
1249  real(RP), intent(inout) :: flux (ka,ia,ja)
1250  real(RP), intent(in) :: mom (ka,ia,ja)
1251  real(RP), intent(in) :: val (ka,ia,ja)
1252  real(RP), intent(in) :: DENS (ka,ia,ja)
1253  real(RP), intent(in) :: GSQRT (ka,ia,ja)
1254  real(RP), intent(in) :: MAPF ( ia,ja,2)
1255  real(RP), intent(in) :: num_diff(ka,ia,ja)
1256  real(RP), intent(in) :: CDZ (ka)
1257  integer, intent(in) :: IIS, IIE, JJS, JJE
1258 
1259  real(RP) :: vel
1260  integer :: k, i, j
1261  !---------------------------------------------------------------------------
1262 
1263  ! note that y-index is added by -1
1264 
1265  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
1266  !$omp private(vel) &
1267  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,mom,val,DENS,flux,GSQRT,MAPF,num_diff)
1268  do j = jjs, jje+1
1269  do i = iis, iie
1270  do k = ks, ke
1271 #ifdef DEBUG
1272  call check( __line__, mom(k,i ,j) )
1273  call check( __line__, mom(k,i,j-1) )
1274 
1275  call check( __line__, val(k,i,j-1) )
1276  call check( __line__, val(k,i,j) )
1277 
1278 #endif
1279  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i,j-1) ) ) &
1280  / ( dens(k,i,j) )
1281  flux(k,i,j-1) = gsqrt(k,i,j) / mapf(i,j,+1) * vel &
1282  * ( f2 * ( val(k,i,j)+val(k,i,j-1) ) ) &
1283  + gsqrt(k,i,j) * num_diff(k,i,j)
1284  enddo
1285  enddo
1286  enddo
1287 #ifdef DEBUG
1288  k = iundef; i = iundef; j = iundef
1289 #endif
1290 
1291  return
1292  end subroutine atmos_dyn_fvm_fluxy_xvz_cd2
1293 
1294 
1295 
1296 
1297 
1298 
1299 
1301 
1302 !--
1303 ! vi:set readonly sw=4 ts=8
1304 !
1305 !Local Variables:
1306 !mode: f90
1307 !buffer-read-only: t
1308 !End:
1309 !
1310 !++
module DEBUG
Definition: scale_debug.F90:11
subroutine, public atmos_dyn_fvm_fluxz_xyz_cd2(flux, mflx, val, GSQRT, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation z-flux at XYZ
subroutine, public atmos_dyn_fvm_fluxz_xyw_cd2(flux, mom, val, DENS, GSQRT, J33G, num_diff, CDZ, FDZ, dtrk, IIS, IIE, JJS, JJE)
calculation z-flux at XYW
module scale_atmos_dyn_fvm_flux_cd2
integer, public ia
of whole cells: x, local, with HALO
subroutine, public atmos_dyn_fvm_fluxy_xyz_cd2(flux, mflx, val, GSQRT, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation Y-flux at XYZ
subroutine, public atmos_dyn_fvm_fluxj23_xyw_cd2(flux, mom, val, DENS, GSQRT, J23G, MAPF, CDZ, IIS, IIE, JJS, JJE)
calculation J23-flux at XYW
subroutine, public atmos_dyn_fvm_fluxx_uyz_cd2(flux, mom, val, DENS, GSQRT, MAPF, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation X-flux at UY
integer, public ja
of whole cells: y, local, with HALO
subroutine, public atmos_dyn_fvm_fluxy_uyz_cd2(flux, mom, val, DENS, GSQRT, MAPF, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation Y-flux at UY
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:56
subroutine, public atmos_dyn_fvm_fluxj13_xvz_cd2(flux, mom, val, DENS, GSQRT, J13G, MAPF, CDZ, IIS, IIE, JJS, JJE)
calculation J13-flux at XVZ
subroutine, public atmos_dyn_fvm_fluxz_xvz_cd2(flux, mom, val, DENS, GSQRT, J33G, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation z-flux at XV
subroutine, public atmos_dyn_fvm_fluxz_uyz_cd2(flux, mom, val, DENS, GSQRT, J33G, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation z-flux at UY
subroutine, public atmos_dyn_fvm_fluxy_xyw_cd2(flux, mom, val, DENS, GSQRT, MAPF, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation Y-flux at XYW
real(rp), public const_undef
Definition: scale_const.F90:41
subroutine, public atmos_dyn_fvm_fluxx_xvz_cd2(flux, mom, val, DENS, GSQRT, MAPF, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation X-flux at XV
module TRACER
module Index
Definition: scale_index.F90:11
module atmosphere / grid / cartesC index
integer, public ke
end point of inner domain: z, local
subroutine, public atmos_dyn_fvm_fluxj13_uyz_cd2(flux, mom, val, DENS, GSQRT, J13G, MAPF, CDZ, IIS, IIE, JJS, JJE)
calculation J13-flux at UYZ
module PROCESS
Definition: scale_prc.F90:11
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:38
integer, public ks
start point of inner domain: z, local
subroutine, public atmos_dyn_fvm_fluxj23_xvz_cd2(flux, mom, val, DENS, GSQRT, J23G, MAPF, CDZ, IIS, IIE, JJS, JJE)
calculation J23-flux at XVZ
module CONSTANT
Definition: scale_const.F90:11
module profiler
Definition: scale_prof.F90:11
real(rp), public const_eps
small number
Definition: scale_const.F90:33
module PRECISION
subroutine, public atmos_dyn_fvm_fluxj13_xyw_cd2(flux, mom, val, DENS, GSQRT, J13G, MAPF, CDZ, IIS, IIE, JJS, JJE)
calculation J13-flux at XYW
integer, public ka
of whole cells: z, local, with HALO
subroutine, public atmos_dyn_fvm_fluxy_xvz_cd2(flux, mom, val, DENS, GSQRT, MAPF, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation Y-flux at XV
module STDIO
Definition: scale_io.F90:10
subroutine, public atmos_dyn_fvm_fluxx_xyw_cd2(flux, mom, val, DENS, GSQRT, MAPF, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation X-flux at XYW
subroutine, public atmos_dyn_fvm_fluxx_xyz_cd2(flux, mflx, val, GSQRT, num_diff, CDZ, IIS, IIE, JJS, JJE)
calculation X-flux at XYZ
subroutine, public atmos_dyn_fvm_flux_valuew_z_cd2(valW, mflx, val, GSQRT, CDZ)
value at XYW
subroutine, public atmos_dyn_fvm_fluxj23_uyz_cd2(flux, mom, val, DENS, GSQRT, J23G, MAPF, CDZ, IIS, IIE, JJS, JJE)
calculation J23-flux at UYZ