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