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