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