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