SCALE-RM
scale_atmos_dyn_fvm_flux_ud1.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 :: f1 = 0.5_rp
87 
88 
89 
90 
91 
92 
93 contains
94  !-----------------------------------------------------------------------------
96 !OCL SERIAL
98  valW, &
99  mflx, val, GSQRT, &
100  CDZ )
101  implicit none
102 
103  real(RP), intent(out) :: valW (ka)
104  real(RP), intent(in) :: mflx (ka)
105  real(RP), intent(in) :: val (ka)
106  real(RP), intent(in) :: GSQRT(ka)
107  real(RP), intent(in) :: CDZ (ka)
108 
109  integer :: k
110  !---------------------------------------------------------------------------
111 
112  do k = ks, ke-1
113 #ifdef DEBUG
114  call check( __line__, mflx(k) )
115 
116  call check( __line__, val(k) )
117  call check( __line__, val(k+1) )
118 
119 #endif
120  valw(k) = f1 * ( val(k+1)+val(k) ) - sign(f1,mflx(k)) * ( val(k+1)-val(k) )
121  enddo
122 #ifdef DEBUG
123  k = iundef
124 #endif
125 
126 #ifdef DEBUG
127 
128 #endif
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 
139  CDZ, &
140  IIS, IIE, JJS, JJE )
141  implicit none
142 
143  real(RP), intent(out) :: flux (ka,ia,ja)
144  real(RP), intent(in) :: mflx (ka,ia,ja)
145  real(RP), intent(in) :: val (ka,ia,ja)
146  real(RP), intent(in) :: GSQRT (ka,ia,ja)
147 
148  real(RP), intent(in) :: CDZ (ka)
149  integer, intent(in) :: IIS, IIE, JJS, JJE
150 
151  real(RP) :: vel
152  integer :: k, i, j
153  !---------------------------------------------------------------------------
154 
155  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
156  do j = jjs, jje
157  do i = iis, iie
158  do k = ks, ke-1
159 #ifdef DEBUG
160  call check( __line__, mflx(k,i,j) )
161 
162  call check( __line__, val(k,i,j) )
163  call check( __line__, val(k+1,i,j) )
164 
165 #endif
166  vel = mflx(k,i,j)
167  flux(k,i,j) = vel &
168  * ( f1 * ( val(k+1,i,j)+val(k,i,j) ) - sign(f1,vel) * ( val(k+1,i,j)-val(k,i,j) ) )
169  enddo
170  enddo
171  enddo
172 #ifdef DEBUG
173  k = iundef; i = iundef; j = iundef
174 #endif
175 
176  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
177  do j = jjs, jje
178  do i = iis, iie
179 #ifdef DEBUG
180 
181 #endif
182  flux(ks-1,i,j) = 0.0_rp
183 
184  flux(ke ,i,j) = 0.0_rp
185  enddo
186  enddo
187 #ifdef DEBUG
188  k = iundef; i = iundef; j = iundef
189 #endif
190 
191  return
192  end subroutine atmos_dyn_fvm_fluxz_xyz_ud1
193 
194  !-----------------------------------------------------------------------------
196  subroutine atmos_dyn_fvm_fluxx_xyz_ud1( &
197  flux, &
198  mflx, val, GSQRT, &
199 
200  CDZ, &
201  IIS, IIE, JJS, JJE )
202  implicit none
203 
204  real(RP), intent(out) :: flux (ka,ia,ja)
205  real(RP), intent(in) :: mflx (ka,ia,ja)
206  real(RP), intent(in) :: val (ka,ia,ja)
207  real(RP), intent(in) :: GSQRT (ka,ia,ja)
208 
209  real(RP), intent(in) :: CDZ(ka)
210  integer, intent(in) :: IIS, IIE, JJS, JJE
211 
212  real(RP) :: vel
213  integer :: k, i, j
214  !---------------------------------------------------------------------------
215 
216  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
217  do j = jjs, jje
218  do i = iis-1, iie
219  do k = ks, ke
220 #ifdef DEBUG
221  call check( __line__, mflx(k,i,j) )
222 
223  call check( __line__, val(k,i,j) )
224  call check( __line__, val(k,i+1,j) )
225 
226 #endif
227  vel = mflx(k,i,j)
228  flux(k,i,j) = vel &
229  * ( f1 * ( val(k,i+1,j)+val(k,i,j) ) - sign(f1,vel) * ( val(k,i+1,j)-val(k,i,j) ) )
230  enddo
231  enddo
232  enddo
233 #ifdef DEBUG
234  k = iundef; i = iundef; j = iundef
235 #endif
236 
237  return
238  end subroutine atmos_dyn_fvm_fluxx_xyz_ud1
239 
240  !-----------------------------------------------------------------------------
242  subroutine atmos_dyn_fvm_fluxy_xyz_ud1( &
243  flux, &
244  mflx, val, GSQRT, &
245 
246  CDZ, &
247  IIS, IIE, JJS, JJE )
248  implicit none
249 
250  real(RP), intent(out) :: flux (ka,ia,ja)
251  real(RP), intent(in) :: mflx (ka,ia,ja)
252  real(RP), intent(in) :: val (ka,ia,ja)
253  real(RP), intent(in) :: GSQRT (ka,ia,ja)
254 
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-1, jje
264  do i = iis, 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,j+1) )
271 
272 #endif
273  vel = mflx(k,i,j)
274  flux(k,i,j) = vel &
275  * ( f1 * ( val(k,i,j+1)+val(k,i,j) ) - sign(f1,vel) * ( val(k,i,j+1)-val(k,i,j) ) )
276  enddo
277  enddo
278  enddo
279 #ifdef DEBUG
280  k = iundef; i = iundef; j = iundef
281 #endif
282 
283  return
284  end subroutine atmos_dyn_fvm_fluxy_xyz_ud1
285 
286 
287  !-----------------------------------------------------------------------------
289  subroutine atmos_dyn_fvm_fluxz_xyw_ud1( &
290  flux, &
291  mom, val, DENS, &
292  GSQRT, J33G, &
293 
294  CDZ, FDZ, &
295  dtrk, &
296  IIS, IIE, JJS, JJE )
297  implicit none
298 
299  real(RP), intent(out) :: flux (ka,ia,ja)
300  real(RP), intent(in) :: mom (ka,ia,ja)
301  real(RP), intent(in) :: val (ka,ia,ja)
302  real(RP), intent(in) :: DENS (ka,ia,ja)
303  real(RP), intent(in) :: GSQRT (ka,ia,ja)
304  real(RP), intent(in) :: J33G
305 
306  real(RP), intent(in) :: CDZ (ka)
307  real(RP), intent(in) :: FDZ (ka-1)
308  real(RP), intent(in) :: dtrk
309  integer, intent(in) :: IIS, IIE, JJS, JJE
310 
311  real(RP) :: vel
312  real(RP) :: sw
313  integer :: k, i, j
314  !---------------------------------------------------------------------------
315 
316  ! note than z-index is added by -1
317 
318  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
319  do j = jjs, jje
320  do i = iis, iie
321  do k = ks+1, ke-1
322 #ifdef DEBUG
323  call check( __line__, mom(k-1,i,j) )
324  call check( __line__, mom(k ,i,j) )
325 
326  call check( __line__, val(k-1,i,j) )
327  call check( __line__, val(k,i,j) )
328 
329 #endif
330  vel = ( 0.5_rp * ( mom(k-1,i,j) &
331  + mom(k,i,j) ) ) &
332  / dens(k,i,j)
333  flux(k-1,i,j) = j33g * vel &
334  * ( f1 * ( val(k,i,j)+val(k-1,i,j) ) - sign(f1,vel) * ( val(k,i,j)-val(k-1,i,j) ) )
335  enddo
336  enddo
337  enddo
338 #ifdef DEBUG
339  k = iundef; i = iundef; j = iundef
340 #endif
341 
342  !$omp parallel do private(i,j,sw) OMP_SCHEDULE_ collapse(2)
343  do j = jjs, jje
344  do i = iis, iie
345 #ifdef DEBUG
346 
347 #endif
348  flux(ks-1,i,j) = 0.0_rp ! k = KS
349 
350  ! if w>0; min(f,w*dz/dt)
351  ! else ; max(f,w*dz/dt) = -min(-f,-w*dz/dt)
352  sw = sign( 1.0_rp, mom(ks,i,j) )
353  flux(ks ,i,j) = sw * min( sw*flux(ks,i,j), sw*val(ks,i,j)*gsqrt(ks,i,j)*fdz(ks)/dtrk )
354 
355 
356  flux(ke-1,i,j) = 0.0_rp ! k = KE
357  flux(ke ,i,j) = 0.0_rp ! k = KE+1
358  enddo
359  enddo
360 
361  return
362  end subroutine atmos_dyn_fvm_fluxz_xyw_ud1
363 
364 
365  !-----------------------------------------------------------------------------
367  subroutine atmos_dyn_fvm_fluxj13_xyw_ud1( &
368  flux, &
369  mom, val, DENS, &
370  GSQRT, J13G, MAPF, &
371  CDZ, &
372  IIS, IIE, JJS, JJE )
373  implicit none
374 
375  real(RP), intent(out) :: flux (ka,ia,ja)
376  real(RP), intent(in) :: mom (ka,ia,ja)
377  real(RP), intent(in) :: val (ka,ia,ja)
378  real(RP), intent(in) :: DENS (ka,ia,ja)
379  real(RP), intent(in) :: GSQRT (ka,ia,ja)
380  real(RP), intent(in) :: J13G (ka,ia,ja)
381  real(RP), intent(in) :: MAPF ( ia,ja,2)
382  real(RP), intent(in) :: CDZ (ka)
383  integer, intent(in) :: IIS, IIE, JJS, JJE
384 
385  real(RP) :: vel
386  integer :: k, i, j
387  !---------------------------------------------------------------------------
388 
389  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
390  do j = jjs, jje
391  do i = iis, iie
392  do k = ks+1, ke-1
393  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i-1,j) ) ) &
394  / dens(k,i,j)
395  flux(k-1,i,j) = j13g(k,i,j) / mapf(i,j,+2) * vel &
396  * ( f1 * ( val(k,i,j)+val(k-1,i,j) ) - sign(f1,vel) * ( val(k,i,j)-val(k-1,i,j) ) )
397  enddo
398  enddo
399  enddo
400 
401  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
402  do j = jjs, jje
403  do i = iis, iie
404  flux(ks-1,i,j) = 0.0_rp
405 
406  flux(ke-1,i,j) = 0.0_rp
407  enddo
408  enddo
409 
410  return
411  end subroutine atmos_dyn_fvm_fluxj13_xyw_ud1
412 
413  !-----------------------------------------------------------------------------
415  subroutine atmos_dyn_fvm_fluxj23_xyw_ud1( &
416  flux, &
417  mom, val, DENS, &
418  GSQRT, J23G, MAPF, &
419  CDZ, &
420  IIS, IIE, JJS, JJE )
421  implicit none
422 
423  real(RP), intent(out) :: flux (ka,ia,ja)
424  real(RP), intent(in) :: mom (ka,ia,ja)
425  real(RP), intent(in) :: val (ka,ia,ja)
426  real(RP), intent(in) :: DENS (ka,ia,ja)
427  real(RP), intent(in) :: GSQRT (ka,ia,ja)
428  real(RP), intent(in) :: J23G (ka,ia,ja)
429  real(RP), intent(in) :: MAPF ( ia,ja,2)
430  real(RP), intent(in) :: CDZ (ka)
431  integer, intent(in) :: IIS, IIE, JJS, JJE
432 
433  real(RP) :: vel
434  integer :: k, i, j
435  !---------------------------------------------------------------------------
436 
437  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
438  do j = jjs, jje
439  do i = iis, iie
440  do k = ks+1, ke-1
441  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i,j-1) ) ) &
442  / dens(k,i,j)
443  flux(k-1,i,j) = j23g(k,i,j) / mapf(i,j,+1) * vel &
444  * ( f1 * ( val(k,i,j)+val(k-1,i,j) ) - sign(f1,vel) * ( val(k,i,j)-val(k-1,i,j) ) )
445  enddo
446  enddo
447  enddo
448 
449  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
450  do j = jjs, jje
451  do i = iis, iie
452  flux(ks-1,i,j) = 0.0_rp
453 
454  flux(ke-1,i,j) = 0.0_rp
455  enddo
456  enddo
457 
458  return
459  end subroutine atmos_dyn_fvm_fluxj23_xyw_ud1
460 
461 
462  !-----------------------------------------------------------------------------
464  subroutine atmos_dyn_fvm_fluxx_xyw_ud1( &
465  flux, &
466  mom, val, DENS, &
467  GSQRT, MAPF, &
468 
469  CDZ, &
470  IIS, IIE, JJS, JJE )
471  implicit none
472 
473  real(RP), intent(out) :: flux (ka,ia,ja)
474  real(RP), intent(in) :: mom (ka,ia,ja)
475  real(RP), intent(in) :: val (ka,ia,ja)
476  real(RP), intent(in) :: DENS (ka,ia,ja)
477  real(RP), intent(in) :: GSQRT (ka,ia,ja)
478  real(RP), intent(in) :: MAPF ( ia,ja,2)
479 
480  real(RP), intent(in) :: CDZ (ka)
481  integer, intent(in) :: IIS, IIE, JJS, JJE
482 
483  real(RP) :: vel
484  integer :: k, i, j
485  !---------------------------------------------------------------------------
486 
487  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
488  do j = jjs, jje
489  do i = iis-1, iie
490  do k = ks, ke-1
491 #ifdef DEBUG
492  call check( __line__, mom(k ,i,j) )
493  call check( __line__, mom(k+1,i,j) )
494 
495  call check( __line__, val(k,i,j) )
496  call check( __line__, val(k,i+1,j) )
497 
498 #endif
499  vel = ( f2h(k,1,i_uyz) &
500  * mom(k+1,i,j) &
501  + f2h(k,2,i_uyz) &
502  * mom(k,i,j) ) &
503  / ( f2h(k,1,i_xyz) &
504  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i+1,j) ) &
505  + f2h(k,2,i_xyz) &
506  * 0.5_rp * ( dens(k,i,j)+dens(k,i+1,j) ) )
507  flux(k,i,j) = gsqrt(k,i,j) / mapf(i,j,+2) * vel &
508  * ( f1 * ( val(k,i+1,j)+val(k,i,j) ) - sign(f1,vel) * ( val(k,i+1,j)-val(k,i,j) ) )
509  enddo
510  enddo
511  enddo
512 #ifdef DEBUG
513  k = iundef; i = iundef; j = iundef
514 #endif
515 
516  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
517  do j = jjs, jje
518  do i = iis-1, iie
519  flux(ke,i,j) = 0.0_rp
520  enddo
521  enddo
522 #ifdef DEBUG
523  k = iundef; i = iundef; j = iundef
524 #endif
525 
526  return
527  end subroutine atmos_dyn_fvm_fluxx_xyw_ud1
528 
529  !-----------------------------------------------------------------------------
531  subroutine atmos_dyn_fvm_fluxy_xyw_ud1( &
532  flux, &
533  mom, val, DENS, &
534  GSQRT, MAPF, &
535 
536  CDZ, &
537  IIS, IIE, JJS, JJE )
538  implicit none
539 
540  real(RP), intent(out) :: flux (ka,ia,ja)
541  real(RP), intent(in) :: mom (ka,ia,ja)
542  real(RP), intent(in) :: val (ka,ia,ja)
543  real(RP), intent(in) :: DENS (ka,ia,ja)
544  real(RP), intent(in) :: GSQRT (ka,ia,ja)
545  real(RP), intent(in) :: MAPF ( ia,ja,2)
546 
547  real(RP), intent(in) :: CDZ (ka)
548  integer, intent(in) :: IIS, IIE, JJS, JJE
549 
550  real(RP) :: vel
551  integer :: k, i, j
552  !---------------------------------------------------------------------------
553 
554  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
555  do j = jjs-1, jje
556  do i = iis, iie
557  do k = ks, ke-1
558 #ifdef DEBUG
559  call check( __line__, mom(k ,i,j) )
560  call check( __line__, mom(k+1,i,j) )
561 
562  call check( __line__, val(k,i,j) )
563  call check( __line__, val(k,i,j+1) )
564 
565 #endif
566  vel = ( f2h(k,1,i_xvz) &
567  * mom(k+1,i,j) &
568  + f2h(k,2,i_xvz) &
569  * mom(k,i,j) ) &
570  / ( f2h(k,1,i_xyz) &
571  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i,j+1) ) &
572  + f2h(k,2,i_xyz) &
573  * 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) ) )
574  flux(k,i,j) = gsqrt(k,i,j) / mapf(i,j,+1) * vel &
575  * ( f1 * ( val(k,i,j+1)+val(k,i,j) ) - sign(f1,vel) * ( val(k,i,j+1)-val(k,i,j) ) )
576  enddo
577  enddo
578  enddo
579 #ifdef DEBUG
580  k = iundef; i = iundef; j = iundef
581 #endif
582 
583  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
584  do j = jjs-1, jje
585  do i = iis, iie
586  flux(ke,i,j) = 0.0_rp
587  enddo
588  enddo
589 #ifdef DEBUG
590  k = iundef; i = iundef; j = iundef
591 #endif
592 
593  return
594  end subroutine atmos_dyn_fvm_fluxy_xyw_ud1
595 
596 
597  !-----------------------------------------------------------------------------
599  subroutine atmos_dyn_fvm_fluxz_uyz_ud1( &
600  flux, &
601  mom, val, DENS, &
602  GSQRT, J33G, &
603 
604  CDZ, &
605  IIS, IIE, JJS, JJE )
606  implicit none
607 
608  real(RP), intent(out) :: flux (ka,ia,ja)
609  real(RP), intent(in) :: mom (ka,ia,ja)
610  real(RP), intent(in) :: val (ka,ia,ja)
611  real(RP), intent(in) :: DENS (ka,ia,ja)
612  real(RP), intent(in) :: GSQRT (ka,ia,ja)
613  real(RP), intent(in) :: J33G
614 
615  real(RP), intent(in) :: CDZ (ka)
616  integer, intent(in) :: IIS, IIE, JJS, JJE
617 
618  real(RP) :: vel
619  integer :: k, i, j
620  !---------------------------------------------------------------------------
621 
622  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
623  do j = jjs, jje
624  do i = iis, iie
625  do k = ks, ke-1
626 #ifdef DEBUG
627  call check( __line__, mom(k,i,j) )
628  call check( __line__, mom(k,i+1,j) )
629 
630  call check( __line__, val(k,i,j) )
631  call check( __line__, val(k+1,i,j) )
632 
633 #endif
634  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i+1,j) ) ) &
635  / ( f2h(k,1,i_xyz) &
636  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i+1,j) ) &
637  + f2h(k,2,i_xyz) &
638  * 0.5_rp * ( dens(k,i,j)+dens(k,i+1,j) ) )
639  flux(k,i,j) = j33g * vel &
640  * ( f1 * ( val(k+1,i,j)+val(k,i,j) ) - sign(f1,vel) * ( val(k+1,i,j)-val(k,i,j) ) )
641  enddo
642  enddo
643  enddo
644 #ifdef DEBUG
645  k = iundef; i = iundef; j = iundef
646 #endif
647 
648  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
649  do j = jjs, jje
650  do i = iis, iie
651 #ifdef DEBUG
652 
653 #endif
654  flux(ks-1,i,j) = 0.0_rp
655 
656  flux(ke,i,j) = 0.0_rp
657  enddo
658  enddo
659 #ifdef DEBUG
660  k = iundef; i = iundef; j = iundef
661 #endif
662 
663  return
664  end subroutine atmos_dyn_fvm_fluxz_uyz_ud1
665 
666  !-----------------------------------------------------------------------------
668  subroutine atmos_dyn_fvm_fluxj13_uyz_ud1( &
669  flux, &
670  mom, val, DENS, &
671  GSQRT, J13G, MAPF, &
672  CDZ, &
673  IIS, IIE, JJS, JJE )
674  implicit none
675 
676  real(RP), intent(out) :: flux (ka,ia,ja)
677  real(RP), intent(in) :: mom (ka,ia,ja)
678  real(RP), intent(in) :: val (ka,ia,ja)
679  real(RP), intent(in) :: DENS (ka,ia,ja)
680  real(RP), intent(in) :: GSQRT (ka,ia,ja)
681  real(RP), intent(in) :: J13G (ka,ia,ja)
682  real(RP), intent(in) :: MAPF ( ia,ja,2)
683  real(RP), intent(in) :: CDZ (ka)
684  integer, intent(in) :: IIS, IIE, JJS, JJE
685 
686  real(RP) :: vel
687  integer :: k, i, j
688  !---------------------------------------------------------------------------
689 
690  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
691  do j = jjs, jje
692  do i = iis, iie
693  do k = ks, ke-1
694  vel = ( f2h(k,1,i_uyz) &
695  * mom(k+1,i,j) &
696  + f2h(k,2,i_uyz) &
697  * mom(k,i,j) ) &
698  / ( f2h(k,1,i_xyz) &
699  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i+1,j) ) &
700  + f2h(k,2,i_xyz) &
701  * 0.5_rp * ( dens(k,i,j)+dens(k,i+1,j) ) )
702  flux(k,i,j) = j13g(k,i,j) / mapf(i,j,+2) * vel &
703  * ( f1 * ( val(k+1,i,j)+val(k,i,j) ) - sign(f1,vel) * ( val(k+1,i,j)-val(k,i,j) ) )
704  enddo
705  enddo
706  enddo
707 
708  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
709  do j = jjs, jje
710  do i = iis, iie
711  flux(ks-1,i,j) = 0.0_rp
712 
713  flux(ke ,i,j) = 0.0_rp
714  enddo
715  enddo
716 
717  return
718  end subroutine atmos_dyn_fvm_fluxj13_uyz_ud1
719 
720  !-----------------------------------------------------------------------------
722  subroutine atmos_dyn_fvm_fluxj23_uyz_ud1( &
723  flux, &
724  mom, val, DENS, &
725  GSQRT, J23G, MAPF, &
726  CDZ, &
727  IIS, IIE, JJS, JJE )
728  implicit none
729 
730  real(RP), intent(out) :: flux (ka,ia,ja)
731  real(RP), intent(in) :: mom (ka,ia,ja)
732  real(RP), intent(in) :: val (ka,ia,ja)
733  real(RP), intent(in) :: DENS (ka,ia,ja)
734  real(RP), intent(in) :: GSQRT (ka,ia,ja)
735  real(RP), intent(in) :: J23G (ka,ia,ja)
736  real(RP), intent(in) :: MAPF ( ia,ja,2)
737  real(RP), intent(in) :: CDZ (ka)
738  integer, intent(in) :: IIS, IIE, JJS, JJE
739 
740  real(RP) :: vel
741  integer :: k, i, j
742  !---------------------------------------------------------------------------
743 
744  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
745  do j = jjs, jje
746  do i = iis, iie
747  do k = ks, ke-1
748  vel = ( f2h(k,1,i_xvz) &
749  * 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) ) &
750  + f2h(k,2,i_xvz) &
751  * 0.25_rp * ( mom(k,i,j)+mom(k,i+1,j)+mom(k,i,j-1)+mom(k,i+1,j-1) ) ) &
752  / ( f2h(k,1,i_xyz) &
753  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i+1,j) ) &
754  + f2h(k,2,i_xyz) &
755  * 0.5_rp * ( dens(k,i,j)+dens(k,i+1,j) ) )
756  flux(k,i,j) = j23g(k,i,j) / mapf(i,j,+1) * vel &
757  * ( f1 * ( val(k+1,i,j)+val(k,i,j) ) - sign(f1,vel) * ( val(k+1,i,j)-val(k,i,j) ) )
758  enddo
759  enddo
760  enddo
761 
762  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
763  do j = jjs, jje
764  do i = iis, iie
765  flux(ks-1,i,j) = 0.0_rp
766 
767  flux(ke ,i,j) = 0.0_rp
768  enddo
769  enddo
770 
771  return
772  end subroutine atmos_dyn_fvm_fluxj23_uyz_ud1
773 
774  !-----------------------------------------------------------------------------
776  subroutine atmos_dyn_fvm_fluxx_uyz_ud1( &
777  flux, &
778  mom, val, DENS, &
779  GSQRT, MAPF, &
780 
781  CDZ, &
782  IIS, IIE, JJS, JJE )
783  implicit none
784 
785  real(RP), intent(out) :: flux (ka,ia,ja)
786  real(RP), intent(in) :: mom (ka,ia,ja)
787  real(RP), intent(in) :: val (ka,ia,ja)
788  real(RP), intent(in) :: DENS (ka,ia,ja)
789  real(RP), intent(in) :: GSQRT (ka,ia,ja)
790  real(RP), intent(in) :: MAPF ( ia,ja,2)
791 
792  real(RP), intent(in) :: CDZ (ka)
793  integer, intent(in) :: IIS, IIE, JJS, JJE
794 
795  real(RP) :: vel
796  integer :: k, i, j
797  !---------------------------------------------------------------------------
798 
799  ! note that x-index is added by -1
800 
801  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
802  do j = jjs, jje
803  do i = iis, iie+1
804  do k = ks, ke
805 #ifdef DEBUG
806  call check( __line__, mom(k,i ,j) )
807  call check( __line__, mom(k,i-1,j) )
808 
809  call check( __line__, val(k,i-1,j) )
810  call check( __line__, val(k,i,j) )
811 
812 #endif
813  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i-1,j) ) ) &
814  / ( dens(k,i,j) )
815  flux(k,i-1,j) = gsqrt(k,i,j) / mapf(i,j,+2) * vel &
816  * ( f1 * ( val(k,i,j)+val(k,i-1,j) ) - sign(f1,vel) * ( val(k,i,j)-val(k,i-1,j) ) )
817  enddo
818  enddo
819  enddo
820 #ifdef DEBUG
821  k = iundef; i = iundef; j = iundef
822 #endif
823 
824  return
825  end subroutine atmos_dyn_fvm_fluxx_uyz_ud1
826 
827  !-----------------------------------------------------------------------------
829  subroutine atmos_dyn_fvm_fluxy_uyz_ud1( &
830  flux, &
831  mom, val, DENS, &
832  GSQRT, MAPF, &
833 
834  CDZ, &
835  IIS, IIE, JJS, JJE )
836  implicit none
837 
838  real(RP), intent(out) :: flux (ka,ia,ja)
839  real(RP), intent(in) :: mom (ka,ia,ja)
840  real(RP), intent(in) :: val (ka,ia,ja)
841  real(RP), intent(in) :: DENS (ka,ia,ja)
842  real(RP), intent(in) :: GSQRT (ka,ia,ja)
843  real(RP), intent(in) :: MAPF ( ia,ja,2)
844 
845  real(RP), intent(in) :: CDZ (ka)
846  integer, intent(in) :: IIS, IIE, JJS, JJE
847 
848  real(RP) :: vel
849  integer :: k, i, j
850  !---------------------------------------------------------------------------
851 
852  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
853  do j = jjs-1, jje
854  do i = iis, iie
855  do k = ks, ke
856 #ifdef DEBUG
857  call check( __line__, mom(k,i ,j) )
858  call check( __line__, mom(k,i-1,j) )
859 
860  call check( __line__, val(k,i,j) )
861  call check( __line__, val(k,i,j+1) )
862 
863 #endif
864  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i+1,j) ) ) &
865  / ( 0.25_rp * ( dens(k,i,j)+dens(k,i+1,j)+dens(k,i,j+1)+dens(k,i+1,j+1) ) )
866  flux(k,i,j) = gsqrt(k,i,j) / mapf(i,j,+1) * vel &
867  * ( f1 * ( val(k,i,j+1)+val(k,i,j) ) - sign(f1,vel) * ( val(k,i,j+1)-val(k,i,j) ) )
868  enddo
869  enddo
870  enddo
871 #ifdef DEBUG
872  k = iundef; i = iundef; j = iundef
873 #endif
874 
875  return
876  end subroutine atmos_dyn_fvm_fluxy_uyz_ud1
877 
878 
879 
880  !-----------------------------------------------------------------------------
882  subroutine atmos_dyn_fvm_fluxz_xvz_ud1( &
883  flux, &
884  mom, val, DENS, &
885  GSQRT, J33G, &
886 
887  CDZ, &
888  IIS, IIE, JJS, JJE )
889  implicit none
890 
891  real(RP), intent(out) :: flux (ka,ia,ja)
892  real(RP), intent(in) :: mom (ka,ia,ja)
893  real(RP), intent(in) :: val (ka,ia,ja)
894  real(RP), intent(in) :: DENS (ka,ia,ja)
895  real(RP), intent(in) :: GSQRT (ka,ia,ja)
896  real(RP), intent(in) :: J33G
897 
898  real(RP), intent(in) :: CDZ (ka)
899  integer, intent(in) :: IIS, IIE, JJS, JJE
900 
901  real(RP) :: vel
902  integer :: k, i, j
903  !---------------------------------------------------------------------------
904 
905  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
906  do j = jjs, jje
907  do i = iis, iie
908  do k = ks, ke-1
909 #ifdef DEBUG
910  call check( __line__, mom(k,i,j) )
911  call check( __line__, mom(k,i,j+1) )
912 
913  call check( __line__, val(k,i,j) )
914  call check( __line__, val(k+1,i,j) )
915 
916 #endif
917  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i,j+1) ) ) &
918  / ( f2h(k,1,i_xyz) &
919  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i,j+1) ) &
920  + f2h(k,2,i_xyz) &
921  * 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) ) )
922  flux(k,i,j) = j33g * vel &
923  * ( f1 * ( val(k+1,i,j)+val(k,i,j) ) - sign(f1,vel) * ( val(k+1,i,j)-val(k,i,j) ) )
924  enddo
925  enddo
926  enddo
927 #ifdef DEBUG
928  k = iundef; i = iundef; j = iundef
929 #endif
930 
931  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
932  do j = jjs, jje
933  do i = iis, iie
934 #ifdef DEBUG
935 
936 #endif
937  flux(ks-1,i,j) = 0.0_rp
938 
939  flux(ke,i,j) = 0.0_rp
940  enddo
941  enddo
942 #ifdef DEBUG
943  k = iundef; i = iundef; j = iundef
944 #endif
945 
946  return
947  end subroutine atmos_dyn_fvm_fluxz_xvz_ud1
948 
949  !-----------------------------------------------------------------------------
951  subroutine atmos_dyn_fvm_fluxj13_xvz_ud1( &
952  flux, &
953  mom, val, DENS, &
954  GSQRT, J13G, MAPF, &
955  CDZ, &
956  IIS, IIE, JJS, JJE )
957  implicit none
958 
959  real(RP), intent(out) :: flux (ka,ia,ja)
960  real(RP), intent(in) :: mom (ka,ia,ja)
961  real(RP), intent(in) :: val (ka,ia,ja)
962  real(RP), intent(in) :: DENS (ka,ia,ja)
963  real(RP), intent(in) :: GSQRT (ka,ia,ja)
964  real(RP), intent(in) :: J13G (ka,ia,ja)
965  real(RP), intent(in) :: MAPF ( ia,ja,2)
966  real(RP), intent(in) :: CDZ (ka)
967  integer, intent(in) :: IIS, IIE, JJS, JJE
968 
969  real(RP) :: vel
970  integer :: k, i, j
971  !---------------------------------------------------------------------------
972 
973  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
974  do j = jjs, jje
975  do i = iis, iie
976  do k = ks, ke-1
977  vel = ( f2h(k,1,i_uyz) &
978  * 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) ) &
979  + f2h(k,2,i_uyz) &
980  * 0.25_rp * ( mom(k,i,j)+mom(k,i-1,j)+mom(k,i,j+1)+mom(k,i-1,j+1) ) ) &
981  / ( f2h(k,1,i_xyz) &
982  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i,j+1) ) &
983  + f2h(k,2,i_xyz) &
984  * 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) ) )
985  flux(k,i,j) = j13g(k,i,j) / mapf(i,j,+2) * vel &
986  * ( f1 * ( val(k+1,i,j)+val(k,i,j) ) - sign(f1,vel) * ( val(k+1,i,j)-val(k,i,j) ) )
987  enddo
988  enddo
989  enddo
990 
991  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
992  do j = jjs, jje
993  do i = iis, iie
994  flux(ks-1,i,j) = 0.0_rp
995 
996  flux(ke ,i,j) = 0.0_rp
997  enddo
998  enddo
999 
1000  return
1001  end subroutine atmos_dyn_fvm_fluxj13_xvz_ud1
1002 
1003  !-----------------------------------------------------------------------------
1005  subroutine atmos_dyn_fvm_fluxj23_xvz_ud1( &
1006  flux, &
1007  mom, val, DENS, &
1008  GSQRT, J23G, MAPF, &
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) :: J23G (ka,ia,ja)
1019  real(RP), intent(in) :: MAPF ( ia,ja,2)
1020  real(RP), intent(in) :: CDZ (ka)
1021  integer, intent(in) :: IIS, IIE, JJS, JJE
1022 
1023  real(RP) :: vel
1024  integer :: k, i, j
1025  !---------------------------------------------------------------------------
1026 
1027  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1028  do j = jjs, jje
1029  do i = iis, iie
1030  do k = ks, ke-1
1031  vel = ( f2h(k,1,i_xvz) &
1032  * mom(k+1,i,j) &
1033  + f2h(k,2,i_xvz) &
1034  * mom(k,i,j) ) &
1035  / ( f2h(k,1,i_xyz) &
1036  * 0.5_rp * ( dens(k+1,i,j)+dens(k+1,i,j+1) ) &
1037  + f2h(k,2,i_xyz) &
1038  * 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) ) )
1039  flux(k,i,j) = j23g(k,i,j) / mapf(i,j,+1) * vel &
1040  * ( f1 * ( val(k+1,i,j)+val(k,i,j) ) - sign(f1,vel) * ( val(k+1,i,j)-val(k,i,j) ) )
1041  enddo
1042  enddo
1043  enddo
1044 
1045  !$omp parallel do private(i,j) OMP_SCHEDULE_ collapse(2)
1046  do j = jjs, jje
1047  do i = iis, iie
1048  flux(ks-1,i,j) = 0.0_rp
1049 
1050  flux(ke ,i,j) = 0.0_rp
1051  enddo
1052  enddo
1053 
1054  return
1055  end subroutine atmos_dyn_fvm_fluxj23_xvz_ud1
1056 
1057  !-----------------------------------------------------------------------------
1059  subroutine atmos_dyn_fvm_fluxx_xvz_ud1( &
1060  flux, &
1061  mom, val, DENS, &
1062  GSQRT, MAPF, &
1063 
1064  CDZ, &
1065  IIS, IIE, JJS, JJE )
1066  implicit none
1067 
1068  real(RP), intent(out) :: flux (ka,ia,ja)
1069  real(RP), intent(in) :: mom (ka,ia,ja)
1070  real(RP), intent(in) :: val (ka,ia,ja)
1071  real(RP), intent(in) :: DENS (ka,ia,ja)
1072  real(RP), intent(in) :: GSQRT (ka,ia,ja)
1073  real(RP), intent(in) :: MAPF ( ia,ja,2)
1074 
1075  real(RP), intent(in) :: CDZ (ka)
1076  integer, intent(in) :: IIS, IIE, JJS, JJE
1077 
1078  real(RP) :: vel
1079  integer :: k, i, j
1080  !---------------------------------------------------------------------------
1081 
1082  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1083  do j = jjs, jje
1084  do i = iis-1, iie
1085  do k = ks, ke
1086 #ifdef DEBUG
1087  call check( __line__, mom(k,i ,j) )
1088  call check( __line__, mom(k,i,j-1) )
1089 
1090  call check( __line__, val(k,i,j) )
1091  call check( __line__, val(k,i+1,j) )
1092 
1093 #endif
1094  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i,j+1) ) ) &
1095  / ( 0.25_rp * ( dens(k,i,j)+dens(k,i+1,j)+dens(k,i,j+1)+dens(k,i+1,j+1) ) )
1096  flux(k,i,j) = gsqrt(k,i,j) / mapf(i,j,+2) * vel &
1097  * ( f1 * ( val(k,i+1,j)+val(k,i,j) ) - sign(f1,vel) * ( val(k,i+1,j)-val(k,i,j) ) )
1098  enddo
1099  enddo
1100  enddo
1101 #ifdef DEBUG
1102  k = iundef; i = iundef; j = iundef
1103 #endif
1104 
1105  return
1106  end subroutine atmos_dyn_fvm_fluxx_xvz_ud1
1107 
1108  !-----------------------------------------------------------------------------
1110  subroutine atmos_dyn_fvm_fluxy_xvz_ud1( &
1111  flux, &
1112  mom, val, DENS, &
1113  GSQRT, MAPF, &
1114 
1115  CDZ, &
1116  IIS, IIE, JJS, JJE )
1117  implicit none
1118 
1119  real(RP), intent(out) :: flux (ka,ia,ja)
1120  real(RP), intent(in) :: mom (ka,ia,ja)
1121  real(RP), intent(in) :: val (ka,ia,ja)
1122  real(RP), intent(in) :: DENS (ka,ia,ja)
1123  real(RP), intent(in) :: GSQRT (ka,ia,ja)
1124  real(RP), intent(in) :: MAPF ( ia,ja,2)
1125 
1126  real(RP), intent(in) :: CDZ (ka)
1127  integer, intent(in) :: IIS, IIE, JJS, JJE
1128 
1129  real(RP) :: vel
1130  integer :: k, i, j
1131  !---------------------------------------------------------------------------
1132 
1133  ! note that y-index is added by -1
1134 
1135  !$omp parallel do private(i,j,k) OMP_SCHEDULE_ collapse(2)
1136  do j = jjs, jje+1
1137  do i = iis, iie
1138  do k = ks, ke
1139 #ifdef DEBUG
1140  call check( __line__, mom(k,i ,j) )
1141  call check( __line__, mom(k,i,j-1) )
1142 
1143  call check( __line__, val(k,i,j-1) )
1144  call check( __line__, val(k,i,j) )
1145 
1146 #endif
1147  vel = ( 0.5_rp * ( mom(k,i,j)+mom(k,i,j-1) ) ) &
1148  / ( dens(k,i,j) )
1149  flux(k,i,j-1) = gsqrt(k,i,j) / mapf(i,j,+1) * vel &
1150  * ( f1 * ( val(k,i,j)+val(k,i,j-1) ) - sign(f1,vel) * ( val(k,i,j)-val(k,i,j-1) ) )
1151  enddo
1152  enddo
1153  enddo
1154 #ifdef DEBUG
1155  k = iundef; i = iundef; j = iundef
1156 #endif
1157 
1158  return
1159  end subroutine atmos_dyn_fvm_fluxy_xvz_ud1
1160 
1161 
1162 
1163 
1164 
1165 
1166 
1168 
1169 !--
1170 ! vi:set readonly sw=4 ts=8
1171 !
1172 !Local Variables:
1173 !mode: f90
1174 !buffer-read-only: t
1175 !End:
1176 !
1177 !++
module DEBUG
Definition: scale_debug.F90:13
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_fluxz_xyw_ud1(flux, mom, val, DENS, GSQRT, J33G, CDZ, FDZ, dtrk, IIS, IIE, JJS, JJE)
calculation z-flux at XYW
subroutine, public atmos_dyn_fvm_fluxy_xvz_ud1(flux, mom, val, DENS, GSQRT, MAPF, CDZ, IIS, IIE, JJS, JJE)
calculation Y-flux at XV
subroutine, public atmos_dyn_fvm_fluxy_xyz_ud1(flux, mflx, val, GSQRT, CDZ, IIS, IIE, JJS, JJE)
calculation Y-flux at XYZ
subroutine, public atmos_dyn_fvm_fluxz_uyz_ud1(flux, mom, val, DENS, GSQRT, J33G, CDZ, IIS, IIE, JJS, JJE)
calculation z-flux at UY
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
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
real(rp), public const_undef
Definition: scale_const.F90:43
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 grid index
module TRACER
module Index
Definition: scale_index.F90:14
integer, public ia
of x whole cells (local, with HALO)
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_fluxx_xyw_ud1(flux, mom, val, DENS, GSQRT, MAPF, CDZ, IIS, IIE, JJS, JJE)
calculation X-flux at XYW
subroutine, public atmos_dyn_fvm_fluxy_uyz_ud1(flux, mom, val, DENS, GSQRT, MAPF, CDZ, IIS, IIE, JJS, JJE)
calculation Y-flux at UY
integer, public ka
of z whole cells (local, with HALO)
subroutine, public atmos_dyn_fvm_fluxx_uyz_ud1(flux, mom, val, DENS, GSQRT, MAPF, CDZ, IIS, IIE, JJS, JJE)
calculation X-flux at UY
subroutine, public atmos_dyn_fvm_fluxx_xvz_ud1(flux, mom, val, DENS, GSQRT, MAPF, CDZ, IIS, IIE, JJS, JJE)
calculation X-flux at XV
subroutine, public atmos_dyn_fvm_fluxz_xyz_ud1(flux, mflx, val, GSQRT, CDZ, IIS, IIE, JJS, JJE)
calculation z-flux at XYZ
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:40
module PROCESS
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
module profiler
Definition: scale_prof.F90:10
subroutine, public atmos_dyn_fvm_fluxy_xyw_ud1(flux, mom, val, DENS, GSQRT, MAPF, CDZ, IIS, IIE, JJS, JJE)
calculation Y-flux at XYW
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
subroutine, public atmos_dyn_fvm_fluxx_xyz_ud1(flux, mflx, val, GSQRT, CDZ, IIS, IIE, JJS, JJE)
calculation X-flux at XYZ
subroutine, public atmos_dyn_fvm_flux_valuew_z_ud1(valW, mflx, val, GSQRT, CDZ)
value at XYW
subroutine, public atmos_dyn_fvm_fluxz_xvz_ud1(flux, mom, val, DENS, GSQRT, J33G, CDZ, IIS, IIE, JJS, JJE)
calculation z-flux at XV
module scale_atmos_dyn_fvm_flux_ud1
integer, public ja
of y whole cells (local, with HALO)