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