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