SCALE-RM
scale_atmos_phy_tb_common.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
11 #include "scalelib.h"
13  !-----------------------------------------------------------------------------
14  !
15  !++ used modules
16  !
17  use scale_precision
18  use scale_io
19  use scale_prof
21  use scale_tracer
22 
23 #ifdef DEBUG
24  use scale_debug, only: &
25  check
26  use scale_const, only: &
27  undef => const_undef, &
28  iundef => const_undef2
29 #endif
30  !-----------------------------------------------------------------------------
31  implicit none
32  private
33  !-----------------------------------------------------------------------------
34  !
35  !++ Public procedure
36  !
43 
44  !-----------------------------------------------------------------------------
45  !
46  !++ Public parameters & variables
47  !
48  !-----------------------------------------------------------------------------
49  !
50  !++ Private procedure
51  !
52  !-----------------------------------------------------------------------------
53  !
54  !++ Private parameters & variables
55  !
56  !-----------------------------------------------------------------------------
57 contains
58  !-----------------------------------------------------------------------------
60  S33_C, S11_C, S22_C, &
61  S31_C, S12_C, S23_C, &
62  S12_Z, S23_X, S31_Y, &
63  S2, &
64  DENS, MOMZ, MOMX, MOMY, &
65  GSQRT, J13G, J23G, J33G, MAPF )
67  fdz => atmos_grid_cartesc_fdz, &
68  fdx => atmos_grid_cartesc_fdx, &
69  fdy => atmos_grid_cartesc_fdy, &
70  rcdz => atmos_grid_cartesc_rcdz, &
71  rcdx => atmos_grid_cartesc_rcdx, &
72  rcdy => atmos_grid_cartesc_rcdy, &
73  rfdz => atmos_grid_cartesc_rfdz, &
74  rfdx => atmos_grid_cartesc_rfdx, &
76  implicit none
77 
78  real(rp), intent(out) :: s33_c (ka,ia,ja) ! (cell center)
79  real(rp), intent(out) :: s11_c (ka,ia,ja)
80  real(rp), intent(out) :: s22_c (ka,ia,ja)
81  real(rp), intent(out) :: s31_c (ka,ia,ja)
82  real(rp), intent(out) :: s12_c (ka,ia,ja)
83  real(rp), intent(out) :: s23_c (ka,ia,ja)
84  real(rp), intent(out) :: s12_z (ka,ia,ja) ! (z edge or x-y plane)
85  real(rp), intent(out) :: s23_x (ka,ia,ja) ! (x edge or y-z plane)
86  real(rp), intent(out) :: s31_y (ka,ia,ja) ! (y edge or z-x plane)
87  real(rp), intent(out) :: s2 (ka,ia,ja) ! |S|^2
88 
89  real(rp), intent(in) :: dens (ka,ia,ja)
90  real(rp), intent(in) :: momz (ka,ia,ja)
91  real(rp), intent(in) :: momx (ka,ia,ja)
92  real(rp), intent(in) :: momy (ka,ia,ja)
93 
94  real(rp), intent(in) :: gsqrt (ka,ia,ja,7)
95  real(rp), intent(in) :: j13g (ka,ia,ja,7)
96  real(rp), intent(in) :: j23g (ka,ia,ja,7)
97  real(rp), intent(in) :: j33g
98  real(rp), intent(in) :: mapf (ia,ja,2,4)
99 
100  ! velocity
101  real(rp) :: velz_c (ka,ia,ja)
102  real(rp) :: velz_xy(ka,ia,ja)
103  real(rp) :: velx_c (ka,ia,ja)
104  real(rp) :: velx_yz(ka,ia,ja)
105  real(rp) :: vely_c (ka,ia,ja)
106  real(rp) :: vely_zx(ka,ia,ja)
107 
108  ! work space
109  real(rp) :: work_v(ka,ia,ja) ! work space (vertex)
110  real(rp) :: work_z(ka,ia,ja) ! (z edge or x-y plane)
111  real(rp) :: work_x(ka,ia,ja) ! (x edge or y-z plane)
112  real(rp) :: work_y(ka,ia,ja) ! (y edge or z-x plane)
113 
114  integer :: iis, iie, jjs, jje
115  integer :: k, i, j
116 
117  !$acc data copyout(S33_C, S11_C, S22_C, S31_C, S12_C, S23_C, S12_Z, S23_X, S31_y, S2) &
118  !$acc copyin(DENS, MOMZ, MOMX, MOMY, GSQRT, J13G, J23G, MAPF, &
119  !$acc FDZ, FDX, FDY, RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY) &
120  !$acc create(VELZ_C, VELZ_XY, VELX_C, VELX_YZ, VELY_C, VELY_ZX, &
121  !$acc WORK_V, WORK_Z, WORK_X, WORK_Y)
122 
123 #ifdef DEBUG
124  !$acc kernels
125  s33_c(:,:,:) = undef
126  s11_c(:,:,:) = undef
127  s22_c(:,:,:) = undef
128  s31_c(:,:,:) = undef
129  s12_c(:,:,:) = undef
130  s23_c(:,:,:) = undef
131  s12_z(:,:,:) = undef
132  s23_x(:,:,:) = undef
133  s31_y(:,:,:) = undef
134  s2(:,:,:) = undef
135 
136  velz_c(:,:,:) = undef
137  velz_xy(:,:,:) = undef
138  velx_c(:,:,:) = undef
139  velx_yz(:,:,:) = undef
140  vely_c(:,:,:) = undef
141  vely_zx(:,:,:) = undef
142 
143  work_v(:,:,:) = undef
144  work_z(:,:,:) = undef
145  work_x(:,:,:) = undef
146  work_y(:,:,:) = undef
147  !$end kernels
148 #endif
149 
150  ! momentum -> velocity
151  !$omp parallel do
152  !$acc kernels
153  do j = js-1, je+1
154  do i = is-1, ie+1
155  do k = ks, ke-1
156 #ifdef DEBUG
157  call check( __line__, momz(k,i,j) )
158  call check( __line__, dens(k+1,i,j) )
159  call check( __line__, dens(k,i,j) )
160 #endif
161  velz_xy(k,i,j) = 2.0_rp * momz(k,i,j) / ( dens(k+1,i,j)+dens(k,i,j) )
162  enddo
163  enddo
164  enddo
165  !$acc end kernels
166 #ifdef DEBUG
167  i = iundef; j = iundef; k = iundef
168 #endif
169  !$omp parallel do
170  !$acc kernels
171  do j = js-1, je+1
172  do i = is-1, ie+1
173  velz_xy(ke,i,j) = 0.0_rp
174  enddo
175  enddo
176  !$acc end kernels
177 #ifdef DEBUG
178  i = iundef; j = iundef; k = iundef
179 #endif
180  !$omp parallel do
181  !$acc kernels
182  do j = js-2, je+2
183  do i = is-2, ie+2
184  do k = ks+1, ke
185 #ifdef DEBUG
186  call check( __line__, momz(k,i,j) )
187  call check( __line__, momz(k-1,i,j) )
188  call check( __line__, dens(k,i,j) )
189 #endif
190  velz_c(k,i,j) = 0.5_rp * ( momz(k,i,j) + momz(k-1,i,j) ) / dens(k,i,j)
191  enddo
192  enddo
193  enddo
194  !$acc end kernels
195 #ifdef DEBUG
196  i = iundef; j = iundef; k = iundef
197 #endif
198  !$omp parallel do
199  !$acc kernels
200  do j = js-2, je+2
201  do i = is-2, ie+2
202 #ifdef DEBUG
203  call check( __line__, momz(ks,i,j) )
204  call check( __line__, dens(ks,i,j) )
205 #endif
206  velz_c(ks,i,j) = 0.5_rp * momz(ks,i,j) / dens(ks,i,j) ! MOMZ(KS-1,i,j) = 0
207  enddo
208  enddo
209  !$acc end kernels
210 #ifdef DEBUG
211  i = iundef; j = iundef; k = iundef
212 #endif
213 
214  !$omp parallel do
215  !$acc kernels
216  do j = js-1, je+1
217  do i = is-2, ie+1
218  do k = ks, ke
219 #ifdef DEBUG
220  call check( __line__, momx(k,i,j) )
221  call check( __line__, dens(k,i+1,j) )
222  call check( __line__, dens(k,i,j) )
223 #endif
224  velx_yz(k,i,j) = 2.0_rp * momx(k,i,j) / ( dens(k,i+1,j)+dens(k,i,j) )
225  enddo
226  enddo
227  enddo
228  !$acc end kernels
229 #ifdef DEBUG
230  i = iundef; j = iundef; k = iundef
231 #endif
232  !$omp parallel do
233  !$acc kernels
234  do j = js-1, je+1
235  do i = is-2, ie+1
236  velx_yz(ke+1,i,j) = 0.0_rp
237  enddo
238  enddo
239  !$acc end kernels
240 #ifdef DEBUG
241  i = iundef; j = iundef; k = iundef
242 #endif
243  !$omp parallel do
244  !$acc kernels
245  do j = js-2, je+2
246  do i = is-1, ie+2
247  do k = ks, ke
248 #ifdef DEBUG
249  call check( __line__, momx(k,i,j) )
250  call check( __line__, momx(k,i-1,j) )
251  call check( __line__, dens(k,i,j) )
252 #endif
253  velx_c(k,i,j) = 0.5_rp * ( momx(k,i,j) + momx(k,i-1,j) ) / dens(k,i,j)
254  enddo
255  enddo
256  enddo
257  !$acc end kernels
258 #ifdef DEBUG
259  i = iundef; j = iundef; k = iundef
260 #endif
261  !$omp parallel do default(none) &
262  !$omp shared(JS,JE,IS,IE,KS,KE,MOMY,DENS,VELY_ZX) &
263  !$omp private(i,j,k) OMP_SCHEDULE_ collapse(2)
264  !$acc kernels
265  do j = js-2, je+1
266  do i = is-1, ie+1
267  do k = ks, ke
268 #ifdef DEBUG
269  call check( __line__, momy(k,i,j) )
270  call check( __line__, dens(k,i,j+1) )
271  call check( __line__, dens(k,i,j) )
272 #endif
273  vely_zx(k,i,j) = 2.0_rp * momy(k,i,j) / ( dens(k,i,j+1)+dens(k,i,j) )
274  enddo
275  enddo
276  enddo
277  !$acc end kernels
278 #ifdef DEBUG
279  i = iundef; j = iundef; k = iundef
280 #endif
281  !$omp parallel do
282  !$acc kernels
283  do j = js-2, je+1
284  do i = is-1, ie+1
285  vely_zx(ke+1,i,j) = 0.0_rp
286  enddo
287  enddo
288  !$acc end kernels
289 #ifdef DEBUG
290  i = iundef; j = iundef; k = iundef
291 #endif
292  !$omp parallel do
293  !$acc kernels
294  do j = js-1, je+2
295  do i = is-2, ie+2
296  do k = ks, ke
297 #ifdef DEBUG
298  call check( __line__, momy(k,i,j) )
299  call check( __line__, momy(k,i,j-1) )
300  call check( __line__, dens(k,i,j) )
301 #endif
302  vely_c(k,i,j) = 0.5_rp * ( momy(k,i,j) + momy(k,i,j-1) ) / dens(k,i,j)
303  enddo
304  enddo
305  enddo
306  !$acc end kernels
307 #ifdef DEBUG
308  i = iundef; j = iundef; k = iundef
309 #endif
310 
311  do jjs = js, je, jblock
312  jje = jjs+jblock-1
313  do iis = is, ie, iblock
314  iie = iis+iblock-1
315 
316 #ifdef DEBUG
317  work_z(:,:,:) = undef; work_x(:,:,:) = undef; work_y(:,:,:) = undef; work_v(:,:,:) = undef
318 #endif
319  ! w
320  ! (x-y plane; x,y,w)
321  ! WORK_Z = VELZ_XY
322  ! (y-z plane; u,y,z)
323  !$omp parallel do
324  !$acc kernels
325  do j = jjs-1, jje+1
326  do i = iis-1, iie+1
327  do k = ks, ke
328 #ifdef DEBUG
329  call check( __line__, velz_c(k,i+1,j) )
330  call check( __line__, velz_c(k,i,j) )
331 #endif
332  work_x(k,i,j) = 0.5_rp * ( velz_c(k,i+1,j) + velz_c(k,i,j) )
333  enddo
334  enddo
335  enddo
336  !$acc end kernels
337 #ifdef DEBUG
338  i = iundef; j = iundef; k = iundef
339 #endif
340  !$omp parallel do
341  !$acc kernels
342  do j = jjs-1, jje+1
343  do i = iis-1, iie+1
344  work_x(ke+1,i,j) = 0.0_rp
345  enddo
346  enddo
347  !$acc end kernels
348 #ifdef DEBUG
349  i = iundef; j = iundef; k = iundef
350 #endif
351  ! (z-x plane; x,v,z)
352  !$omp parallel do
353  !$acc kernels
354  do j = jjs-1, jje+1
355  do i = iis-1, iie+1
356  do k = ks, ke
357 #ifdef DEBUG
358  call check( __line__, velz_c(k,i,j+1) )
359  call check( __line__, velz_c(k,i,j) )
360 #endif
361  work_y(k,i,j) = 0.5_rp * ( velz_c(k,i,j+1) + velz_c(k,i,j) )
362  enddo
363  enddo
364  enddo
365  !$acc end kernels
366 #ifdef DEBUG
367  i = iundef; j = iundef; k = iundef
368 #endif
369  !$omp parallel do
370  !$acc kernels
371  do j = jjs-1, jje+1
372  do i = iis-1, iie+1
373  work_y(ke+1,i,j) = 0.0_rp
374  enddo
375  enddo
376  !$acc end kernels
377 #ifdef DEBUG
378  i = iundef; j = iundef; k = iundef
379 #endif
380 
381  ! dw/dz
382  ! (cell center; x,y,z)
383  !$omp parallel do
384  !$acc kernels
385  do j = jjs-1, jje+1
386  do i = iis-1, iie+1
387  do k = ks+1, ke
388 #ifdef DEBUG
389  call check( __line__, velz_xy(k,i,j) )
390  call check( __line__, velz_xy(k-1,i,j) )
391  call check( __line__, rcdz(k) )
392 #endif
393  s33_c(k,i,j) = ( velz_xy(k,i,j) - velz_xy(k-1,i,j) ) * rcdz(k) &
394  * j33g / gsqrt(k,i,j,i_xyz)
395  enddo
396  enddo
397  enddo
398  !$acc end kernels
399 #ifdef DEBUG
400  i = iundef; j = iundef; k = iundef
401 #endif
402  !$omp parallel do
403  !$acc kernels
404  do j = jjs-1, jje+1
405  do i = iis-1, iie+1
406 #ifdef DEBUG
407  call check( __line__, velz_xy(ks,i,j) )
408  call check( __line__, gsqrt(ks,i,j,i_xyz) )
409  call check( __line__, rcdz(ks) )
410 #endif
411  s33_c(ks,i,j) = velz_xy(ks,i,j) * rcdz(ks) & ! VELZ_XY(KS-1,i,j) == 0
412  * j33g / gsqrt(ks,i,j,i_xyz)
413  enddo
414  enddo
415  !$acc end kernels
416 #ifdef DEBUG
417  i = iundef; j = iundef; k = iundef
418 #endif
419 
420  ! 1/2 * dw/dx
421  ! (cell center; x,y,z)
422  !$omp parallel do
423  !$acc kernels
424 !OCL LOOP_NOINTERCHANGE
425  do j = jjs-1, jje+1
426  do i = iis-1, iie+1
427  do k = ks+1, ke-1
428 #ifdef DEBUG
429  call check( __line__, velz_c(k,i+1,j) )
430  call check( __line__, velz_c(k,i-1,j) )
431  call check( __line__, gsqrt(k,i+1,j,i_xyz) )
432  call check( __line__, gsqrt(k,i-1,j,i_xyz) )
433  call check( __line__, velz_xy(k,i,j) )
434  call check( __line__, velz_xy(k-1,i,j) )
435  call check( __line__, j13g(k,i,j,i_xyw) )
436  call check( __line__, j13g(k-1,i,j,i_xyw) )
437  call check( __line__, fdx(i) )
438  call check( __line__, fdx(i-1) )
439 #endif
440  s31_c(k,i,j) = 0.5_rp * ( &
441  ( gsqrt(k,i+1,j,i_xyz)*velz_c(k,i+1,j) - gsqrt(k,i-1,j,i_xyz)*velz_c(k,i-1,j) ) / ( fdx(i) + fdx(i-1) ) &
442  + ( j13g(k,i,j,i_xyw)*velz_xy(k,i,j) - j13g(k-1,i,j,i_xyw)*velz_xy(k-1,i,j) ) * rcdz(k) &
443  ) * mapf(i,j,1,i_xy)
444 
445  enddo
446  enddo
447  enddo
448  !$acc end kernels
449 #ifdef DEBUG
450  i = iundef; j = iundef; k = iundef
451 #endif
452  !$omp parallel do
453  !$acc kernels
454  do j = jjs-1, jje+1
455  do i = iis-1, iie+1
456 #ifdef DEBUG
457  call check( __line__, velz_c(ks,i+1,j) )
458  call check( __line__, velz_c(ks,i-1,j) )
459  call check( __line__, gsqrt(ks,i+1,j,i_xyz) )
460  call check( __line__, gsqrt(ks,i-1,j,i_xyz) )
461  call check( __line__, velz_xy(ks,i,j) )
462  call check( __line__, j13g(ks,i,j,i_xyw) )
463  call check( __line__, velz_c(ke,i+1,j) )
464  call check( __line__, velz_c(ke,i-1,j) )
465  call check( __line__, gsqrt(ke,i+1,j,i_xyz) )
466  call check( __line__, gsqrt(ke,i-1,j,i_xyz) )
467  call check( __line__, velz_xy(ke,i,j) )
468  call check( __line__, j13g(ke,i,j,i_xyw) )
469  call check( __line__, fdx(i) )
470  call check( __line__, fdx(i-1) )
471 #endif
472  s31_c(ks,i,j) = 0.5_rp * ( &
473  ( gsqrt(ks,i+1,j,i_xyz)*velz_c(ks,i+1,j) - gsqrt(ks,i-1,j,i_xyz)*velz_c(ks,i-1,j) ) / ( fdx(i) + fdx(i-1) ) &
474  + ( j13g(ks,i,j,i_xyw)*velz_xy(ks,i,j) ) * rcdz(ks) &
475  ) * mapf(i,j,1,i_xy)
476  s31_c(ke,i,j) = 0.5_rp * ( &
477  ( gsqrt(ke,i+1,j,i_xyz)*velz_c(ke,i+1,j) - gsqrt(ke,i-1,j,i_xyz)*velz_c(ke,i-1,j) ) / ( fdx(i) + fdx(i-1) ) &
478  - ( j13g(ke-1,i,j,i_xyw)*velz_xy(ke-1,i,j) ) * rcdz(ke) &
479  ) * mapf(i,j,1,i_xy)
480  enddo
481  enddo
482  !$acc end kernels
483 #ifdef DEBUG
484  i = iundef; j = iundef; k = iundef
485 #endif
486 
487  ! (y edge, u,y,w)
488  !$omp parallel do
489  !$acc kernels
490  do j = jjs , jje
491  do i = iis-1, iie
492  do k = ks, ke-1
493 #ifdef DEBUG
494  call check( __line__, velz_xy(k,i+1,j) )
495  call check( __line__, velz_xy(k,i,j) )
496  call check( __line__, rfdx(i) )
497 #endif
498  s31_y(k,i,j) = 0.5_rp * ( &
499  ( gsqrt(k,i+1,j,i_xyw)*velz_xy(k,i+1,j) - gsqrt(k,i,j,i_xyw)*velz_xy(k,i,j) ) * rfdx(i) &
500  + ( j13g(k+1,i,j,i_uyz)*work_x(k+1,i,j) - j13g(k,i,j,i_uyz)*work_x(k,i,j)) * rfdz(k) &
501  ) * mapf(i,j,1,i_uy)
502  enddo
503  enddo
504  enddo
505  !$acc end kernels
506 #ifdef DEBUG
507  i = iundef; j = iundef; k = iundef
508 #endif
509 
510  ! 1/2 * dw/dy
511  ! (cell center; x,y,z)
512  !$omp parallel do
513  !$acc kernels
514  do j = jjs-1, jje+1
515  do i = iis-1, iie+1
516  do k = ks+1, ke-1
517 #ifdef DEBUG
518  call check( __line__, velz_c(k,i,j+1) )
519  call check( __line__, velz_c(k,i,j-1) )
520  call check( __line__, gsqrt(k,i,j+1,i_xyz) )
521  call check( __line__, gsqrt(k,i,j-1,i_xyz) )
522  call check( __line__, velz_xy(k,i,j) )
523  call check( __line__, velz_xy(k-1,i,j) )
524  call check( __line__, j23g(k,i,j,i_xyw) )
525  call check( __line__, j23g(k-1,i,j,i_xyw) )
526  call check( __line__, fdy(j) )
527  call check( __line__, fdy(j-1) )
528 #endif
529  s23_c(k,i,j) = 0.5_rp * ( &
530  ( gsqrt(k,i,j+1,i_xyz)*velz_c(k,i,j+1) - gsqrt(k,i,j-1,i_xyz)*velz_c(k,i,j-1) ) / ( fdy(j) + fdy(j-1) ) &
531  + ( j23g(k,i,j,i_xyw)*velz_xy(k,i,j) - j23g(k-1,i,j,i_xyw)*velz_xy(k-1,i,j) ) * rcdz(k) &
532  ) * mapf(i,j,2,i_xy)
533  enddo
534  enddo
535  enddo
536  !$acc end kernels
537 #ifdef DEBUG
538  i = iundef; j = iundef; k = iundef
539 #endif
540  !$omp parallel do
541  !$acc kernels
542  do j = jjs-1, jje+1
543  do i = iis-1, iie+1
544 #ifdef DEBUG
545  call check( __line__, velz_c(ks,i,j+1) )
546  call check( __line__, velz_c(ks,i,j-1) )
547  call check( __line__, gsqrt(ks,i,j+1,i_xyz) )
548  call check( __line__, gsqrt(ks,i,j-1,i_xyz) )
549  call check( __line__, velz_xy(ks,i,j) )
550  call check( __line__, j23g(ks,i,j,i_xyw) )
551  call check( __line__, velz_c(ke,i,j+1) )
552  call check( __line__, velz_c(ke,i,j-1) )
553  call check( __line__, gsqrt(ke,i,j+1,i_xyz) )
554  call check( __line__, gsqrt(ke,i,j-1,i_xyz) )
555  call check( __line__, velz_xy(ke,i,j) )
556  call check( __line__, j23g(ke,i,j,i_xyw) )
557  call check( __line__, fdy(j) )
558  call check( __line__, fdy(j-1) )
559 #endif
560  s23_c(ks,i,j) = 0.5_rp * ( &
561  ( gsqrt(ks,i,j+1,i_xyz)*velz_c(ks,i,j+1) - gsqrt(ks,i,j-1,i_xyz)*velz_c(ks,i,j-1) ) / ( fdy(j) + fdy(j-1) ) &
562  + ( j23g(ks,i,j,i_xyw)*velz_xy(ks,i,j) ) * rcdz(ks) &
563  ) * mapf(i,j,2,i_xy)
564  s23_c(ke,i,j) = 0.5_rp * ( &
565  ( gsqrt(ke,i,j+1,i_xyz)*velz_c(ke,i,j+1) - gsqrt(ke,i,j-1,i_xyz)*velz_c(ke,i,j-1) ) / ( fdy(j) + fdy(j-1) ) &
566  - ( j23g(ke-1,i,j,i_xyw)*velz_xy(ke-1,i,j) ) * rcdz(ke) &
567  ) * mapf(i,j,2,i_xy)
568  enddo
569  enddo
570  !$acc end kernels
571 #ifdef DEBUG
572  i = iundef; j = iundef; k = iundef
573 #endif
574 
575  ! (x edge; x,v,w)
576  !$omp parallel do
577  !$acc kernels
578  do j = jjs-1, jje
579  do i = iis , iie
580  do k = ks, ke-1
581 #ifdef DEBUG
582  call check( __line__, velz_xy(k,i,j+1) )
583  call check( __line__, velz_xy(k,i,j) )
584  call check( __line__, rfdy(j) )
585 #endif
586  s23_x(k,i,j) = 0.5_rp * ( &
587  ( gsqrt(k,i,j+1,i_xyw)*velz_xy(k,i,j+1) - gsqrt(k,i,j,i_xyw)*velz_xy(k,i,j) ) * rfdy(j) &
588  + ( j23g(k+1,i,j,i_xvz)*work_y(k+1,i,j) - j23g(k,i,j,i_xvz)*work_y(k,i,j) ) * rfdz(k) &
589  ) * mapf(i,j,2,i_xv)
590  enddo
591  enddo
592  enddo
593  !$acc end kernels
594 #ifdef DEBUG
595  i = iundef; j = iundef; k = iundef
596 #endif
597 
598 #ifdef DEBUG
599  work_z(:,:,:) = undef; work_x(:,:,:) = undef; work_y(:,:,:) = undef; work_v(:,:,:) = undef
600 #endif
601  ! u
602  ! (x-y plane; x,y,w)
603  !$omp parallel do
604  !$acc kernels
605  do j = jjs-1, jje+1
606  do i = iis-1, iie+1
607  do k = ks, ke-1
608 #ifdef DEBUG
609  call check( __line__, velx_c(k+1,i,j) )
610  call check( __line__, velx_c(k,i,j) )
611 #endif
612  work_z(k,i,j) = 0.5_rp * ( velx_c(k+1,i,j) + velx_c(k,i,j) )
613  enddo
614  enddo
615  enddo
616  !$acc end kernels
617 #ifdef DEBUG
618  i = iundef; j = iundef; k = iundef
619 #endif
620  ! (y-z plane; u,y,z)
621  ! WORK_X = VELX_YZ
622  ! (z-x plane; x,v,z)
623  !$omp parallel do
624  !$acc kernels
625  do j = jjs-1, jje
626  do i = iis-1, iie+1
627  do k = ks, ke
628 #ifdef DEBUG
629  call check( __line__, velx_c(k,i,j+1) )
630  call check( __line__, velx_c(k,i,j) )
631 #endif
632  work_y(k,i,j) = 0.5_rp * ( velx_c(k,i,j+1) + velx_c(k,i,j) )
633  enddo
634  enddo
635  enddo
636  !$acc end kernels
637 #ifdef DEBUG
638  i = iundef; j = iundef; k = iundef
639 #endif
640  ! (vertex; u,v,w)
641  !$omp parallel do
642  !$acc kernels
643  do j = jjs-1, jje
644  do i = iis-1, iie
645  do k = ks, ke-1
646 #ifdef DEBUG
647  call check( __line__, velx_yz(k,i,j) )
648  call check( __line__, velx_yz(k,i,j+1) )
649  call check( __line__, velx_yz(k+1,i,j) )
650  call check( __line__, velx_yz(k+1,i,j+1) )
651  call check( __line__, j23g(k ,i,j ,i_uvz) )
652  call check( __line__, j23g(k+1,i,j ,i_uvz) )
653  call check( __line__, j23g(k ,i,j+1,i_uvz) )
654  call check( __line__, j23g(k+1,i,j+1,i_uvz) )
655 #endif
656  work_v(k,i,j) = 0.25_rp &
657  * ( j23g(k ,i,j ,i_uyz)*velx_yz(k ,i,j ) &
658  + j23g(k+1,i,j ,i_uyz)*velx_yz(k+1,i,j ) &
659  + j23g(k ,i,j+1,i_uyz)*velx_yz(k ,i,j+1) &
660  + j23g(k+1,i,j+1,i_uyz)*velx_yz(k+1,i,j+1) )
661  enddo
662  enddo
663  enddo
664  !$acc end kernels
665 #ifdef DEBUG
666  i = iundef; j = iundef; k = iundef
667 #endif
668 
669  ! du/dx
670  ! (cell center; x,y,z)
671  !$omp parallel do
672  !$acc kernels
673  do j = jjs-1, jje+1
674  do i = iis-1, iie+1
675  do k = ks+1, ke-1
676 #ifdef DEBUG
677  call check( __line__, velx_yz(k,i,j) )
678  call check( __line__, velx_yz(k,i-1,j) )
679  call check( __line__, gsqrt(k,i,j,i_uyz) )
680  call check( __line__, gsqrt(k,i-1,j,i_uyz) )
681  call check( __line__, work_z(k,i,j) )
682  call check( __line__, work_z(k-1,i,j) )
683  call check( __line__, j13g(k,i,j,i_xyw) )
684  call check( __line__, j13g(k-1,i,j,i_xyw) )
685  call check( __line__, gsqrt(k,i,j,i_xyz) )
686  call check( __line__, rcdx(i) )
687 #endif
688  s11_c(k,i,j) = ( &
689  ( gsqrt(k,i,j,i_uyz)*velx_yz(k,i,j) - gsqrt(k,i-1,j,i_uyz)*velx_yz(k,i-1,j) ) * rcdx(i) &
690  + ( j13g(k,i,j,i_xyw)*work_z(k,i,j) - j13g(k-1,i,j,i_xyw)*work_z(k-1,i,j) ) * rcdz(k) &
691  ) * mapf(i,j,1,i_xy) / gsqrt(k,i,j,i_xyz)
692  enddo
693  enddo
694  enddo
695  !$acc end kernels
696 #ifdef DEBUG
697  i = iundef; j = iundef; k = iundef
698 #endif
699  !$omp parallel do
700  !$acc kernels
701  do j = jjs-1, jje+1
702  do i = iis-1, iie+1
703 #ifdef DEBUG
704  call check( __line__, velx_yz(ks,i,j) )
705  call check( __line__, velx_yz(ks,i-1,j) )
706  call check( __line__, gsqrt(ks,i,j,i_uyz) )
707  call check( __line__, gsqrt(ks,i-1,j,i_uyz) )
708  call check( __line__, velx_c(ks+1,i,j) )
709  call check( __line__, velx_c(ks,i,j) )
710  call check( __line__, j13g(ks+1,i,j,i_xyz) )
711  call check( __line__, j13g(ks,i,j,i_xyz) )
712  call check( __line__, gsqrt(ks,i,j,i_xyz) )
713  call check( __line__, velx_yz(ke,i,j) )
714  call check( __line__, velx_yz(ke,i-1,j) )
715  call check( __line__, gsqrt(ke,i,j,i_uyz) )
716  call check( __line__, gsqrt(ke,i-1,j,i_uyz) )
717  call check( __line__, velx_c(ke,i,j) )
718  call check( __line__, velx_c(ke-1,i,j) )
719  call check( __line__, j13g(ke,i,j,i_xyz) )
720  call check( __line__, j13g(ke-1,i,j,i_xyz) )
721  call check( __line__, gsqrt(ke,i,j,i_xyz) )
722  call check( __line__, rcdx(i) )
723 #endif
724  s11_c(ks,i,j) = ( &
725  ( gsqrt(ks,i,j,i_uyz)*velx_yz(ks,i,j) - gsqrt(ks,i-1,j,i_uyz)*velx_yz(ks,i-1,j) ) * rcdx(i) &
726  + ( j13g(ks+1,i,j,i_xyz)*velx_c(ks+1,i,j) - j13g(ks,i,j,i_xyz)*velx_c(ks,i,j) ) * rfdz(ks) &
727  ) * mapf(i,j,1,i_xy) / gsqrt(ks,i,j,i_xyz)
728  s11_c(ke,i,j) = ( &
729  ( gsqrt(ke,i,j,i_uyz)*velx_yz(ke,i,j) - gsqrt(ke,i-1,j,i_uyz)*velx_yz(ke,i-1,j) ) * rcdx(i) &
730  + ( j13g(ke,i,j,i_xyz)*velx_c(ke,i,j) - j13g(ke-1,i,j,i_xyz)*velx_c(ke-1,i,j) ) * rfdz(ke-1) &
731  ) * mapf(i,j,1,i_xy) / gsqrt(ke,i,j,i_xyz)
732  enddo
733  enddo
734  !$acc end kernels
735 #ifdef DEBUG
736  i = iundef; j = iundef; k = iundef
737 #endif
738 
739  ! 1/2 * du/dz
740  ! (cell center; x,y,z)
741  !$omp parallel do
742  !$acc kernels
743  do j = jjs-1, jje+1
744  do i = iis-1, iie+1
745  do k = ks+1, ke-1
746 #ifdef DEBUG
747  call check( __line__, s31_c(k,i,j) )
748  call check( __line__, velx_c(k+1,i,j) )
749  call check( __line__, velx_c(k-1,i,j) )
750  call check( __line__, fdz(k) )
751  call check( __line__, fdz(k-1) )
752 #endif
753  s31_c(k,i,j) = ( s31_c(k,i,j) & ! dw/dx
754  + 0.5_rp * ( velx_c(k+1,i,j) - velx_c(k-1,i,j) ) * j33g / ( fdz(k) + fdz(k-1) ) &
755  ) / gsqrt(k,i,j,i_xyz)
756  enddo
757  enddo
758  enddo
759  !$acc end kernels
760 #ifdef DEBUG
761  i = iundef; j = iundef; k = iundef
762 #endif
763  !$omp parallel do
764  !$acc kernels
765  do j = jjs-1, jje+1
766  do i = iis-1, iie+1
767 #ifdef DEBUG
768  call check( __line__, s31_c(ks,i,j) )
769  call check( __line__, velx_c(ks+1,i,j) )
770  call check( __line__, velx_c(ks,i,j) )
771  call check( __line__, rfdz(ks) )
772  call check( __line__, s31_c(ke,i,j) )
773  call check( __line__, velx_c(ke,i,j) )
774  call check( __line__, velx_c(ke-1,i,j) )
775  call check( __line__, rfdz(ke-1) )
776 #endif
777  s31_c(ks,i,j) = ( s31_c(ks,i,j) &
778  + 0.5_rp * ( velx_c(ks+1,i,j) - velx_c(ks,i,j) ) * j33g * rfdz(ks) &
779  ) / gsqrt(ks,i,j,i_xyz)
780  s31_c(ke,i,j) = ( s31_c(ke,i,j) &
781  + 0.5_rp * ( velx_c(ke,i,j) - velx_c(ke-1,i,j) ) * j33g * rfdz(ke-1) &
782  ) / gsqrt(ke,i,j,i_xyz)
783  enddo
784  enddo
785  !$acc end kernels
786 #ifdef DEBUG
787  i = iundef; j = iundef; k = iundef
788 #endif
789  ! (y edge; u,y,w)
790  !$omp parallel do
791  !$acc kernels
792  do j = jjs , jje
793  do i = iis-1, iie
794  do k = ks, ke-1
795 #ifdef DEBUG
796  call check( __line__, s31_y(k,i,j) )
797  call check( __line__, velx_yz(k+1,i,j) )
798  call check( __line__, velx_yz(k,i,j) )
799  call check( __line__, rfdz(k) )
800 #endif
801  s31_y(k,i,j) = ( s31_y(k,i,j) & ! dw/dx
802  + 0.5_rp * ( velx_yz(k+1,i,j) - velx_yz(k,i,j) ) * j33g * rfdz(k) &
803  ) / gsqrt(k,i,j,i_uyw)
804  enddo
805  enddo
806  enddo
807  !$acc end kernels
808 #ifdef DEBUG
809  i = iundef; j = iundef; k = iundef
810 #endif
811 
812  ! 1/2 * du/dy
813  ! (cell center; x,y,z)
814  !$omp parallel do
815  !$acc kernels
816  do j = jjs-1, jje+1
817  do i = iis-1, iie+1
818  do k = ks+1, ke-1
819 #ifdef DEBUG
820  call check( __line__, velx_c(k,i,j+1) )
821  call check( __line__, velx_c(k,i,j-1) )
822  call check( __line__, gsqrt(k,i,j+1,i_xyz) )
823  call check( __line__, gsqrt(k,i,j-1,i_xyz) )
824  call check( __line__, work_z(k,i,j) )
825  call check( __line__, work_z(k-1,i,j) )
826  call check( __line__, j23g(k,i,j,i_xyw) )
827  call check( __line__, j23g(k-1,i,j,i_xyw) )
828 #endif
829  s12_c(k,i,j) = 0.5_rp * ( &
830  ( gsqrt(k,i,j+1,i_xyz)*velx_c(k,i,j+1) - gsqrt(k,i,j-1,i_xyz)*velx_c(k,i,j-1) ) / ( fdy(j) + fdy(j-1) ) &
831  + ( j23g(k,i,j,i_xyw)*work_z(k,i,j) - j23g(k-1,i,j,i_xyw)*work_z(k-1,i,j) ) * rcdz(k) &
832  ) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz)
833  enddo
834  enddo
835  enddo
836  !$acc end kernels
837 #ifdef DEBUG
838  i = iundef; j = iundef; k = iundef
839 #endif
840  !$omp parallel do
841  !$acc kernels
842  do j = jjs-1, jje+1
843  do i = iis-1, iie+1
844 #ifdef DEBUG
845  call check( __line__, velx_c(ks,i,j+1) )
846  call check( __line__, velx_c(ks,i,j-1) )
847  call check( __line__, gsqrt(ks,i,j+1,i_xyz) )
848  call check( __line__, gsqrt(ks,i,j-1,i_xyz) )
849  call check( __line__, velx_c(ks+1,i,j) )
850  call check( __line__, velx_c(ks,i,j) )
851  call check( __line__, j23g(ks+1,i,j,i_xyz) )
852  call check( __line__, j23g(ks,i,j,i_xyz) )
853  call check( __line__, gsqrt(ks,i,j,i_xyz) )
854  call check( __line__, velx_c(ke,i,j+1) )
855  call check( __line__, velx_c(ke,i,j-1) )
856  call check( __line__, gsqrt(ke,i,j+1,i_xyz) )
857  call check( __line__, gsqrt(ke,i,j-1,i_xyz) )
858  call check( __line__, velx_c(ke,i,j) )
859  call check( __line__, velx_c(ke-1,i,j) )
860  call check( __line__, j23g(ke,i,j,i_xyz) )
861  call check( __line__, j23g(ke-1,i,j,i_xyz) )
862  call check( __line__, gsqrt(ke,i,j,i_xyz) )
863  call check( __line__, fdy(j) )
864  call check( __line__, fdy(j-1) )
865 #endif
866  s12_c(ks,i,j) = 0.5_rp * ( &
867  ( gsqrt(ks,i,j+1,i_xyz)*velx_c(ks,i,j+1) - gsqrt(ks,i,j-1,i_xyz)*velx_c(ks,i,j-1) ) / ( fdy(j) + fdy(j-1) ) &
868  + ( j23g(ks+1,i,j,i_xyz)*velx_c(ks+1,i,j) - j23g(ks,i,j,i_xyz)*velx_c(ks,i,j) ) * rfdz(ks) &
869  ) * mapf(i,j,2,i_xy) / gsqrt(ks,i,j,i_xyz)
870  s12_c(ke,i,j) = 0.5_rp * ( &
871  ( gsqrt(ke,i,j+1,i_xyz)*velx_c(ke,i,j+1) - gsqrt(ke,i,j-1,i_xyz)*velx_c(ke,i,j-1) ) / ( fdy(j) + fdy(j-1) ) &
872  + ( j23g(ke,i,j,i_xyz)*velx_c(ke,i,j) - j23g(ke-1,i,j,i_xyz)*velx_c(ke-1,i,j) ) * rfdz(ke-1) &
873  ) * mapf(i,j,2,i_xy) / gsqrt(ke,i,j,i_xyz)
874  enddo
875  enddo
876  !$acc end kernels
877 #ifdef DEBUG
878  i = iundef; j = iundef; k = iundef
879 #endif
880 
881  ! (z edge; u,v,z)
882  !$omp parallel do
883  !$acc kernels
884  do j = jjs-1, jje
885  do i = iis-1, iie
886  do k = ks+1, ke-1
887 #ifdef DEBUG
888  call check( __line__, velx_yz(k,i,j+1) )
889  call check( __line__, velx_yz(k,i,j) )
890  call check( __line__, work_v(k,i,j) )
891  call check( __line__, work_v(k-1,i,j) )
892  call check( __line__, rfdy(j) )
893 #endif
894  s12_z(k,i,j) = 0.5_rp * ( &
895  ( gsqrt(k,i,j+1,i_uyz)*velx_yz(k,i,j+1) - gsqrt(k,i,j,i_uyz)*velx_yz(k,i,j) ) * rfdy(j) &
896  + ( work_v(k,i,j) - work_v(k-1,i,j) ) * rcdz(k) &
897  ) * mapf(i,j,2,i_uv)
898  enddo
899  enddo
900  enddo
901  !$acc end kernels
902 #ifdef DEBUG
903  i = iundef; j = iundef; k = iundef
904 #endif
905  !$omp parallel do
906  !$acc kernels
907  do j = jjs-1, jje
908  do i = iis-1, iie
909 #ifdef DEBUG
910  call check( __line__, velx_yz(ks,i,j+1) )
911  call check( __line__, velx_yz(ks,i,j) )
912  call check( __line__, velx_yz(ks+1,i,j) )
913  call check( __line__, velx_yz(ks+1,i,j+1) )
914  call check( __line__, j23g(ks+1,i,j,i_uvz) )
915  call check( __line__, j23g(ks ,i,j,i_uvz) )
916  call check( __line__, velx_yz(ke,i,j+1) )
917  call check( __line__, velx_yz(ke,i,j) )
918  call check( __line__, velx_yz(ke-1,i,j) )
919  call check( __line__, velx_yz(ke-1,i,j+1) )
920  call check( __line__, j23g(ke ,i,j,i_uvz) )
921  call check( __line__, j23g(ke-1,i,j,i_uvz) )
922 #endif
923  s12_z(ks,i,j) = 0.25_rp * ( &
924  ( gsqrt(ks,i,j+1,i_uyz)*velx_yz(ks,i,j+1) - gsqrt(ks,i,j,i_uyz)*velx_yz(ks,i,j) ) * rfdy(j) &
925  + ( j23g(ks+1,i,j,i_uvz) * ( velx_yz(ks+1,i,j) + velx_yz(ks+1,i,j+1) ) &
926  - j23g(ks ,i,j,i_uvz) * ( velx_yz(ks ,i,j) + velx_yz(ks ,i,j+1) ) ) * rfdz(ks) &
927  ) * mapf(i,j,2,i_uv)
928  s12_z(ke,i,j) = 0.25_rp * ( &
929  ( gsqrt(ke,i,j+1,i_uyz)*velx_yz(ke,i,j+1) - gsqrt(ke,i,j,i_uyz)*velx_yz(ke,i,j) ) * rfdy(j) &
930  + ( j23g(ke ,i,j,i_uvz) * ( velx_yz(ke ,i,j) + velx_yz(ke ,i,j+1) ) &
931  - j23g(ke-1,i,j,i_uvz) * ( velx_yz(ke-1,i,j) + velx_yz(ke-1,i,j+1) ) ) * rfdz(ke-1) &
932  ) * mapf(i,j,2,i_uv)
933  enddo
934  enddo
935  !$acc end kernels
936 #ifdef DEBUG
937  i = iundef; j = iundef; k = iundef
938 #endif
939 
940 #ifdef DEBUG
941  work_z(:,:,:) = undef; work_x(:,:,:) = undef; work_y(:,:,:) = undef; work_v(:,:,:) = undef
942 #endif
943  ! v
944  ! (x-y plane; x,y,w)
945  !$omp parallel do
946  !$acc kernels
947  do j = jjs-1, jje+1
948  do i = iis-1, iie+1
949  do k = ks, ke-1
950 #ifdef DEBUG
951  call check( __line__, vely_c(k+1,i,j) )
952  call check( __line__, vely_c(k,i,j) )
953 #endif
954  work_z(k,i,j) = 0.5_rp * ( vely_c(k+1,i,j) + vely_c(k,i,j) )
955  enddo
956  enddo
957  enddo
958  !$acc end kernels
959 #ifdef DEBUG
960  i = iundef; j = iundef; k = iundef
961 #endif
962  ! (y-z plane; u,y,z)
963  !$omp parallel do
964  !$acc kernels
965  do j = jjs-1, jje+1
966  do i = iis-1, iie
967  do k = ks, ke
968 #ifdef DEBUG
969  call check( __line__, vely_c(k,i+1,j) )
970  call check( __line__, vely_c(k,i,j) )
971 #endif
972  work_x(k,i,j) = 0.5_rp * ( velz_c(k,i+1,j) + velz_c(k,i,j) )
973  enddo
974  enddo
975  enddo
976  !$acc end kernels
977 #ifdef DEBUG
978  i = iundef; j = iundef; k = iundef
979 #endif
980  ! (z-x plane; x,v,z)
981  ! WORK_Y = VELY_ZX
982  ! (vertex; u,v,w)
983  !$omp parallel do
984  !$acc kernels
985  do j = jjs-1, jje
986  do i = iis-1, iie
987  do k = ks, ke-1
988 #ifdef DEBUG
989  call check( __line__, vely_zx(k,i,j) )
990  call check( __line__, vely_zx(k+1,i,j) )
991  call check( __line__, vely_zx(k,i+1,j) )
992  call check( __line__, vely_zx(k+1,i+1,j) )
993 #endif
994  work_v(k,i,j) = 0.25_rp &
995  * ( j13g(k ,i ,j,i_xvz)*vely_zx(k ,i ,j) &
996  + j13g(k+1,i ,j,i_xvz)*vely_zx(k+1,i ,j) &
997  + j13g(k ,i+1,j,i_xvz)*vely_zx(k ,i+1,j) &
998  + j13g(k+1,i+1,j,i_xvz)*vely_zx(k+1,i+1,j) )
999  enddo
1000  enddo
1001  enddo
1002  !$acc end kernels
1003 #ifdef DEBUG
1004  i = iundef; j = iundef; k = iundef
1005 #endif
1006 
1007  ! dv/dy
1008  ! (cell center; x,y,z)
1009  !$omp parallel do
1010  !$acc kernels
1011  do j = jjs-1, jje+1
1012  do i = iis-1, iie+1
1013  do k = ks+1, ke-1
1014 #ifdef DEBUG
1015  call check( __line__, vely_zx(k,i,j) )
1016  call check( __line__, vely_zx(k,i,j-1) )
1017  call check( __line__, gsqrt(k,i,j,i_xvz) )
1018  call check( __line__, gsqrt(k,i,j-1,i_xvz) )
1019  call check( __line__, work_z(k,i,j) )
1020  call check( __line__, work_z(k-1,i,j) )
1021  call check( __line__, j23g(k,i,j,i_xyw) )
1022  call check( __line__, j23g(k-1,i,j,i_xyw) )
1023  call check( __line__, rcdy(j) )
1024 #endif
1025  s22_c(k,i,j) = ( &
1026  ( gsqrt(k,i,j,i_xvz)*vely_zx(k,i,j) - gsqrt(k,i,j-1,i_xvz)*vely_zx(k,i,j-1) ) * rcdy(j) &
1027  + ( j23g(k,i,j,i_xyw)*work_z(k,i,j) - j23g(k-1,i,j,i_xyw)*work_z(k-1,i,j) ) * rcdz(k) &
1028  ) * mapf(i,j,2,i_xy) / gsqrt(k,i,j,i_xyz)
1029  enddo
1030  enddo
1031  enddo
1032  !$acc end kernels
1033 #ifdef DEBUG
1034  i = iundef; j = iundef; k = iundef
1035 #endif
1036  !$omp parallel do
1037  !$acc kernels
1038  do j = jjs-1, jje+1
1039  do i = iis-1, iie+1
1040 #ifdef DEBUG
1041  call check( __line__, vely_zx(ks,i,j) )
1042  call check( __line__, vely_zx(ks,i,j-1) )
1043  call check( __line__, gsqrt(ks,i,j,i_xvz) )
1044  call check( __line__, gsqrt(ks,i,j-1,i_xvz) )
1045  call check( __line__, vely_c(ks+1,i,j) )
1046  call check( __line__, vely_c(ks,i,j) )
1047  call check( __line__, j23g(ks+1,i,j,i_xyz) )
1048  call check( __line__, j23g(ks,i,j,i_xyz) )
1049  call check( __line__, rcdy(j) )
1050  call check( __line__, vely_zx(ke,i,j) )
1051  call check( __line__, vely_zx(ke,i,j-1) )
1052  call check( __line__, gsqrt(ke,i,j,i_xvz) )
1053  call check( __line__, gsqrt(ke,i,j-1,i_xvz) )
1054  call check( __line__, vely_c(ke,i,j) )
1055  call check( __line__, vely_c(ke-1,i,j) )
1056  call check( __line__, j23g(ke,i,j,i_xyz) )
1057  call check( __line__, j23g(ke-1,i,j,i_xyz) )
1058 #endif
1059  s22_c(ks,i,j) = ( &
1060  ( gsqrt(ks,i,j,i_xvz)*vely_zx(ks,i,j) - gsqrt(ks,i,j-1,i_xvz)*vely_zx(ks,i,j-1) ) * rcdy(j) &
1061  + ( j23g(ks+1,i,j,i_xyz)*vely_c(ks+1,i,j) - j23g(ks,i,j,i_xyz)*vely_c(ks,i,j) ) * rfdz(ks) &
1062  ) * mapf(i,j,2,i_xy) / gsqrt(ks,i,j,i_xyz)
1063  s22_c(ke,i,j) = ( &
1064  ( gsqrt(ke,i,j,i_xvz)*vely_zx(ke,i,j) - gsqrt(ke,i,j-1,i_xvz)*vely_zx(ke,i,j-1) ) * rcdy(j) &
1065  + ( j23g(ke,i,j,i_xyz)*vely_c(ke,i,j) - j23g(ke-1,i,j,i_xyz)*vely_c(ke-1,i,j) ) * rfdz(ke-1) &
1066  ) * mapf(i,j,2,i_xy) / gsqrt(ke,i,j,i_xyz)
1067  enddo
1068  enddo
1069  !$acc end kernels
1070 #ifdef DEBUG
1071  i = iundef; j = iundef; k = iundef
1072 #endif
1073 
1074  ! 1/2 * dv/dx
1075  ! (cell center; x,y,z)
1076  !$omp parallel do
1077  !$acc kernels
1078  do j = jjs-1, jje+1
1079  do i = iis-1, iie+1
1080  do k = ks+1, ke-1
1081 #ifdef DEBUG
1082  call check( __line__, s12_c(k,i,j) )
1083  call check( __line__, vely_c(k,i+1,j) )
1084  call check( __line__, vely_c(k,i-1,j) )
1085  call check( __line__, gsqrt(k,i+1,j,i_xyz) )
1086  call check( __line__, gsqrt(k,i-1,j,i_xyz) )
1087  call check( __line__, work_z(k,i,j) )
1088  call check( __line__, work_z(k-1,i,j) )
1089  call check( __line__, j13g(k,i,j,i_xyw) )
1090  call check( __line__, j13g(k-1,i,j,i_xyw) )
1091  call check( __line__, gsqrt(k,i,j,i_xyz) )
1092  call check( __line__, fdx(i) )
1093  call check( __line__, fdx(i-1) )
1094 #endif
1095  s12_c(k,i,j) = ( s12_c(k,i,j) & ! du/dy
1096  + 0.5_rp * ( &
1097  ( gsqrt(k,i+1,j,i_xyz)*vely_c(k,i+1,j) - gsqrt(k,i-1,j,i_xyz)*vely_c(k,i-1,j) ) / ( fdx(i) + fdx(i-1) ) &
1098  + ( j13g(k,i,j,i_xyw)*work_z(k,i,j) - j13g(k-1,i,j,i_xyw)*work_z(k-1,i,j) ) * rcdz(k) ) * mapf(i,j,1,i_xy) &
1099  ) / gsqrt(k,i,j,i_xyz)
1100  enddo
1101  enddo
1102  enddo
1103  !$acc end kernels
1104 #ifdef DEBUG
1105  i = iundef; j = iundef; k = iundef
1106 #endif
1107  !$omp parallel do
1108  !$acc kernels
1109  do j = jjs-1, jje+1
1110  do i = iis-1, iie+1
1111 #ifdef DEBUG
1112  call check( __line__, s12_c(ks,i,j) )
1113  call check( __line__, vely_c(ks,i+1,j) )
1114  call check( __line__, vely_c(ks,i-1,j) )
1115  call check( __line__, gsqrt(ks,i+1,j,i_xyz) )
1116  call check( __line__, gsqrt(ks,i-1,j,i_xyz) )
1117  call check( __line__, vely_c(ks+1,i,j) )
1118  call check( __line__, vely_c(ks,i,j) )
1119  call check( __line__, j13g(ks+1,i,j,i_xyz) )
1120  call check( __line__, j13g(ks,i,j,i_xyz) )
1121  call check( __line__, gsqrt(ks,i,j,i_xyz) )
1122  call check( __line__, s12_c(ke,i,j) )
1123  call check( __line__, vely_c(ke,i+1,j) )
1124  call check( __line__, vely_c(ke,i-1,j) )
1125  call check( __line__, gsqrt(ke,i+1,j,i_xyz) )
1126  call check( __line__, gsqrt(ke,i-1,j,i_xyz) )
1127  call check( __line__, vely_c(ke,i,j) )
1128  call check( __line__, vely_c(ke-1,i,j) )
1129  call check( __line__, j13g(ke,i,j,i_xyz) )
1130  call check( __line__, j13g(ke-1,i,j,i_xyz) )
1131  call check( __line__, gsqrt(ke,i,j,i_xyz) )
1132  call check( __line__, fdx(i) )
1133  call check( __line__, fdx(i-1) )
1134 #endif
1135  s12_c(ks,i,j) = ( s12_c(ks,i,j) & ! du/dy
1136  + 0.5_rp * ( &
1137  ( gsqrt(ks,i+1,j,i_xyz)*vely_c(ks,i+1,j) - gsqrt(ks,i-1,j,i_xyz)*vely_c(ks,i-1,j) ) / ( fdx(i) + fdx(i-1) ) &
1138  + ( j13g(ks+1,i,j,i_xyz)*vely_c(ks+1,i,j) - j13g(ks,i,j,i_xyz)*vely_c(ks,i,j) ) * rfdz(ks) ) &
1139  * mapf(i,j,1,i_xy) &
1140  ) / gsqrt(ks,i,j,i_xyz)
1141  s12_c(ke,i,j) = ( s12_c(ke,i,j) & ! du/dy
1142  + 0.5_rp * ( &
1143  ( gsqrt(ke,i+1,j,i_xyz)*vely_c(ke,i+1,j) - gsqrt(ke,i-1,j,i_xyz)*vely_c(ke,i-1,j) ) / ( fdx(i) + fdx(i-1) ) &
1144  + ( j13g(ke,i,j,i_xyz)*vely_c(ke,i,j) - j13g(ke-1,i,j,i_xyz)*vely_c(ke-1,i,j) ) * rfdz(ke-1) ) &
1145  * mapf(i,j,1,i_xy) &
1146  ) / gsqrt(ke,i,j,i_xyz)
1147  enddo
1148  enddo
1149  !$acc end kernels
1150 #ifdef DEBUG
1151  i = iundef; j = iundef; k = iundef
1152 #endif
1153  ! (z edge; u,v,z)
1154  !$omp parallel do
1155  !$acc kernels
1156  do j = jjs-1, jje
1157  do i = iis-1, iie
1158  do k = ks+1, ke-1
1159 #ifdef DEBUG
1160  call check( __line__, s12_z(k,i,j) )
1161  call check( __line__, vely_zx(k,i+1,j) )
1162  call check( __line__, vely_zx(k,i,j) )
1163  call check( __line__, work_v(k,i,j) )
1164  call check( __line__, work_v(k-1,i,j) )
1165  call check( __line__, rfdx(i) )
1166 #endif
1167  s12_z(k,i,j) = ( s12_z(k,i,j) &
1168  + 0.5_rp * ( &
1169  ( gsqrt(k,i+1,j,i_xvz)*vely_zx(k,i+1,j) - gsqrt(k,i,j,i_xvz)*vely_zx(k,i,j) ) * rfdx(i) &
1170  + ( work_v(k,i,j) - work_v(k-1,i,j) ) * rcdz(k) ) * mapf(i,j,1,i_uv) &
1171  ) / gsqrt(k,i,j,i_uvz)
1172  enddo
1173  enddo
1174  enddo
1175  !$acc end kernels
1176 #ifdef DEBUG
1177  i = iundef; j = iundef; k = iundef
1178 #endif
1179  !$omp parallel do
1180  !$acc kernels
1181  do j = jjs-1, jje
1182  do i = iis-1, iie
1183 #ifdef DEBUG
1184  call check( __line__, s12_z(ks,i,j) )
1185  call check( __line__, vely_zx(ks,i+1,j) )
1186  call check( __line__, vely_zx(ks,i,j) )
1187  call check( __line__, vely_zx(ks+1,i,j) )
1188  call check( __line__, vely_zx(ks+1,i+1,j) )
1189  call check( __line__, s12_z(ke,i,j) )
1190  call check( __line__, vely_zx(ke,i+1,j) )
1191  call check( __line__, vely_zx(ke,i,j) )
1192  call check( __line__, vely_zx(ke-1,i,j) )
1193  call check( __line__, vely_zx(ke-1,i+1,j) )
1194  call check( __line__, rfdx(i) )
1195 #endif
1196  s12_z(ks,i,j) = ( s12_z(ks,i,j) &
1197  + 0.5_rp * ( &
1198  ( gsqrt(ks,i+1,j,i_xvz)*vely_zx(ks,i+1,j) - gsqrt(ks,i,j,i_xvz)*vely_zx(ks,i,j) ) * rfdx(i) &
1199  + ( j13g(ks+1,i,j,i_uvz) * ( vely_zx(ks+1,i,j) + vely_zx(ks+1,i+1,j) ) &
1200  - j13g(ks ,i,j,i_uvz) * ( vely_zx(ks ,i,j) + vely_zx(ks ,i+1,j) ) ) * rfdz(ks) ) * mapf(i,j,1,i_uv) &
1201  ) / gsqrt(ks,i,j,i_uvz)
1202  s12_z(ke,i,j) = ( s12_z(ke,i,j) &
1203  + 0.5_rp * ( &
1204  ( gsqrt(ke,i+1,j,i_xvz)*vely_zx(ke,i+1,j) - gsqrt(ke,i,j,i_xvz)*vely_zx(ke,i,j) ) * rfdx(i) &
1205  + ( j13g(ke ,i,j,i_uvz) * ( vely_zx(ke ,i,j) + vely_zx(ke ,i+1,j) ) &
1206  - j13g(ke-1,i,j,i_uvz) * ( vely_zx(ke-1,i,j) + vely_zx(ke-1,i+1,j) ) ) * rfdz(ke-1) ) * mapf(i,j,1,i_uv) &
1207  ) / gsqrt(ke,i,j,i_uvz)
1208  enddo
1209  enddo
1210  !$acc end kernels
1211 #ifdef DEBUG
1212  i = iundef; j = iundef; k = iundef
1213 #endif
1214 
1215  ! 1/2 * dv/dz
1216  ! (cell center; x,y,z)
1217  !$omp parallel do
1218  !$acc kernels
1219  do j = jjs-1, jje+1
1220  do i = iis-1, iie+1
1221  do k = ks+1, ke-1
1222 #ifdef DEBUG
1223  call check( __line__, s23_c(k,i,j) )
1224  call check( __line__, vely_c(k+1,i,j) )
1225  call check( __line__, vely_c(k-1,i,j) )
1226  call check( __line__, fdz(k) )
1227  call check( __line__, fdz(k-1) )
1228 #endif
1229  s23_c(k,i,j) = ( s23_c(k,i,j) & ! dw/dy
1230  + 0.5_rp * ( vely_c(k+1,i,j) - vely_c(k-1,i,j) ) * j33g / ( fdz(k) + fdz(k-1) ) &
1231  ) / gsqrt(k,i,j,i_xyz)
1232  enddo
1233  enddo
1234  enddo
1235  !$acc end kernels
1236 #ifdef DEBUG
1237  i = iundef; j = iundef; k = iundef
1238 #endif
1239  !$omp parallel do
1240  !$acc kernels
1241  do j = jjs-1, jje+1
1242  do i = iis-1, iie+1
1243 #ifdef DEBUG
1244  call check( __line__, s23_c(ks,i,j) )
1245  call check( __line__, vely_c(ks+1,i,j) )
1246  call check( __line__, vely_c(ks,i,j) )
1247  call check( __line__, rfdz(ks) )
1248  call check( __line__, s23_c(ke,i,j) )
1249  call check( __line__, vely_c(ke,i,j) )
1250  call check( __line__, vely_c(ke-1,i,j) )
1251  call check( __line__, rfdz(ke-1) )
1252 #endif
1253  s23_c(ks,i,j) = ( s23_c(ks,i,j) &
1254  + 0.5_rp * ( vely_c(ks+1,i,j) - vely_c(ks,i,j) ) * j33g * rfdz(ks) &
1255  ) / gsqrt(ks,i,j,i_xyz)
1256  s23_c(ke,i,j) = ( s23_c(ke,i,j) &
1257  + 0.5_rp * ( vely_c(ke,i,j) - vely_c(ke-1,i,j) ) * j33g * rfdz(ke-1) &
1258  ) / gsqrt(ke,i,j,i_xyz)
1259  enddo
1260  enddo
1261  !$acc end kernels
1262 #ifdef DEBUG
1263  i = iundef; j = iundef; k = iundef
1264 #endif
1265 
1266  ! (x edge; x,v,w)
1267  !$omp parallel do
1268  !$acc kernels
1269  do j = jjs-1, jje
1270  do i = iis , iie
1271  do k = ks, ke-1
1272 #ifdef DEBUG
1273  call check( __line__, s23_x(k,i,j) )
1274  call check( __line__, vely_zx(k+1,i,j) )
1275  call check( __line__, vely_zx(k,i,j) )
1276  call check( __line__, rfdz(k) )
1277 #endif
1278  s23_x(k,i,j) = ( s23_x(k,i,j) &
1279  + 0.5_rp * ( vely_zx(k+1,i,j) - vely_zx(k,i,j) ) * j33g * rfdz(k) &
1280  ) / gsqrt(k,i,j,i_xvw)
1281  enddo
1282  enddo
1283  enddo
1284  !$acc end kernels
1285 #ifdef DEBUG
1286  i = iundef; j = iundef; k = iundef
1287 #endif
1288 
1289 
1290  ! |S|^2 = 2*Sij*Sij
1291 #ifdef DEBUG
1292  work_z(:,:,:) = undef; work_x(:,:,:) = undef; work_y(:,:,:) = undef
1293 #endif
1294  ! (cell center)
1295  !$omp parallel do default(none) &
1296  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,S11_C,S22_C,S33_C,S31_C,S12_C,S23_C,S2) &
1297  !$omp private(i,j,k) OMP_SCHEDULE_ collapse(2)
1298  !$acc kernels
1299  do j = jjs-1, jje+1
1300  do i = iis-1, iie+1
1301  do k = ks, ke
1302 #ifdef DEBUG
1303  call check( __line__, s11_c(k,i,j) )
1304  call check( __line__, s22_c(k,i,j) )
1305  call check( __line__, s33_c(k,i,j) )
1306  call check( __line__, s31_c(k,i,j) )
1307  call check( __line__, s12_c(k,i,j) )
1308  call check( __line__, s23_c(k,i,j) )
1309 #endif
1310  s2(k,i,j) = 2.0_rp * ( s11_c(k,i,j)**2 + s22_c(k,i,j)**2 + s33_c(k,i,j)**2 ) &
1311  + 4.0_rp * ( s31_c(k,i,j)**2 + s12_c(k,i,j)**2 + s23_c(k,i,j)**2 )
1312  enddo
1313  enddo
1314  enddo
1315  !$acc end kernels
1316 #ifdef DEBUG
1317  i = iundef; j = iundef; k = iundef
1318 #endif
1319 
1320  enddo
1321  enddo
1322 
1323  !$acc end data
1324 
1325  return
1326  end subroutine atmos_phy_tb_calc_strain_tensor
1327 
1328  !-----------------------------------------------------------------------------
1329  subroutine atmos_phy_tb_calc_flux_phi( &
1330  qflx_phi, &
1331  DENS, PHI, Kh, FACT, &
1332  GSQRT, J13G, J23G, J33G, MAPF, &
1333  horizontal, &
1334  implicit, &
1335  a, b, c, dt, &
1336  IIS, IIE, JJS, JJE )
1338  rcdz => atmos_grid_cartesc_rcdz, &
1339  rfdz => atmos_grid_cartesc_rfdz, &
1340  rfdx => atmos_grid_cartesc_rfdx, &
1341  rfdy => atmos_grid_cartesc_rfdy
1342  use scale_matrix, only: &
1343  matrix_solver_tridiagonal
1344  implicit none
1345 
1346  real(rp), intent(inout) :: qflx_phi(ka,ia,ja,3)
1347  real(rp), intent(in) :: dens (ka,ia,ja)
1348  real(rp), intent(in) :: phi (ka,ia,ja)
1349  real(rp), intent(in) :: kh (ka,ia,ja)
1350  real(rp), intent(in) :: fact
1351  real(rp), intent(in) :: gsqrt (ka,ia,ja,7)
1352  real(rp), intent(in) :: j13g (ka,ia,ja,7)
1353  real(rp), intent(in) :: j23g (ka,ia,ja,7)
1354  real(rp), intent(in) :: j33g
1355  real(rp), intent(in) :: mapf (ia,ja,2,4)
1356  logical, intent(in) :: horizontal
1357  logical, intent(in) :: implicit
1358  real(rp), intent(in) :: a (ka,ia,ja)
1359  real(rp), intent(in) :: b (ka,ia,ja)
1360  real(rp), intent(in) :: c (ka,ia,ja)
1361  real(dp), intent(in) :: dt
1362  integer, intent(in) :: iis
1363  integer, intent(in) :: iie
1364  integer, intent(in) :: jjs
1365  integer, intent(in) :: jje
1366 
1367  real(rp) :: tend(ka,ia,ja)
1368  real(rp) :: tend2(ka,ia,ja)
1369  real(rp) :: d(ka)
1370 
1371  integer :: k, i, j
1372 
1373  !$acc data copy(qflx_phi) &
1374  !$acc copyin(DENS, PHI, Kh, GSQRT, J13G, J23G, MAPF, a, b, c, &
1375  !$acc RCDZ, RFDZ, RFDX, RFDY) &
1376  !$acc create(TEND,TEND2)
1377 
1378 
1379  ! (x-y plane; x,y,w)
1380  if ( horizontal ) then
1381  !$omp parallel workshare
1382  !$acc kernels
1383 !XFILL
1384  qflx_phi(:,:,:,zdir) = 0.0_rp
1385  !$acc end kernels
1386  !$omp end parallel workshare
1387  else
1388  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
1389  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,DENS,Kh,FACT,PHI,qflx_phi,GSQRT,I_XYW,RFDZ,J33G)
1390  !$acc kernels
1391  do j = jjs, jje
1392  do i = iis, iie
1393  do k = ks, ke-1
1394 #ifdef DEBUG
1395  call check( __line__, dens(k,i,j) )
1396  call check( __line__, dens(k+1,i,j) )
1397  call check( __line__, kh(k,i,j) )
1398  call check( __line__, kh(k+1,i,j) )
1399  call check( __line__, phi(k+1,i,j) )
1400  call check( __line__, phi(k,i,j) )
1401  call check( __line__, rfdz(k) )
1402 #endif
1403  qflx_phi(k,i,j,zdir) = - 0.25_rp & ! 1/2/2
1404  * ( dens(k,i,j)+dens(k+1,i,j) ) &
1405  * ( kh(k,i,j) + kh(k+1,i,j) ) * fact &
1406  * ( phi(k+1,i,j)-phi(k,i,j) ) * rfdz(k) * j33g &
1407  / gsqrt(k,i,j,i_xyw)
1408  enddo
1409  enddo
1410  enddo
1411  !$acc end kernels
1412 #ifdef DEBUG
1413  i = iundef; j = iundef; k = iundef
1414 #endif
1415  !$omp parallel do
1416  !$acc kernels
1417  do j = jjs, jje
1418  do i = iis, iie
1419  qflx_phi(ks-1,i,j,zdir) = 0.0_rp
1420  qflx_phi(ke ,i,j,zdir) = 0.0_rp
1421  enddo
1422  enddo
1423  !$acc end kernels
1424 #ifdef DEBUG
1425  i = iundef; j = iundef; k = iundef
1426 #endif
1427  end if
1428 
1429  ! (y-z plane; u,y,z)
1430  !$omp parallel do default(none) private(i,j,k) OMP_SCHEDULE_ collapse(2) &
1431  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,DENS,Kh,FACT,PHI,qflx_phi,GSQRT,I_XYZ,RFDX,J13G,I_UYZ) &
1432  !$omp shared(RCDZ)
1433  !$acc kernels
1434  do j = jjs, jje
1435  do i = iis-1, iie
1436  do k = ks+1, ke-1
1437 #ifdef DEBUG
1438  call check( __line__, dens(k,i,j) )
1439  call check( __line__, dens(k,i+1,j) )
1440  call check( __line__, kh(k,i,j) )
1441  call check( __line__, kh(k,i+1,j) )
1442  call check( __line__, phi(k,i+1,j) )
1443  call check( __line__, phi(k,i,j) )
1444  call check( __line__, rfdx(i) )
1445 #endif
1446  qflx_phi(k,i,j,xdir) = - 0.25_rp & ! 1/2/2
1447  * ( dens(k,i,j) + dens(k,i+1,j) ) &
1448  * ( kh(k,i,j) + kh(k,i+1,j) ) * fact &
1449  * ( &
1450  ( gsqrt(k,i+1,j,i_xyz) * phi(k,i+1,j) &
1451  - gsqrt(k,i ,j,i_xyz) * phi(k,i ,j) ) * rfdx(i) &
1452  + ( j13g(k ,i,j,i_uyz) * ( phi(k+1,i+1,j)+phi(k+1,i,j)+phi(k ,i+1,j)+phi(k ,i,j) ) &
1453  - j13g(k-1,i,j,i_uyz) * ( phi(k ,i+1,j)+phi(k ,i,j)+phi(k-1,i+1,j)+phi(k-1,i,j) ) &
1454  ) * 0.25_rp * rcdz(k) &
1455  ) / gsqrt(k,i,j,i_uyz)
1456  enddo
1457  enddo
1458  enddo
1459  !$acc end kernels
1460 #ifdef DEBUG
1461  i = iundef; j = iundef; k = iundef
1462 #endif
1463  !$omp parallel do
1464  !$acc kernels
1465  do j = jjs, jje
1466  do i = iis-1, iie
1467 #ifdef DEBUG
1468  call check( __line__, dens(ks,i,j) )
1469  call check( __line__, dens(ks,i+1,j) )
1470  call check( __line__, kh(ks,i,j) )
1471  call check( __line__, kh(ks,i+1,j) )
1472  call check( __line__, phi(ks,i+1,j) )
1473  call check( __line__, phi(ks,i,j) )
1474  call check( __line__, rfdx(i) )
1475 #endif
1476  qflx_phi(ks,i,j,xdir) = - 0.25_rp & ! 1/2/2
1477  * ( dens(ks,i,j) + dens(ks,i+1,j) ) &
1478  * ( kh(ks,i,j) + kh(ks,i+1,j) ) * fact &
1479  * ( &
1480  ( gsqrt(ks,i+1,j,i_xyz) * phi(ks,i+1,j) &
1481  - gsqrt(ks,i ,j,i_xyz) * phi(ks,i ,j) ) * rfdx(i) &
1482  + ( j13g(ks+1,i,j,i_uyz) * ( phi(ks+1,i+1,j)+phi(ks+1,i,j) ) &
1483  - j13g(ks ,i,j,i_uyz) * ( phi(ks ,i+1,j)+phi(ks ,i,j) ) &
1484  ) * 0.5_rp * rfdz(ks) &
1485  ) * mapf(i,j,1,i_uy) / gsqrt(ks,i,j,i_uyz)
1486  qflx_phi(ke,i,j,xdir) = - 0.25_rp & ! 1/2/2
1487  * ( dens(ke,i,j) + dens(ke,i+1,j) ) &
1488  * ( kh(ke,i,j) + kh(ke,i+1,j) ) * fact &
1489  * ( &
1490  ( gsqrt(ke,i+1,j,i_xyz) * phi(ke,i+1,j) &
1491  - gsqrt(ke,i ,j,i_xyz) * phi(ke,i ,j) ) * rfdx(i) &
1492  + ( j13g(ke ,i,j,i_uyz) * ( phi(ke ,i+1,j)+phi(ke ,i,j) ) &
1493  - j13g(ke-1,i,j,i_uyz) * ( phi(ke-1,i+1,j)+phi(ke-1,i,j) ) &
1494  ) * 0.5_rp * rfdz(ke-1) &
1495  ) * mapf(i,j,1,i_uy) / gsqrt(ke,i,j,i_uyz)
1496  enddo
1497  enddo
1498  !$acc end kernels
1499 #ifdef DEBUG
1500  i = iundef; j = iundef; k = iundef
1501 #endif
1502  ! (z-x plane; x,v,z)
1503  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
1504  !$omp private(i,j,k) &
1505  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,Kh,FACT,PHI,RFDY,DENS,qflx_phi,GSQRT,I_XYZ,J23G,I_XVZ,RCDZ) &
1506  !$omp shared(MAPF,I_XV)
1507  !$acc kernels
1508  do j = jjs-1, jje
1509  do i = iis, iie
1510  do k = ks+1, ke-1
1511 #ifdef DEBUG
1512  call check( __line__, dens(k,i,j) )
1513  call check( __line__, dens(k,i,j+1) )
1514  call check( __line__, kh(k,i,j) )
1515  call check( __line__, kh(k,i,j+1) )
1516  call check( __line__, phi(k,i,j+1) )
1517  call check( __line__, phi(k,i,j) )
1518  call check( __line__, rfdy(j) )
1519 #endif
1520  qflx_phi(k,i,j,ydir) = - 0.25_rp &
1521  * ( dens(k,i,j) + dens(k,i,j+1) ) &
1522  * ( kh(k,i,j) + kh(k,i,j+1) ) * fact &
1523  * ( &
1524  ( gsqrt(k,i,j+1,i_xyz) * phi(k,i,j+1) &
1525  - gsqrt(k,i,j ,i_xyz) * phi(k,i,j ) ) * rfdy(j) &
1526  + ( j23g(k ,i,j,i_xvz) * ( phi(k+1,i,j+1)+phi(k+1,i,j)+phi(k ,i,j+1)+phi(k ,i,j) ) &
1527  - j23g(k-1,i,j,i_xvz) * ( phi(k ,i,j+1)+phi(k ,i,j)+phi(k-1,i,j+1)+phi(k-1,i,j) ) &
1528  ) * 0.25_rp * rcdz(k) &
1529  ) * mapf(i,j,2,i_xv) / gsqrt(k,i,j,i_xvz)
1530  enddo
1531  enddo
1532  enddo
1533  !$acc end kernels
1534 #ifdef DEBUG
1535  i = iundef; j = iundef; k = iundef
1536 #endif
1537  !$omp parallel do
1538  !$acc kernels
1539  do j = jjs-1, jje
1540  do i = iis, iie
1541 #ifdef DEBUG
1542  call check( __line__, dens(ks,i,j) )
1543  call check( __line__, dens(ks,i,j+1) )
1544  call check( __line__, kh(ks,i,j) )
1545  call check( __line__, kh(ks,i,j+1) )
1546  call check( __line__, phi(ks,i,j+1) )
1547  call check( __line__, phi(ks,i,j) )
1548  call check( __line__, rfdy(j) )
1549 #endif
1550  qflx_phi(ks,i,j,ydir) = - 0.25_rp &
1551  * ( dens(ks,i,j) + dens(ks,i,j+1) ) &
1552  * ( kh(ks,i,j) + kh(ks,i,j+1) ) * fact &
1553  * ( &
1554  ( gsqrt(ks,i,j+1,i_xyz) * phi(ks,i,j+1) &
1555  - gsqrt(ks,i,j ,i_xyz) * phi(ks,i,j ) ) * rfdy(j) &
1556  + ( j23g(ks+1,i,j,i_xvz) * ( phi(ks+1,i,j+1)+phi(ks+1,i,j) ) &
1557  - j23g(ks ,i,j,i_xvz) * ( phi(ks ,i,j+1)+phi(ks ,i,j) ) &
1558  ) * 0.5_rp * rfdz(ks) &
1559  ) * mapf(i,j,2,i_xv) / gsqrt(ks,i,j,i_xvz)
1560  qflx_phi(ke,i,j,ydir) = - 0.25_rp &
1561  * ( dens(ke,i,j) + dens(ke,i,j+1) ) &
1562  * ( kh(ke,i,j) + kh(ke,i,j+1) ) * fact &
1563  * ( &
1564  ( gsqrt(ke,i,j+1,i_xyz) * phi(ke,i,j+1) &
1565  - gsqrt(ke,i,j ,i_xyz) * phi(ke,i,j ) ) * rfdy(j) &
1566  + ( j23g(ke ,i,j,i_xvz) * ( phi(ke ,i,j+1)+phi(ke ,i,j) ) &
1567  - j23g(ke-1,i,j,i_xvz) * ( phi(ke-1,i,j+1)+phi(ke-1,i,j) ) &
1568  ) * 0.5_rp * rfdz(ke-1) &
1569  ) * mapf(i,j,2,i_xv) / gsqrt(ke,i,j,i_xvz)
1570  enddo
1571  enddo
1572  !$acc end kernels
1573 #ifdef DEBUG
1574  i = iundef; j = iundef; k = iundef
1575 #endif
1576 
1577  if ( (.not. horizontal) .and. implicit ) then
1578  call atmos_phy_tb_calc_tend_phi( tend, & ! (out)
1579  qflx_phi, & ! (in)
1580  gsqrt, j13g, j23g, j33g, mapf, & ! (in)
1581  iis, iie, jjs, jje ) ! (in)
1582 
1583  call matrix_solver_tridiagonal( ka, ks, ke, ia, iis, iie, ja, jjs, jje, &
1584  a(:,:,:), b(:,:,:), c(:,:,:), tend(:,:,:), &
1585  tend2(:,:,:) )
1586 
1587  !$omp parallel do collapse(2)
1588  !$acc kernels
1589  do j = jjs, jje
1590  do i = iis, iie
1591  do k = ks, ke-1
1592  qflx_phi(k,i,j,zdir) = qflx_phi(k,i,j,zdir) &
1593  - 0.25_rp & ! 1/2/2
1594  * ( dens(k,i,j)+dens(k+1,i,j) ) &
1595  * ( kh(k,i,j) + kh(k+1,i,j) ) * fact &
1596  * dt * ( tend2(k+1,i,j)-tend2(k,i,j) ) * rfdz(k) * j33g &
1597  / gsqrt(k,i,j,i_xyw)
1598  end do
1599  end do
1600  end do
1601  !$acc end kernels
1602 
1603  end if
1604 
1605  !$acc end data
1606 
1607  return
1608  end subroutine atmos_phy_tb_calc_flux_phi
1609 
1610  !-----------------------------------------------------------------------------
1611  subroutine atmos_phy_tb_calc_tend_momz( &
1612  MOMZ_t_TB, &
1613  QFLX_MOMZ, &
1614  GSQRT, J13G, J23G, J33G, MAPF, &
1615  IIS, IIE, JJS, JJE )
1617  rcdx => atmos_grid_cartesc_rcdx, &
1618  rcdy => atmos_grid_cartesc_rcdy, &
1619  rfdz => atmos_grid_cartesc_rfdz
1620  implicit none
1621 
1622  real(rp), intent(out) :: momz_t_tb(ka,ia,ja)
1623 
1624  real(rp), intent(in) :: qflx_momz(ka,ia,ja,3)
1625  real(rp), intent(in) :: gsqrt(ka,ia,ja,7)
1626  real(rp), intent(in) :: j13g(ka,ia,ja,7)
1627  real(rp), intent(in) :: j23g(ka,ia,ja,7)
1628  real(rp), intent(in) :: j33g
1629  real(rp), intent(in) :: mapf(ia,ja,2,4)
1630  integer , intent(in) :: iis
1631  integer , intent(in) :: iie
1632  integer , intent(in) :: jjs
1633  integer , intent(in) :: jje
1634 
1635  real(rp) :: fluxz(ka)
1636  integer :: k, i, j
1637 
1638  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
1639  !$omp private(i,j,k, &
1640  !$omp fluxZ) &
1641  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,I_XYZ,I_XYW,I_UYW,I_XVW,I_XY,I_XV,I_UY, &
1642  !$omp MOMZ_t_TB,QFLX_MOMZ,GSQRT,J13G,J23G,J33G,MAPF,RFDZ,RCDX,RCDY)
1643  !$acc kernels
1644  do j = jjs, jje
1645  !$acc loop private(fluxZ)
1646  do i = iis, iie
1647  do k = ks+1, ke-1
1648  fluxz(k) = ( ( qflx_momz(k ,i,j,xdir) + qflx_momz(k ,i-1,j,xdir) &
1649  + qflx_momz(k-1,i,j,xdir) + qflx_momz(k-1,i-1,j,xdir) ) * j13g(k,i,j,i_xyz) * mapf(i,j,1,i_xy) &
1650  + ( qflx_momz(k ,i,j,ydir) + qflx_momz(k ,i,j-1,ydir) &
1651  + qflx_momz(k-1,i,j,ydir) + qflx_momz(k-1,i,j-1,ydir) ) * j23g(k,i,j,i_xyz) * mapf(i,j,2,i_xy) &
1652  ) * 0.25_rp &
1653  + j33g * qflx_momz(k,i,j,zdir)
1654  end do
1655  fluxz(ks) = 0.0_rp
1656  fluxz(ke) = 0.0_rp
1657  do k = ks, ke-1
1658  momz_t_tb(k,i,j) = &
1659  - ( ( ( gsqrt(k,i ,j,i_uyw) * qflx_momz(k,i ,j,xdir) / mapf(i ,j,2,i_uy) &
1660  - gsqrt(k,i-1,j,i_uyw) * qflx_momz(k,i-1,j,xdir) / mapf(i-1,j,2,i_uy) ) * rcdx(i) &
1661  + ( gsqrt(k,i,j ,i_xvw) * qflx_momz(k,i,j ,ydir) / mapf(i,j ,1,i_xv) &
1662  - gsqrt(k,i,j-1,i_xvw) * qflx_momz(k,i,j-1,ydir) / mapf(i,j-1,1,i_xv) ) * rcdy(j) &
1663  ) * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) &
1664  + ( fluxz(k+1) - fluxz(k) ) * rfdz(k) &
1665  ) / gsqrt(k,i,j,i_xyw)
1666  enddo
1667  enddo
1668  enddo
1669  !$acc end kernels
1670 
1671  return
1672  end subroutine atmos_phy_tb_calc_tend_momz
1673 
1674  subroutine atmos_phy_tb_calc_tend_momx( &
1675  MOMX_t_TB, &
1676  QFLX_MOMX, &
1677  GSQRT, J13G, J23G, J33G, MAPF, &
1678  IIS, IIE, JJS, JJE )
1680  rcdz => atmos_grid_cartesc_rcdz, &
1681  rcdy => atmos_grid_cartesc_rcdy, &
1682  rfdx => atmos_grid_cartesc_rfdx
1683  implicit none
1684  real(rp), intent(out) :: momx_t_tb(ka,ia,ja)
1685 
1686  real(rp), intent(in) :: qflx_momx(ka,ia,ja,3)
1687  real(rp), intent(in) :: gsqrt(ka,ia,ja,7)
1688  real(rp), intent(in) :: mapf(ia,ja,2,4)
1689  real(rp), intent(in) :: j13g(ka,ia,ja,7)
1690  real(rp), intent(in) :: j23g(ka,ia,ja,7)
1691  real(rp), intent(in) :: j33g
1692  integer , intent(in) :: iis
1693  integer , intent(in) :: iie
1694  integer , intent(in) :: jjs
1695  integer , intent(in) :: jje
1696 
1697  real(rp) :: fluxz(ka)
1698  integer :: k, i, j
1699 
1700  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
1701  !$omp private(i,j,k, &
1702  !$omp fluxZ) &
1703  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,I_XYZ,I_UYZ,I_UVZ,I_UYW,I_UY,I_XY,I_UV, &
1704  !$omp MOMX_t_TB,QFLX_MOMX,GSQRT,J13G,J23G,J33G,MAPF,RCDZ,RFDX,RCDY)
1705  !$acc kernels
1706  do j = jjs, jje
1707  !$acc loop private(fluxZ)
1708  do i = iis, iie
1709  do k = ks, ke-1
1710  fluxz(k) = ( ( qflx_momx(k+1,i+1,j,xdir) + qflx_momx(k+1,i,j ,xdir) &
1711  + qflx_momx(k ,i+1,j,xdir) + qflx_momx(k ,i,j ,xdir) ) * j13g(k,i,j,i_uyw) * mapf(i,j,1,i_uy) &
1712  + ( qflx_momx(k+1,i ,j,ydir) + qflx_momx(k+1,i,j-1,ydir) &
1713  + qflx_momx(k ,i ,j,ydir) + qflx_momx(k ,i,j-1,ydir) ) * j23g(k,i,j,i_uyw) * mapf(i,j,2,i_uy) &
1714  ) * 0.25_rp &
1715  + j33g * qflx_momx(k,i,j,zdir)
1716  end do
1717  fluxz(ks-1) = 0.0_rp
1718  fluxz(ke ) = 0.0_rp
1719  do k = ks, ke
1720  momx_t_tb(k,i,j) = &
1721  - ( ( ( gsqrt(k,i+1,j,i_xyz) * qflx_momx(k,i+1,j,xdir) / mapf(i+1,j ,2,i_xy) &
1722  - gsqrt(k,i ,j,i_xyz) * qflx_momx(k,i ,j,xdir) / mapf(i ,j ,2,i_xy) ) * rfdx(i) &
1723  + ( gsqrt(k,i,j ,i_uvz) * qflx_momx(k,i,j ,ydir) / mapf(i ,j ,1,i_uv) &
1724  - gsqrt(k,i,j-1,i_uvz) * qflx_momx(k,i,j-1,ydir) / mapf(i ,j-1,1,i_uv) ) * rcdy(j) &
1725  ) * mapf(i,j,1,i_uy) * mapf(i,j,2,i_uy) &
1726  + ( fluxz(k) - fluxz(k-1) ) * rcdz(k) &
1727  ) / gsqrt(k,i,j,i_uyz)
1728  enddo
1729  enddo
1730  enddo
1731  !$acc end kernels
1732 
1733  return
1734  end subroutine atmos_phy_tb_calc_tend_momx
1735 
1736  subroutine atmos_phy_tb_calc_tend_momy( &
1737  MOMY_t_TB, &
1738  QFLX_MOMY, &
1739  GSQRT, J13G, J23G, J33G, MAPF, &
1740  IIS, IIE, JJS, JJE )
1742  rcdz => atmos_grid_cartesc_rcdz, &
1743  rcdx => atmos_grid_cartesc_rcdx, &
1744  rfdy => atmos_grid_cartesc_rfdy
1745  implicit none
1746 
1747  real(rp), intent(out) :: momy_t_tb(ka,ia,ja)
1748 
1749  real(rp), intent(in) :: qflx_momy(ka,ia,ja,3)
1750  real(rp), intent(in) :: gsqrt(ka,ia,ja,7)
1751  real(rp), intent(in) :: mapf(ia,ja,2,4)
1752  real(rp), intent(in) :: j13g(ka,ia,ja,7)
1753  real(rp), intent(in) :: j23g(ka,ia,ja,7)
1754  real(rp), intent(in) :: j33g
1755  integer , intent(in) :: iis
1756  integer , intent(in) :: iie
1757  integer , intent(in) :: jjs
1758  integer , intent(in) :: jje
1759 
1760  real(rp) :: fluxz(ka)
1761  integer :: k, i, j
1762 
1763  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
1764  !$omp private(i,j,k, &
1765  !$omp fluxZ) &
1766  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,I_XVZ,I_UVZ,I_XYZ,I_XVW,I_XV,I_UV,I_XY, &
1767  !$omp MOMY_t_TB,QFLX_MOMY,GSQRT,J13G,J23G,J33G,MAPF,RCDZ,RCDX,RFDY)
1768  !$acc kernels
1769  do j = jjs, jje
1770  !$acc loop private(fluxZ)
1771  do i = iis, iie
1772  do k = ks, ke-1
1773  fluxz(k) = ( ( qflx_momy(k+1,i,j ,xdir) + qflx_momy(k+1,i-1,j,xdir) &
1774  + qflx_momy(k ,i,j ,xdir) + qflx_momy(k ,i-1,j,xdir) ) * j13g(k,i,j,i_xvw) * mapf(i,j,1,i_xv) &
1775  + ( qflx_momy(k+1,i,j+1,ydir) + qflx_momy(k+1,i ,j,ydir) &
1776  + qflx_momy(k ,i,j+1,ydir) + qflx_momy(k ,i ,j,ydir) ) * j23g(k,i,j,i_xvw) * mapf(i,j,2,i_xv) &
1777  ) * 0.25_rp &
1778  + j33g * qflx_momy(k,i,j,zdir)
1779  end do
1780  fluxz(ks-1) = 0.0_rp
1781  fluxz(ke ) = 0.0_rp
1782  do k = ks, ke
1783  momy_t_tb(k,i,j) = &
1784  - ( ( ( gsqrt(k,i ,j ,i_uvz) * qflx_momy(k,i ,j,xdir) / mapf(i ,j,2,i_uv) &
1785  - gsqrt(k,i-1,j ,i_uvz) * qflx_momy(k,i-1,j,xdir) / mapf(i-1,j,2,i_uv) ) * rcdx(i) &
1786  + ( gsqrt(k,i ,j+1,i_xyz) * qflx_momy(k,i,j+1,ydir) / mapf(i,j+1,2,i_xy) &
1787  - gsqrt(k,i ,j ,i_xyz) * qflx_momy(k,i,j ,ydir) / mapf(i,j ,2,i_xy) ) * rfdy(j) &
1788  ) * mapf(i,j,1,i_xv) * mapf(i,j,2,i_xv) &
1789  + ( fluxz(k) - fluxz(k-1) ) * rcdz(k) &
1790  ) / gsqrt(k,i,j,i_xvz)
1791  enddo
1792  enddo
1793  enddo
1794  !$acc end kernels
1795 
1796  return
1797  end subroutine atmos_phy_tb_calc_tend_momy
1798 
1799  subroutine atmos_phy_tb_calc_tend_phi( &
1800  phi_t_TB, &
1801  QFLX_phi, &
1802  GSQRT, J13G, J23G, J33G, MAPF, &
1803  IIS, IIE, JJS, JJE )
1805  rcdz => atmos_grid_cartesc_rcdz, &
1806  rcdx => atmos_grid_cartesc_rcdx, &
1807  rcdy => atmos_grid_cartesc_rcdy, &
1808  fdz => atmos_grid_cartesc_fdz
1809  implicit none
1810 
1811  real(rp), intent(out) :: phi_t_tb(ka,ia,ja)
1812 
1813  real(rp), intent(in) :: qflx_phi(ka,ia,ja,3)
1814  real(rp), intent(in) :: gsqrt(ka,ia,ja,7)
1815  real(rp), intent(in) :: j13g(ka,ia,ja,7)
1816  real(rp), intent(in) :: j23g(ka,ia,ja,7)
1817  real(rp), intent(in) :: j33g
1818  real(rp), intent(in) :: mapf(ia,ja,2,4)
1819  integer , intent(in) :: iis
1820  integer , intent(in) :: iie
1821  integer , intent(in) :: jjs
1822  integer , intent(in) :: jje
1823 
1824  real(rp) :: fluxz(0:ka)
1825  integer :: k, i, j
1826 
1827  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
1828  !$omp private(i,j,k, &
1829  !$omp fluxZ) &
1830  !$omp shared(JJS,JJE,IIS,IIE,KS,KE,I_XYZ,I_UYZ,I_XVZ,I_XYW,I_XY,I_UY,I_XV, &
1831  !$omp phi_t_TB,QFLX_phi,GSQRT,MAPF, &
1832  !$omp J13G,J23G,J33G,RCDZ,RCDX,RCDY)
1833  !$acc kernels
1834  do j = jjs, jje
1835  !$acc loop private(fluxZ)
1836  do i = iis, iie
1837  do k = ks, ke-1
1838  fluxz(k) = ( ( qflx_phi(k+1,i,j,xdir) + qflx_phi(k+1,i-1,j,xdir) &
1839  + qflx_phi(k ,i,j,xdir) + qflx_phi(k ,i-1,j,xdir) ) * j13g(k,i,j,i_xyw) * mapf(i,j,1,i_xy) &
1840  + ( qflx_phi(k+1,i,j,ydir) + qflx_phi(k+1,i,j-1,ydir) &
1841  + qflx_phi(k ,i,j,ydir) + qflx_phi(k ,i,j-1,ydir) ) * j23g(k,i,j,i_xyw) * mapf(i,j,2,i_xy) &
1842  ) * 0.25_rp &
1843  + j33g * qflx_phi(k,i,j,zdir)
1844  end do
1845  fluxz(ks-1) = 0.0_rp
1846  fluxz(ke ) = 0.0_rp
1847  do k = ks, ke
1848  phi_t_tb(k,i,j) = &
1849  - ( ( ( gsqrt(k,i ,j,i_uyz) * qflx_phi(k,i ,j,xdir) / mapf(i ,j,2,i_uy) &
1850  - gsqrt(k,i-1,j,i_uyz) * qflx_phi(k,i-1,j,xdir) / mapf(i-1,j,2,i_uy) ) * rcdx(i) &
1851  + ( gsqrt(k,i,j ,i_xvz) * qflx_phi(k,i,j ,ydir) / mapf(i,j ,1,i_xv) &
1852  - gsqrt(k,i,j-1,i_xvz) * qflx_phi(k,i,j-1,ydir) / mapf(i,j-1,1,i_xv) ) * rcdy(j) &
1853  ) * mapf(i,j,1,i_xy) * mapf(i,j,2,i_xy) &
1854  + ( fluxz(k) - fluxz(k-1) ) * rcdz(k) &
1855  ) / gsqrt(k,i,j,i_xyz)
1856  enddo
1857  enddo
1858  enddo
1859  !$acc end kernels
1860 
1861  return
1862  end subroutine atmos_phy_tb_calc_tend_phi
1863 
1864 end module scale_atmos_phy_tb_common
scale_atmos_grid_cartesc_index::i_uy
integer, public i_uy
Definition: scale_atmos_grid_cartesC_index.F90:100
scale_atmos_grid_cartesc_index::ke
integer, public ke
end point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:52
scale_atmos_grid_cartesc_index::i_xv
integer, public i_xv
Definition: scale_atmos_grid_cartesC_index.F90:101
scale_atmos_grid_cartesc::atmos_grid_cartesc_rfdx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rfdx
reciprocal of face-dx
Definition: scale_atmos_grid_cartesC.F90:68
scale_const::const_undef2
integer, parameter, public const_undef2
undefined value (INT2)
Definition: scale_const.F90:40
scale_atmos_phy_tb_common::atmos_phy_tb_calc_strain_tensor
subroutine, public atmos_phy_tb_calc_strain_tensor(S33_C, S11_C, S22_C, S31_C, S12_C, S23_C, S12_Z, S23_X, S31_Y, S2, DENS, MOMZ, MOMX, MOMY, GSQRT, J13G, J23G, J33G, MAPF)
Definition: scale_atmos_phy_tb_common.F90:66
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_atmos_grid_cartesc_index::ka
integer, public ka
Definition: scale_atmos_grid_cartesC_index.F90:47
scale_atmos_grid_cartesc_index::jblock
integer, public jblock
block size for cache blocking: y
Definition: scale_atmos_grid_cartesC_index.F90:41
scale_atmos_grid_cartesc_index::i_uv
integer, public i_uv
Definition: scale_atmos_grid_cartesC_index.F90:102
scale_atmos_grid_cartesc::atmos_grid_cartesc_rcdx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rcdx
reciprocal of center-dx
Definition: scale_atmos_grid_cartesC.F90:66
scale_atmos_phy_tb_common::atmos_phy_tb_calc_tend_momx
subroutine, public atmos_phy_tb_calc_tend_momx(MOMX_t_TB, QFLX_MOMX, GSQRT, J13G, J23G, J33G, MAPF, IIS, IIE, JJS, JJE)
Definition: scale_atmos_phy_tb_common.F90:1679
scale_atmos_grid_cartesc_index::iblock
integer, public iblock
block size for cache blocking: x
Definition: scale_atmos_grid_cartesC_index.F90:40
scale_atmos_phy_tb_common::atmos_phy_tb_calc_tend_momy
subroutine, public atmos_phy_tb_calc_tend_momy(MOMY_t_TB, QFLX_MOMY, GSQRT, J13G, J23G, J33G, MAPF, IIS, IIE, JJS, JJE)
Definition: scale_atmos_phy_tb_common.F90:1741
scale_atmos_grid_cartesc::atmos_grid_cartesc_fdy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fdy
y-length of grid(j+1) to grid(j) [m]
Definition: scale_atmos_grid_cartesC.F90:64
scale_atmos_grid_cartesc::atmos_grid_cartesc_rcdz
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rcdz
reciprocal of center-dz
Definition: scale_atmos_grid_cartesC.F90:45
scale_atmos_grid_cartesc_index::i_xyz
integer, public i_xyz
Definition: scale_atmos_grid_cartesC_index.F90:91
scale_atmos_grid_cartesc_index::i_xy
integer, public i_xy
Definition: scale_atmos_grid_cartesC_index.F90:99
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_atmos_grid_cartesc_index::ie
integer, public ie
end point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:54
scale_atmos_grid_cartesc_index::i_uyz
integer, public i_uyz
Definition: scale_atmos_grid_cartesC_index.F90:95
scale_io
module STDIO
Definition: scale_io.F90:10
scale_atmos_grid_cartesc::atmos_grid_cartesc_rcdy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rcdy
reciprocal of center-dy
Definition: scale_atmos_grid_cartesC.F90:67
scale_atmos_grid_cartesc_index::i_xvw
integer, public i_xvw
Definition: scale_atmos_grid_cartesC_index.F90:94
scale_tracer::k
real(rp), public k
Definition: scale_tracer.F90:45
scale_atmos_grid_cartesc_index
module atmosphere / grid / cartesC index
Definition: scale_atmos_grid_cartesC_index.F90:12
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
scale_debug::check
subroutine, public check(current_line, v)
Undefined value checker.
Definition: scale_debug.F90:59
scale_atmos_grid_cartesc_index::zdir
integer, parameter, public zdir
Definition: scale_atmos_grid_cartesC_index.F90:32
scale_atmos_grid_cartesc_index::i_xvz
integer, public i_xvz
Definition: scale_atmos_grid_cartesC_index.F90:96
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_atmos_phy_tb_common::atmos_phy_tb_calc_flux_phi
subroutine, public atmos_phy_tb_calc_flux_phi(qflx_phi, DENS, PHI, Kh, FACT, GSQRT, J13G, J23G, J33G, MAPF, horizontal, implicit, a, b, c, dt, IIS, IIE, JJS, JJE)
Definition: scale_atmos_phy_tb_common.F90:1337
scale_atmos_grid_cartesc_index::i_uyw
integer, public i_uyw
Definition: scale_atmos_grid_cartesC_index.F90:93
scale_atmos_grid_cartesc_index::is
integer, public is
start point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:53
scale_precision::dp
integer, parameter, public dp
Definition: scale_precision.F90:32
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
scale_atmos_grid_cartesc::atmos_grid_cartesc_fdx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fdx
x-length of grid(i+1) to grid(i) [m]
Definition: scale_atmos_grid_cartesC.F90:63
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
scale_atmos_phy_tb_common::atmos_phy_tb_calc_tend_phi
subroutine, public atmos_phy_tb_calc_tend_phi(phi_t_TB, QFLX_phi, GSQRT, J13G, J23G, J33G, MAPF, IIS, IIE, JJS, JJE)
Definition: scale_atmos_phy_tb_common.F90:1804
scale_atmos_grid_cartesc_index::xdir
integer, parameter, public xdir
Definition: scale_atmos_grid_cartesC_index.F90:33
scale_atmos_grid_cartesc_index::ks
integer, public ks
start point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:51
scale_debug
module DEBUG
Definition: scale_debug.F90:11
scale_atmos_phy_tb_common
module ATMOSPHERE / Physics Turbulence
Definition: scale_atmos_phy_tb_common.F90:12
scale_atmos_grid_cartesc_index::js
integer, public js
start point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:55
scale_atmos_grid_cartesc_index::ydir
integer, parameter, public ydir
Definition: scale_atmos_grid_cartesC_index.F90:34
scale_matrix
module MATRIX
Definition: scale_matrix.F90:17
scale_atmos_grid_cartesc::atmos_grid_cartesc_rfdz
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rfdz
reciprocal of face-dz
Definition: scale_atmos_grid_cartesC.F90:46
scale_atmos_grid_cartesc
module atmosphere / grid / cartesC
Definition: scale_atmos_grid_cartesC.F90:12
scale_atmos_grid_cartesc::atmos_grid_cartesc_fdz
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fdz
z-length of grid(i+1) to grid(i) [m]
Definition: scale_atmos_grid_cartesC.F90:44
scale_atmos_grid_cartesc_index::i_uvz
integer, public i_uvz
Definition: scale_atmos_grid_cartesC_index.F90:97
scale_const::const_undef
real(rp), public const_undef
Definition: scale_const.F90:43
scale_atmos_phy_tb_common::atmos_phy_tb_calc_tend_momz
subroutine, public atmos_phy_tb_calc_tend_momz(MOMZ_t_TB, QFLX_MOMZ, GSQRT, J13G, J23G, J33G, MAPF, IIS, IIE, JJS, JJE)
Definition: scale_atmos_phy_tb_common.F90:1616
scale_atmos_grid_cartesc::atmos_grid_cartesc_rfdy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rfdy
reciprocal of face-dy
Definition: scale_atmos_grid_cartesC.F90:69
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56
scale_atmos_grid_cartesc_index::i_xyw
integer, public i_xyw
Definition: scale_atmos_grid_cartesC_index.F90:92