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