63 S33_C, S11_C, S22_C, &
64 S31_C, S12_C, S23_C, &
65 S12_Z, S23_X, S31_Y, &
67 DENS, MOMZ, MOMX, MOMY, &
68 GSQRT, J13G, J23G, J33G, MAPF )
93 real(RP),
intent(out) :: S33_C (
ka,
ia,
ja)
94 real(RP),
intent(out) :: S11_C (
ka,
ia,
ja)
95 real(RP),
intent(out) :: S22_C (
ka,
ia,
ja)
96 real(RP),
intent(out) :: S31_C (
ka,
ia,
ja)
97 real(RP),
intent(out) :: S12_C (
ka,
ia,
ja)
98 real(RP),
intent(out) :: S23_C (
ka,
ia,
ja)
99 real(RP),
intent(out) :: S12_Z (
ka,
ia,
ja)
100 real(RP),
intent(out) :: S23_X (
ka,
ia,
ja)
101 real(RP),
intent(out) :: S31_Y (
ka,
ia,
ja)
102 real(RP),
intent(out) :: S2 (
ka,
ia,
ja)
104 real(RP),
intent(in) :: DENS (
ka,
ia,
ja)
105 real(RP),
intent(in) :: MOMZ (
ka,
ia,
ja)
106 real(RP),
intent(in) :: MOMX (
ka,
ia,
ja)
107 real(RP),
intent(in) :: MOMY (
ka,
ia,
ja)
109 real(RP),
intent(in) :: GSQRT (
ka,
ia,
ja,7)
110 real(RP),
intent(in) :: J13G (
ka,
ia,
ja,7)
111 real(RP),
intent(in) :: J23G (
ka,
ia,
ja,7)
112 real(RP),
intent(in) :: J33G
113 real(RP),
intent(in) :: MAPF (
ia,
ja,2,4)
116 real(RP) :: VELZ_C (
ka,
ia,
ja)
117 real(RP) :: VELZ_XY(
ka,
ia,
ja)
118 real(RP) :: VELX_C (
ka,
ia,
ja)
119 real(RP) :: VELX_YZ(
ka,
ia,
ja)
120 real(RP) :: VELY_C (
ka,
ia,
ja)
121 real(RP) :: VELY_ZX(
ka,
ia,
ja)
124 real(RP) :: WORK_V(
ka,
ia,
ja)
125 real(RP) :: WORK_Z(
ka,
ia,
ja)
126 real(RP) :: WORK_X(
ka,
ia,
ja)
127 real(RP) :: WORK_Y(
ka,
ia,
ja)
129 integer :: IIS, IIE, JJS, JJE
144 velz_c(:,:,:) = undef
145 velz_xy(:,:,:) = undef
146 velx_c(:,:,:) = undef
147 velx_yz(:,:,:) = undef
148 vely_c(:,:,:) = undef
149 vely_zx(:,:,:) = undef
151 work_v(:,:,:) = undef
152 work_z(:,:,:) = undef
153 work_x(:,:,:) = undef
154 work_y(:,:,:) = undef
162 call check( __line__, momz(k,i,j) )
163 call check( __line__, dens(k+1,i,j) )
164 call check( __line__, dens(k,i,j) )
166 velz_xy(k,i,j) = 2.0_rp * momz(k,i,j) / ( dens(k+1,i,j)+dens(k,i,j) )
171 i = iundef; j = iundef; k = iundef
175 velz_xy(
ke,i,j) = 0.0_rp
179 i = iundef; j = iundef; k = iundef
185 call check( __line__, momz(k,i,j) )
186 call check( __line__, momz(k-1,i,j) )
187 call check( __line__, dens(k,i,j) )
189 velz_c(k,i,j) = 0.5_rp * ( momz(k,i,j) + momz(k-1,i,j) ) / dens(k,i,j)
194 i = iundef; j = iundef; k = iundef
199 call check( __line__, momz(
ks,i,j) )
200 call check( __line__, dens(
ks,i,j) )
202 velz_c(
ks,i,j) = 0.5_rp * momz(
ks,i,j) / dens(
ks,i,j)
207 i = iundef; j = iundef; k = iundef
214 call check( __line__, momx(k,i,j) )
215 call check( __line__, dens(k,i+1,j) )
216 call check( __line__, dens(k,i,j) )
218 velx_yz(k,i,j) = 2.0_rp * momx(k,i,j) / ( dens(k,i+1,j)+dens(k,i,j) )
223 i = iundef; j = iundef; k = iundef
227 velx_yz(
ke+1,i,j) = 0.0_rp
231 i = iundef; j = iundef; k = iundef
237 call check( __line__, momx(k,i,j) )
238 call check( __line__, momx(k,i-1,j) )
239 call check( __line__, dens(k,i,j) )
241 velx_c(k,i,j) = 0.5_rp * ( momx(k,i,j) + momx(k,i-1,j) ) / dens(k,i,j)
246 i = iundef; j = iundef; k = iundef
253 call check( __line__, momy(k,i,j) )
254 call check( __line__, dens(k,i,j+1) )
255 call check( __line__, dens(k,i,j) )
257 vely_zx(k,i,j) = 2.0_rp * momy(k,i,j) / ( dens(k,i,j+1)+dens(k,i,j) )
262 i = iundef; j = iundef; k = iundef
266 vely_zx(
ke+1,i,j) = 0.0_rp
270 i = iundef; j = iundef; k = iundef
276 call check( __line__, momy(k,i,j) )
277 call check( __line__, momy(k,i,j-1) )
278 call check( __line__, dens(k,i,j) )
280 vely_c(k,i,j) = 0.5_rp * ( momy(k,i,j) + momy(k,i,j-1) ) / dens(k,i,j)
285 i = iundef; j = iundef; k = iundef
294 work_z(:,:,:) = undef; work_x(:,:,:) = undef; work_y(:,:,:) = undef; work_v(:,:,:) = undef
304 call check( __line__, velz_c(k,i+1,j) )
305 call check( __line__, velz_c(k,i,j) )
307 work_x(k,i,j) = 0.5_rp * ( velz_c(k,i+1,j) + velz_c(k,i,j) )
312 i = iundef; j = iundef; k = iundef
316 work_x(
ke+1,i,j) = 0.0_rp
320 i = iundef; j = iundef; k = iundef
327 call check( __line__, velz_c(k,i,j+1) )
328 call check( __line__, velz_c(k,i,j) )
330 work_y(k,i,j) = 0.5_rp * ( velz_c(k,i,j+1) + velz_c(k,i,j) )
335 i = iundef; j = iundef; k = iundef
339 work_y(
ke+1,i,j) = 0.0_rp
343 i = iundef; j = iundef; k = iundef
352 call check( __line__, velz_xy(k,i,j) )
353 call check( __line__, velz_xy(k-1,i,j) )
354 call check( __line__, rcdz(k) )
356 s33_c(k,i,j) = ( velz_xy(k,i,j) - velz_xy(k-1,i,j) ) * rcdz(k) &
357 * j33g / gsqrt(k,i,j,
i_xyz)
362 i = iundef; j = iundef; k = iundef
367 call check( __line__, velz_xy(
ks,i,j) )
369 call check( __line__, rcdz(
ks) )
371 s33_c(
ks,i,j) = velz_xy(
ks,i,j) * rcdz(
ks) &
376 i = iundef; j = iundef; k = iundef
385 call check( __line__, velz_c(k,i+1,j) )
386 call check( __line__, velz_c(k,i-1,j) )
389 call check( __line__, velz_xy(k,i,j) )
390 call check( __line__, velz_xy(k-1,i,j) )
393 call check( __line__, fdx(i) )
394 call check( __line__, fdx(i-1) )
396 s31_c(k,i,j) = 0.5_rp * ( &
397 ( 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) ) &
398 + ( 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) &
405 i = iundef; j = iundef; k = iundef
410 call check( __line__, velz_c(
ks,i+1,j) )
411 call check( __line__, velz_c(
ks,i-1,j) )
414 call check( __line__, velz_xy(
ks,i,j) )
416 call check( __line__, velz_c(
ke,i+1,j) )
417 call check( __line__, velz_c(
ke,i-1,j) )
420 call check( __line__, velz_xy(
ke,i,j) )
422 call check( __line__, fdx(i) )
423 call check( __line__, fdx(i-1) )
425 s31_c(
ks,i,j) = 0.5_rp * ( &
426 ( 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) ) &
427 + ( j13g(
ks,i,j,
i_xyw)*velz_xy(
ks,i,j) ) * rcdz(
ks) &
429 s31_c(
ke,i,j) = 0.5_rp * ( &
430 ( 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) ) &
431 - ( j13g(
ke-1,i,j,
i_xyw)*velz_xy(
ke-1,i,j) ) * rcdz(
ke) &
436 i = iundef; j = iundef; k = iundef
444 call check( __line__, velz_xy(k,i+1,j) )
445 call check( __line__, velz_xy(k,i,j) )
446 call check( __line__, rfdx(i) )
448 s31_y(k,i,j) = 0.5_rp * ( &
449 ( 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) &
450 + ( 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) &
456 i = iundef; j = iundef; k = iundef
465 call check( __line__, velz_c(k,i,j+1) )
466 call check( __line__, velz_c(k,i,j-1) )
469 call check( __line__, velz_xy(k,i,j) )
470 call check( __line__, velz_xy(k-1,i,j) )
473 call check( __line__, fdy(j) )
474 call check( __line__, fdy(j-1) )
476 s23_c(k,i,j) = 0.5_rp * ( &
477 ( 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) ) &
478 + ( 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) &
484 i = iundef; j = iundef; k = iundef
489 call check( __line__, velz_c(
ks,i,j+1) )
490 call check( __line__, velz_c(
ks,i,j-1) )
493 call check( __line__, velz_xy(
ks,i,j) )
495 call check( __line__, velz_c(
ke,i,j+1) )
496 call check( __line__, velz_c(
ke,i,j-1) )
499 call check( __line__, velz_xy(
ke,i,j) )
501 call check( __line__, fdy(j) )
502 call check( __line__, fdy(j-1) )
504 s23_c(
ks,i,j) = 0.5_rp * ( &
505 ( 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) ) &
506 + ( j23g(
ks,i,j,
i_xyw)*velz_xy(
ks,i,j) ) * rcdz(
ks) &
508 s23_c(
ke,i,j) = 0.5_rp * ( &
509 ( 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) ) &
510 - ( j23g(
ke-1,i,j,
i_xyw)*velz_xy(
ke-1,i,j) ) * rcdz(
ke) &
515 i = iundef; j = iundef; k = iundef
523 call check( __line__, velz_xy(k,i,j+1) )
524 call check( __line__, velz_xy(k,i,j) )
525 call check( __line__, rfdy(j) )
527 s23_x(k,i,j) = 0.5_rp * ( &
528 ( 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) &
529 + ( 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) &
535 i = iundef; j = iundef; k = iundef
539 work_z(:,:,:) = undef; work_x(:,:,:) = undef; work_y(:,:,:) = undef; work_v(:,:,:) = undef
547 call check( __line__, velx_c(k+1,i,j) )
548 call check( __line__, velx_c(k,i,j) )
550 work_z(k,i,j) = 0.5_rp * ( velx_c(k+1,i,j) + velx_c(k,i,j) )
555 i = iundef; j = iundef; k = iundef
564 call check( __line__, velx_c(k,i,j+1) )
565 call check( __line__, velx_c(k,i,j) )
567 work_y(k,i,j) = 0.5_rp * ( velx_c(k,i,j+1) + velx_c(k,i,j) )
572 i = iundef; j = iundef; k = iundef
579 call check( __line__, velx_yz(k,i,j) )
580 call check( __line__, velx_yz(k,i,j+1) )
581 call check( __line__, velx_yz(k+1,i,j) )
582 call check( __line__, velx_yz(k+1,i,j+1) )
588 work_v(k,i,j) = 0.25_rp &
589 * ( j23g(k ,i,j ,
i_uyz)*velx_yz(k ,i,j ) &
590 + j23g(k+1,i,j ,
i_uyz)*velx_yz(k+1,i,j ) &
591 + j23g(k ,i,j+1,
i_uyz)*velx_yz(k ,i,j+1) &
592 + j23g(k+1,i,j+1,
i_uyz)*velx_yz(k+1,i,j+1) )
597 i = iundef; j = iundef; k = iundef
606 call check( __line__, velx_yz(k,i,j) )
607 call check( __line__, velx_yz(k,i-1,j) )
610 call check( __line__, work_z(k,i,j) )
611 call check( __line__, work_z(k-1,i,j) )
615 call check( __line__, rcdx(i) )
618 ( 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) &
619 + ( 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) &
620 ) * mapf(i,j,1,
i_xy) / gsqrt(k,i,j,
i_xyz)
625 i = iundef; j = iundef; k = iundef
630 call check( __line__, velx_yz(
ks,i,j) )
631 call check( __line__, velx_yz(
ks,i-1,j) )
634 call check( __line__, velx_c(
ks+1,i,j) )
635 call check( __line__, velx_c(
ks,i,j) )
639 call check( __line__, velx_yz(
ke,i,j) )
640 call check( __line__, velx_yz(
ke,i-1,j) )
643 call check( __line__, velx_c(
ke,i,j) )
644 call check( __line__, velx_c(
ke-1,i,j) )
648 call check( __line__, rcdx(i) )
651 ( 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) &
652 + ( 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) &
655 ( 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) &
656 + ( 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) &
661 i = iundef; j = iundef; k = iundef
670 call check( __line__, s31_c(k,i,j) )
671 call check( __line__, velx_c(k+1,i,j) )
672 call check( __line__, velx_c(k-1,i,j) )
673 call check( __line__, fdz(k) )
674 call check( __line__, fdz(k-1) )
676 s31_c(k,i,j) = ( s31_c(k,i,j) &
677 + 0.5_rp * ( velx_c(k+1,i,j) - velx_c(k-1,i,j) ) * j33g / ( fdz(k) + fdz(k-1) ) &
678 ) / gsqrt(k,i,j,
i_xyz)
683 i = iundef; j = iundef; k = iundef
688 call check( __line__, s31_c(
ks,i,j) )
689 call check( __line__, velx_c(
ks+1,i,j) )
690 call check( __line__, velx_c(
ks,i,j) )
691 call check( __line__, rfdz(
ks) )
692 call check( __line__, s31_c(
ke,i,j) )
693 call check( __line__, velx_c(
ke,i,j) )
694 call check( __line__, velx_c(
ke-1,i,j) )
695 call check( __line__, rfdz(
ke-1) )
697 s31_c(
ks,i,j) = ( s31_c(
ks,i,j) &
698 + 0.5_rp * ( velx_c(
ks+1,i,j) - velx_c(
ks,i,j) ) * j33g * rfdz(
ks) &
700 s31_c(
ke,i,j) = ( s31_c(
ke,i,j) &
701 + 0.5_rp * ( velx_c(
ke,i,j) - velx_c(
ke-1,i,j) ) * j33g * rfdz(
ke-1) &
706 i = iundef; j = iundef; k = iundef
713 call check( __line__, s31_y(k,i,j) )
714 call check( __line__, velx_yz(k+1,i,j) )
715 call check( __line__, velx_yz(k,i,j) )
716 call check( __line__, rfdz(k) )
718 s31_y(k,i,j) = ( s31_y(k,i,j) &
719 + 0.5_rp * ( velx_yz(k+1,i,j) - velx_yz(k,i,j) ) * j33g * rfdz(k) &
720 ) / gsqrt(k,i,j,
i_uyw)
725 i = iundef; j = iundef; k = iundef
734 call check( __line__, velx_c(k,i,j+1) )
735 call check( __line__, velx_c(k,i,j-1) )
738 call check( __line__, work_z(k,i,j) )
739 call check( __line__, work_z(k-1,i,j) )
743 s12_c(k,i,j) = 0.5_rp * ( &
744 ( 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) ) &
745 + ( 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) &
746 ) * mapf(i,j,2,
i_xy) / gsqrt(k,i,j,
i_xyz)
751 i = iundef; j = iundef; k = iundef
756 call check( __line__, velx_c(
ks,i,j+1) )
757 call check( __line__, velx_c(
ks,i,j-1) )
760 call check( __line__, velx_c(
ks+1,i,j) )
761 call check( __line__, velx_c(
ks,i,j) )
765 call check( __line__, velx_c(
ke,i,j+1) )
766 call check( __line__, velx_c(
ke,i,j-1) )
769 call check( __line__, velx_c(
ke,i,j) )
770 call check( __line__, velx_c(
ke-1,i,j) )
774 call check( __line__, fdy(j) )
775 call check( __line__, fdy(j-1) )
777 s12_c(
ks,i,j) = 0.5_rp * ( &
778 ( 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) ) &
779 + ( 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) &
781 s12_c(
ke,i,j) = 0.5_rp * ( &
782 ( 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) ) &
783 + ( 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) &
788 i = iundef; j = iundef; k = iundef
796 call check( __line__, velx_yz(k,i,j+1) )
797 call check( __line__, velx_yz(k,i,j) )
798 call check( __line__, work_v(k,i,j) )
799 call check( __line__, work_v(k-1,i,j) )
800 call check( __line__, rfdy(j) )
802 s12_z(k,i,j) = 0.5_rp * ( &
803 ( 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) &
804 + ( work_v(k,i,j) - work_v(k-1,i,j) ) * rcdz(k) &
810 i = iundef; j = iundef; k = iundef
815 call check( __line__, velx_yz(
ks,i,j+1) )
816 call check( __line__, velx_yz(
ks,i,j) )
817 call check( __line__, velx_yz(
ks+1,i,j) )
818 call check( __line__, velx_yz(
ks+1,i,j+1) )
821 call check( __line__, velx_yz(
ke,i,j+1) )
822 call check( __line__, velx_yz(
ke,i,j) )
823 call check( __line__, velx_yz(
ke-1,i,j) )
824 call check( __line__, velx_yz(
ke-1,i,j+1) )
828 s12_z(
ks,i,j) = 0.25_rp * ( &
829 ( 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) &
830 + ( j23g(
ks+1,i,j,
i_uvz) * ( velx_yz(
ks+1,i,j) + velx_yz(
ks+1,i,j+1) ) &
831 - j23g(
ks ,i,j,
i_uvz) * ( velx_yz(
ks ,i,j) + velx_yz(
ks ,i,j+1) ) ) * rfdz(
ks) &
833 s12_z(
ke,i,j) = 0.25_rp * ( &
834 ( 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) &
835 + ( j23g(
ke ,i,j,
i_uvz) * ( velx_yz(
ke ,i,j) + velx_yz(
ke ,i,j+1) ) &
836 - j23g(
ke-1,i,j,
i_uvz) * ( velx_yz(
ke-1,i,j) + velx_yz(
ke-1,i,j+1) ) ) * rfdz(
ke-1) &
841 i = iundef; j = iundef; k = iundef
845 work_z(:,:,:) = undef; work_x(:,:,:) = undef; work_y(:,:,:) = undef; work_v(:,:,:) = undef
853 call check( __line__, vely_c(k+1,i,j) )
854 call check( __line__, vely_c(k,i,j) )
856 work_z(k,i,j) = 0.5_rp * ( vely_c(k+1,i,j) + vely_c(k,i,j) )
861 i = iundef; j = iundef; k = iundef
868 call check( __line__, vely_c(k,i+1,j) )
869 call check( __line__, vely_c(k,i,j) )
871 work_x(k,i,j) = 0.5_rp * ( velz_c(k,i+1,j) + velz_c(k,i,j) )
876 i = iundef; j = iundef; k = iundef
885 call check( __line__, vely_zx(k,i,j) )
886 call check( __line__, vely_zx(k+1,i,j) )
887 call check( __line__, vely_zx(k,i+1,j) )
888 call check( __line__, vely_zx(k+1,i+1,j) )
890 work_v(k,i,j) = 0.25_rp &
891 * ( j13g(k ,i ,j,
i_xvz)*vely_zx(k ,i ,j) &
892 + j13g(k+1,i ,j,
i_xvz)*vely_zx(k+1,i ,j) &
893 + j13g(k ,i+1,j,
i_xvz)*vely_zx(k ,i+1,j) &
894 + j13g(k+1,i+1,j,
i_xvz)*vely_zx(k+1,i+1,j) )
899 i = iundef; j = iundef; k = iundef
908 call check( __line__, vely_zx(k,i,j) )
909 call check( __line__, vely_zx(k,i,j-1) )
912 call check( __line__, work_z(k,i,j) )
913 call check( __line__, work_z(k-1,i,j) )
916 call check( __line__, rcdy(j) )
919 ( 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) &
920 + ( 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) &
921 ) * mapf(i,j,2,
i_xy) / gsqrt(k,i,j,
i_xyz)
926 i = iundef; j = iundef; k = iundef
931 call check( __line__, vely_zx(
ks,i,j) )
932 call check( __line__, vely_zx(
ks,i,j-1) )
935 call check( __line__, vely_c(
ks+1,i,j) )
936 call check( __line__, vely_c(
ks,i,j) )
939 call check( __line__, rcdy(j) )
940 call check( __line__, vely_zx(
ke,i,j) )
941 call check( __line__, vely_zx(
ke,i,j-1) )
944 call check( __line__, vely_c(
ke,i,j) )
945 call check( __line__, vely_c(
ke-1,i,j) )
950 ( 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) &
951 + ( 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) &
954 ( 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) &
955 + ( 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) &
960 i = iundef; j = iundef; k = iundef
969 call check( __line__, s12_c(k,i,j) )
970 call check( __line__, vely_c(k,i+1,j) )
971 call check( __line__, vely_c(k,i-1,j) )
974 call check( __line__, work_z(k,i,j) )
975 call check( __line__, work_z(k-1,i,j) )
979 call check( __line__, fdx(i) )
980 call check( __line__, fdx(i-1) )
982 s12_c(k,i,j) = ( s12_c(k,i,j) &
984 ( 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) ) &
985 + ( 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) &
986 ) / gsqrt(k,i,j,
i_xyz)
991 i = iundef; j = iundef; k = iundef
996 call check( __line__, s12_c(
ks,i,j) )
997 call check( __line__, vely_c(
ks,i+1,j) )
998 call check( __line__, vely_c(
ks,i-1,j) )
1001 call check( __line__, vely_c(
ks+1,i,j) )
1002 call check( __line__, vely_c(
ks,i,j) )
1006 call check( __line__, s12_c(
ke,i,j) )
1007 call check( __line__, vely_c(
ke,i+1,j) )
1008 call check( __line__, vely_c(
ke,i-1,j) )
1011 call check( __line__, vely_c(
ke,i,j) )
1012 call check( __line__, vely_c(
ke-1,i,j) )
1016 call check( __line__, fdx(i) )
1017 call check( __line__, fdx(i-1) )
1019 s12_c(
ks,i,j) = ( s12_c(
ks,i,j) &
1021 ( 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) ) &
1022 + ( 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) ) &
1023 * mapf(i,j,1,
i_xy) &
1025 s12_c(
ke,i,j) = ( s12_c(
ke,i,j) &
1027 ( 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) ) &
1028 + ( 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) ) &
1029 * mapf(i,j,1,
i_xy) &
1034 i = iundef; j = iundef; k = iundef
1041 call check( __line__, s12_z(k,i,j) )
1042 call check( __line__, vely_zx(k,i+1,j) )
1043 call check( __line__, vely_zx(k,i,j) )
1044 call check( __line__, work_v(k,i,j) )
1045 call check( __line__, work_v(k-1,i,j) )
1046 call check( __line__, rfdx(i) )
1048 s12_z(k,i,j) = ( s12_z(k,i,j) &
1050 ( 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) &
1051 + ( work_v(k,i,j) - work_v(k-1,i,j) ) * rcdz(k) ) * mapf(i,j,1,
i_uv) &
1052 ) / gsqrt(k,i,j,
i_uvz)
1057 i = iundef; j = iundef; k = iundef
1062 call check( __line__, s12_z(
ks,i,j) )
1063 call check( __line__, vely_zx(
ks,i+1,j) )
1064 call check( __line__, vely_zx(
ks,i,j) )
1065 call check( __line__, vely_zx(
ks+1,i,j) )
1066 call check( __line__, vely_zx(
ks+1,i+1,j) )
1067 call check( __line__, s12_z(
ke,i,j) )
1068 call check( __line__, vely_zx(
ke,i+1,j) )
1069 call check( __line__, vely_zx(
ke,i,j) )
1070 call check( __line__, vely_zx(
ke-1,i,j) )
1071 call check( __line__, vely_zx(
ke-1,i+1,j) )
1072 call check( __line__, rfdx(i) )
1074 s12_z(
ks,i,j) = ( s12_z(
ks,i,j) &
1076 ( 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) &
1077 + ( j13g(
ks+1,i,j,
i_uvz) * ( vely_zx(
ks+1,i,j) + vely_zx(
ks+1,i+1,j) ) &
1078 - 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) &
1080 s12_z(
ke,i,j) = ( s12_z(
ke,i,j) &
1082 ( 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) &
1083 + ( j13g(
ke ,i,j,
i_uvz) * ( vely_zx(
ke ,i,j) + vely_zx(
ke ,i+1,j) ) &
1084 - 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) &
1089 i = iundef; j = iundef; k = iundef
1098 call check( __line__, s23_c(k,i,j) )
1099 call check( __line__, vely_c(k+1,i,j) )
1100 call check( __line__, vely_c(k-1,i,j) )
1101 call check( __line__, fdz(k) )
1102 call check( __line__, fdz(k-1) )
1104 s23_c(k,i,j) = ( s23_c(k,i,j) &
1105 + 0.5_rp * ( vely_c(k+1,i,j) - vely_c(k-1,i,j) ) * j33g / ( fdz(k) + fdz(k-1) ) &
1106 ) / gsqrt(k,i,j,
i_xyz)
1111 i = iundef; j = iundef; k = iundef
1116 call check( __line__, s23_c(
ks,i,j) )
1117 call check( __line__, vely_c(
ks+1,i,j) )
1118 call check( __line__, vely_c(
ks,i,j) )
1119 call check( __line__, rfdz(
ks) )
1120 call check( __line__, s23_c(
ke,i,j) )
1121 call check( __line__, vely_c(
ke,i,j) )
1122 call check( __line__, vely_c(
ke-1,i,j) )
1123 call check( __line__, rfdz(
ke-1) )
1125 s23_c(
ks,i,j) = ( s23_c(
ks,i,j) &
1126 + 0.5_rp * ( vely_c(
ks+1,i,j) - vely_c(
ks,i,j) ) * j33g * rfdz(
ks) &
1128 s23_c(
ke,i,j) = ( s23_c(
ke,i,j) &
1129 + 0.5_rp * ( vely_c(
ke,i,j) - vely_c(
ke-1,i,j) ) * j33g * rfdz(
ke-1) &
1134 i = iundef; j = iundef; k = iundef
1142 call check( __line__, s23_x(k,i,j) )
1143 call check( __line__, vely_zx(k+1,i,j) )
1144 call check( __line__, vely_zx(k,i,j) )
1145 call check( __line__, rfdz(k) )
1147 s23_x(k,i,j) = ( s23_x(k,i,j) &
1148 + 0.5_rp * ( vely_zx(k+1,i,j) - vely_zx(k,i,j) ) * j33g * rfdz(k) &
1149 ) / gsqrt(k,i,j,
i_xvw)
1154 i = iundef; j = iundef; k = iundef
1160 work_z(:,:,:) = undef; work_x(:,:,:) = undef; work_y(:,:,:) = undef
1167 call check( __line__, s11_c(k,i,j) )
1168 call check( __line__, s22_c(k,i,j) )
1169 call check( __line__, s33_c(k,i,j) )
1170 call check( __line__, s31_c(k,i,j) )
1171 call check( __line__, s12_c(k,i,j) )
1172 call check( __line__, s23_c(k,i,j) )
1174 s2(k,i,j) = max( 1e-10_rp, &
1175 2.0_rp * ( s11_c(k,i,j)**2 + s22_c(k,i,j)**2 + s33_c(k,i,j)**2 ) &
1176 + 4.0_rp * ( s31_c(k,i,j)**2 + s12_c(k,i,j)**2 + s23_c(k,i,j)**2 ) )
1181 i = iundef; j = iundef; k = iundef
1193 DENS, PHI, Kh, FACT, &
1194 GSQRT, J13G, J23G, J33G, MAPF, &
1197 IIS, IIE, JJS, JJE )
1222 real(RP),
intent(out) :: qflx_phi(
ka,
ia,
ja,3)
1223 real(RP),
intent(in) :: DENS(
ka,
ia,
ja)
1224 real(RP),
intent(in) :: PHI(
ka,
ia,
ja)
1225 real(RP),
intent(in) :: Kh(
ka,
ia,
ja)
1226 real(RP),
intent(in) :: FACT
1227 real(RP),
intent(in) :: GSQRT(
ka,
ia,
ja,7)
1228 real(RP),
intent(in) :: J13G(
ka,
ia,
ja,7)
1229 real(RP),
intent(in) :: J23G(
ka,
ia,
ja,7)
1230 real(RP),
intent(in) :: J33G
1231 real(RP),
intent(in) :: MAPF(
ia,
ja,2,4)
1232 real(RP),
intent(in) :: a(
ka,
ia,
ja)
1233 real(RP),
intent(in) :: b(
ka,
ia,
ja)
1234 real(RP),
intent(in) :: c(
ka,
ia,
ja)
1235 real(DP),
intent(in) :: dt
1236 logical,
intent(in) :: implicit
1237 integer,
intent(in) :: IIS
1238 integer,
intent(in) :: IIE
1239 integer,
intent(in) :: JJS
1240 integer,
intent(in) :: JJE
1242 real(RP) :: TEND(
ka,
ia,
ja)
1252 call check( __line__, dens(k,i,j) )
1253 call check( __line__, dens(k+1,i,j) )
1254 call check( __line__, kh(k,i,j) )
1255 call check( __line__, kh(k+1,i,j) )
1256 call check( __line__, phi(k+1,i,j) )
1257 call check( __line__, phi(k,i,j) )
1258 call check( __line__, rfdz(k) )
1260 qflx_phi(k,i,j,
zdir) = - 0.25_rp &
1261 * ( dens(k,i,j)+dens(k+1,i,j) ) &
1262 * ( kh(k,i,j) + kh(k+1,i,j) ) &
1263 * ( phi(k+1,i,j)-phi(k,i,j) ) * rfdz(k) * j33g &
1264 / gsqrt(k,i,j,
i_xyw)
1269 i = iundef; j = iundef; k = iundef
1273 qflx_phi(
ks-1,i,j,
zdir) = 0.0_rp
1274 qflx_phi(
ke ,i,j,
zdir) = 0.0_rp
1278 i = iundef; j = iundef; k = iundef
1286 call check( __line__, dens(k,i,j) )
1287 call check( __line__, dens(k,i+1,j) )
1288 call check( __line__, kh(k,i,j) )
1289 call check( __line__, kh(k,i+1,j) )
1290 call check( __line__, phi(k,i+1,j) )
1291 call check( __line__, phi(k,i,j) )
1292 call check( __line__, rfdx(i) )
1294 qflx_phi(k,i,j,
xdir) = - 0.25_rp &
1295 * ( dens(k,i,j)+dens(k,i+1,j) ) &
1296 * ( kh(k,i,j) + kh(k,i+1,j) ) &
1298 ( gsqrt(k,i+1,j,
i_xyz) * phi(k,i+1,j) &
1299 - gsqrt(k,i ,j,
i_xyz) * phi(k,i ,j) ) * rfdx(i) &
1300 + ( j13g(k+1,i,j,
i_uyz) * ( phi(k+1,i+1,j)+phi(k+1,i,j) ) &
1301 - j13g(k-1,i,j,
i_uyz) * ( phi(k-1,i+1,j)+phi(k-1,i,j) ) &
1302 ) * 0.5_rp / ( fdz(k)+fdz(k-1) ) &
1303 ) / gsqrt(k,i,j,
i_uyz)
1308 i = iundef; j = iundef; k = iundef
1313 call check( __line__, dens(
ks,i,j) )
1314 call check( __line__, dens(
ks,i+1,j) )
1315 call check( __line__, kh(
ks,i,j) )
1316 call check( __line__, kh(
ks,i+1,j) )
1317 call check( __line__, phi(
ks,i+1,j) )
1318 call check( __line__, phi(
ks,i,j) )
1319 call check( __line__, rfdx(i) )
1321 qflx_phi(
ks,i,j,
xdir) = - 0.25_rp &
1322 * ( dens(
ks,i,j)+dens(
ks,i+1,j) ) &
1323 * ( kh(
ks,i,j) + kh(
ks,i+1,j) ) &
1325 ( gsqrt(
ks,i+1,j,
i_xyz) * phi(
ks,i+1,j) &
1326 - gsqrt(
ks,i ,j,
i_xyz) * phi(
ks,i ,j) ) * rfdx(i) &
1327 + ( j13g(
ks+1,i,j,
i_uyz) * ( phi(
ks+1,i+1,j)+phi(
ks+1,i,j) ) &
1328 - j13g(
ks ,i,j,
i_uyz) * ( phi(
ks ,i+1,j)+phi(
ks ,i,j) ) &
1329 ) * 0.5_rp * rfdz(
ks) &
1331 qflx_phi(
ke,i,j,
xdir) = - 0.25_rp &
1332 * ( dens(
ke,i,j)+dens(
ke,i+1,j) ) &
1333 * ( kh(
ke,i,j) + kh(
ke,i+1,j) ) &
1335 ( gsqrt(
ke,i+1,j,
i_xyz) * phi(
ke,i+1,j) &
1336 - gsqrt(
ke,i ,j,
i_xyz) * phi(
ke,i ,j) ) * rfdx(i) &
1337 + ( j13g(
ke ,i,j,
i_uyz) * ( phi(
ke ,i+1,j)+phi(
ke ,i,j) ) &
1338 - j13g(
ke-1,i,j,
i_uyz) * ( phi(
ke-1,i+1,j)+phi(
ke-1,i,j) ) &
1339 ) * 0.5_rp * rfdz(
ke-1) &
1344 i = iundef; j = iundef; k = iundef
1351 call check( __line__, dens(k,i,j) )
1352 call check( __line__, dens(k,i,j+1) )
1353 call check( __line__, kh(k,i,j) )
1354 call check( __line__, kh(k,i,j+1) )
1355 call check( __line__, phi(k,i,j+1) )
1356 call check( __line__, phi(k,i,j) )
1357 call check( __line__, rfdy(j) )
1359 qflx_phi(k,i,j,
ydir) = - 0.25_rp &
1360 * ( dens(k,i,j)+dens(k,i,j+1) ) &
1361 * ( kh(k,i,j) + kh(k,i,j+1) ) &
1363 ( gsqrt(k,i,j+1,
i_xyz) * phi(k,i,j+1) &
1364 - gsqrt(k,i,j ,
i_xyz) * phi(k,i,j ) ) * rfdy(j) &
1365 + ( j23g(k+1,i,j,
i_xvz) * ( phi(k+1,i,j+1)+phi(k+1,i,j) ) &
1366 - j23g(k-1,i,j,
i_xvz) * ( phi(k-1,i,j+1)+phi(k-1,i,j) ) &
1367 ) * 0.5_rp / ( fdz(k)+fdz(k-1) ) &
1368 ) * mapf(i,j,2,
i_xv) / gsqrt(k,i,j,
i_xvz)
1373 i = iundef; j = iundef; k = iundef
1378 call check( __line__, dens(
ks,i,j) )
1379 call check( __line__, dens(
ks,i,j+1) )
1380 call check( __line__, kh(
ks,i,j) )
1381 call check( __line__, kh(
ks,i,j+1) )
1382 call check( __line__, phi(
ks,i,j+1) )
1383 call check( __line__, phi(
ks,i,j) )
1384 call check( __line__, rfdy(j) )
1386 qflx_phi(
ks,i,j,
ydir) = - 0.25_rp &
1387 * ( dens(
ks,i,j)+dens(
ks,i,j+1) ) &
1388 * ( kh(
ks,i,j) + kh(
ks,i,j+1) ) &
1390 ( gsqrt(
ks,i,j+1,
i_xyz) * phi(
ks,i,j+1) &
1391 - gsqrt(
ks,i,j ,
i_xyz) * phi(
ks,i,j ) ) * rfdy(j) &
1392 + ( j23g(
ks+1,i,j,
i_xvz) * ( phi(
ks+1,i,j+1)+phi(
ks+1,i,j) ) &
1393 - j23g(
ks ,i,j,
i_xvz) * ( phi(
ks ,i,j+1)+phi(
ks ,i,j) ) &
1394 ) * 0.5_rp * rfdz(
ks) &
1396 qflx_phi(
ke,i,j,
ydir) = - 0.25_rp &
1397 * ( dens(
ke,i,j)+dens(
ke,i,j+1) ) &
1398 * ( kh(
ke,i,j) + kh(
ke,i,j+1) ) &
1400 ( gsqrt(
ke,i,j+1,
i_xyz) * phi(
ke,i,j+1) &
1401 - gsqrt(
ke,i,j ,
i_xyz) * phi(
ke,i,j ) ) * rfdy(j) &
1402 + ( j23g(
ke ,i,j,
i_xvz) * ( phi(
ke ,i,j+1)+phi(
ke ,i,j) ) &
1403 - j23g(
ke-1,i,j,
i_xvz) * ( phi(
ke-1,i,j+1)+phi(
ke-1,i,j) ) &
1404 ) * 0.5_rp * rfdz(
ke-1) &
1409 i = iundef; j = iundef; k = iundef
1412 if (
implicit )
then 1415 gsqrt, j13g, j23g, j33g, mapf, &
1416 iis, iie, jjs, jje )
1427 a(:,i,j), b(:,i,j), c(:,i,j), d, &
1431 qflx_phi(k,i,j,
zdir) = qflx_phi(k,i,j,
zdir) &
1433 * ( dens(k,i,j)+dens(k+1,i,j) ) &
1434 * ( kh(k,i,j) + kh(k+1,i,j) ) &
1435 * dt * ( tend(k+1,i,j)-tend(k,i,j) ) * rfdz(k) * j33g &
1436 / gsqrt(k,i,j,
i_xyw)
1454 real(RP),
intent(out) :: phi(
ka)
1455 real(RP),
intent(in) :: a(
ka)
1456 real(RP),
intent(in) :: b(
ka)
1457 real(RP),
intent(in) :: c(
ka)
1458 real(RP),
intent(in) :: d(
ka)
1459 integer,
intent(in) :: KE_TB
1467 do k =
ks+1, ke_tb-1
1468 denom = b(k) + c(k)*e(k-1)
1469 e(k) = - a(k) / denom
1470 f(k) = ( d(k) - c(k)*f(k-1) ) / denom
1474 phi(ke_tb) = ( d(ke_tb) - c(ke_tb)*f(ke_tb-1) ) / ( b(ke_tb) + c(ke_tb)*e(ke_tb-1) )
1476 do k = ke_tb-1,
ks, -1
1477 phi(k) = e(k) * phi(k+1) + f(k)
1493 GSQRT, J13G, J23G, J33G, MAPF, &
1494 IIS, IIE, JJS, JJE )
1509 real(RP),
intent(out) :: MOMZ_t_TB(
ka,
ia,
ja)
1510 real(RP),
intent(in) :: QFLX_MOMZ(
ka,
ia,
ja,3)
1511 real(RP),
intent(in) :: GSQRT(
ka,
ia,
ja,7)
1512 real(RP),
intent(in) :: J13G(
ka,
ia,
ja,7)
1513 real(RP),
intent(in) :: J23G(
ka,
ia,
ja,7)
1514 real(RP),
intent(in) :: J33G
1515 real(RP),
intent(in) :: MAPF(
ia,
ja,2,4)
1516 integer ,
intent(in) :: IIS
1517 integer ,
intent(in) :: IIE
1518 integer ,
intent(in) :: JJS
1519 integer ,
intent(in) :: JJE
1526 momz_t_tb(k,i,j) = &
1527 - ( ( gsqrt(k,i ,j,
i_uyw) * qflx_momz(k,i ,j,
xdir) &
1528 - gsqrt(k,i-1,j,
i_uyw) * qflx_momz(k,i-1,j,
xdir) ) * rcdx(i) * mapf(i,j,1,
i_xy) &
1529 + ( gsqrt(k,i,j ,
i_xvw) * qflx_momz(k,i,j ,
ydir) &
1530 - gsqrt(k,i,j-1,
i_xvw) * qflx_momz(k,i,j-1,
ydir) ) * rcdy(j) &
1531 + ( ( j13g(k+1,i,j,
i_xyz) * ( qflx_momz(k+1,i,j,
xdir) + qflx_momz(k+1,i-1,j,
xdir) ) &
1532 - j13g(k-1,i,j,
i_xyz) * ( qflx_momz(k-1,i,j,
xdir) + qflx_momz(k-1,i-1,j,
xdir) ) &
1533 ) * mapf(i,j,1,
i_xy) &
1534 + ( j23g(k+1,i,j,
i_xyz) * ( qflx_momz(k+1,i,j,
ydir) + qflx_momz(k+1,i,j-1,
ydir) ) &
1535 - j23g(k-1,i,j,
i_xyz) * ( qflx_momz(k-1,i,j,
ydir) + qflx_momz(k-1,i,j-1,
ydir) ) &
1536 ) * mapf(i,j,2,
i_xy) &
1537 ) * 0.5_rp / ( cdz(k+1)+cdz(k) ) &
1538 + j33g * ( qflx_momz(k+1,i,j,
zdir) - qflx_momz(k,i,j,
zdir) ) * rfdz(k) ) &
1539 / gsqrt(k,i,j,
i_xyw)
1546 momz_t_tb(
ks,i,j) = &
1548 - gsqrt(
ks,i-1,j,
i_uyw) * qflx_momz(
ks,i-1,j,
xdir) ) * rcdx(i) * mapf(i,j,1,
i_xy) &
1550 - gsqrt(
ks,i,j-1,
i_xvw) * qflx_momz(
ks,i,j-1,
ydir) ) * rcdy(j) * mapf(i,j,2,
i_xy) &
1553 ) * mapf(i,j,1,
i_xy) &
1556 ) * mapf(i,j,2,
i_xy) &
1557 ) * 0.5_rp * rcdz(
ks+1) &
1558 + j33g * ( qflx_momz(
ks+1,i,j,
zdir) - qflx_momz(
ks,i,j,
zdir) ) * rfdz(
ks) ) &
1560 momz_t_tb(
ke-1,i,j) = &
1562 - gsqrt(
ke-1,i-1,j,
i_uyw) * qflx_momz(
ke-1,i-1,j,
xdir) ) * rcdx(i) * mapf(i,j,1,
i_xy) &
1564 - gsqrt(
ke-1,i,j-1,
i_xvw) * qflx_momz(
ke-1,i,j-1,
ydir) ) * rcdy(j) * mapf(i,j,2,
i_xy) &
1567 ) * mapf(i,j,1,
i_xy) &
1570 ) * mapf(i,j,2,
i_xy) &
1571 ) * 0.5_rp * rcdz(
ke-1) &
1572 + j33g * ( qflx_momz(
ke,i,j,
zdir) - qflx_momz(
ke-1,i,j,
zdir) ) * rfdz(
ke-1) ) &
1583 GSQRT, J13G, J23G, J33G, MAPF, &
1584 IIS, IIE, JJS, JJE )
1598 real(RP),
intent(out) :: MOMX_t_TB(
ka,
ia,
ja)
1599 real(RP),
intent(in) :: QFLX_MOMX(
ka,
ia,
ja,3)
1600 real(RP),
intent(in) :: GSQRT(
ka,
ia,
ja,7)
1601 real(RP),
intent(in) :: MAPF(
ia,
ja,2,4)
1602 real(RP),
intent(in) :: J13G(
ka,
ia,
ja,7)
1603 real(RP),
intent(in) :: J23G(
ka,
ia,
ja,7)
1604 real(RP),
intent(in) :: J33G
1605 integer ,
intent(in) :: IIS
1606 integer ,
intent(in) :: IIE
1607 integer ,
intent(in) :: JJS
1608 integer ,
intent(in) :: JJE
1615 momx_t_tb(k,i,j) = &
1616 - ( ( gsqrt(k,i+1,j,
i_xyz) * qflx_momx(k,i+1,j,
xdir) &
1617 - gsqrt(k,i ,j,
i_xyz) * qflx_momx(k,i ,j,
xdir) ) * rfdx(i) * mapf(i,j,1,
i_uy) &
1618 + ( gsqrt(k,i,j ,
i_uvz) * qflx_momx(k,i,j ,
ydir) &
1619 - gsqrt(k,i,j-1,
i_uvz) * qflx_momx(k,i,j-1,
ydir) ) * rcdy(j) * mapf(i,j,2,
i_uy) &
1620 + ( ( j13g(k+1,i,j,
i_uyw) * ( qflx_momx(k+1,i+1,j,
xdir) + qflx_momx(k+1,i,j ,
xdir) ) &
1621 - j13g(k-1,i,j,
i_uyw) * ( qflx_momx(k-1,i+1,j,
xdir) + qflx_momx(k-1,i,j ,
xdir) ) &
1622 ) * mapf(i,j,1,
i_uy) &
1623 + ( j23g(k+1,i,j,
i_uyw) * ( qflx_momx(k+1,i ,j,
ydir) + qflx_momx(k+1,i,j-1,
ydir) ) &
1624 - j23g(k-1,i,j,
i_uyw) * ( qflx_momx(k-1,i ,j,
ydir) + qflx_momx(k-1,i,j-1,
ydir) ) &
1625 ) * mapf(i,j,2,
i_uy) &
1626 ) * 0.5_rp / ( fdz(k)+fdz(k-1) ) &
1627 + j33g * ( qflx_momx(k,i,j,
zdir) - qflx_momx(k-1,i,j,
zdir) ) * rcdz(k) ) &
1628 / gsqrt(k,i,j,
i_uyz)
1634 momx_t_tb(
ks,i,j) = &
1636 - gsqrt(
ks,i ,j,
i_xyz) * qflx_momx(
ks,i ,j,
xdir) ) * rfdx(i) * mapf(i,j,1,
i_uy) &
1638 - gsqrt(
ks,i,j-1,
i_uvz) * qflx_momx(
ks,i,j-1,
ydir) ) * rcdy(j) &
1641 ) * mapf(i,j,1,
i_uy) &
1644 ) * mapf(i,j,2,
i_uy) &
1645 ) * 0.5_rp * rcdz(
ks) &
1646 + j33g * ( qflx_momx(
ks,i,j,
zdir) ) * rfdz(
ks) ) &
1648 momx_t_tb(
ke,i,j) = &
1650 - gsqrt(
ke,i ,j,
i_xyz) * qflx_momx(
ke,i ,j,
xdir) ) * rfdx(i) * mapf(i,j,1,
i_uy) &
1652 - gsqrt(
ke,i,j-1,
i_uvz) * qflx_momx(
ke,i,j-1,
ydir) ) * rcdy(j) * mapf(i,j,2,
i_uy)&
1655 ) * mapf(i,j,1,
i_uy) &
1658 ) * mapf(i,j,2,
i_uy) &
1659 ) * 0.5_rp * rfdz(
ke-1) &
1660 - j33g * ( qflx_momx(
ke-1,i,j,
zdir) ) * rcdz(
ke) ) &
1671 GSQRT, J13G, J23G, J33G, MAPF, &
1672 IIS, IIE, JJS, JJE )
1686 real(RP),
intent(out) :: MOMY_t_TB(
ka,
ia,
ja)
1687 real(RP),
intent(in) :: QFLX_MOMY(
ka,
ia,
ja,3)
1688 real(RP),
intent(in) :: GSQRT(
ka,
ia,
ja,7)
1689 real(RP),
intent(in) :: MAPF(
ia,
ja,2,4)
1690 real(RP),
intent(in) :: J13G(
ka,
ia,
ja,7)
1691 real(RP),
intent(in) :: J23G(
ka,
ia,
ja,7)
1692 real(RP),
intent(in) :: J33G
1693 integer ,
intent(in) :: IIS
1694 integer ,
intent(in) :: IIE
1695 integer ,
intent(in) :: JJS
1696 integer ,
intent(in) :: JJE
1703 momy_t_tb(k,i,j) = &
1704 - ( ( gsqrt(k,i ,j ,
i_uvz) * qflx_momy(k,i ,j,
xdir) &
1705 - gsqrt(k,i-1,j ,
i_uvz) * qflx_momy(k,i-1,j,
xdir) ) * rcdx(i) * mapf(i,j,1,
i_xv) &
1706 + ( gsqrt(k,i ,j+1,
i_xyz) * qflx_momy(k,i,j+1,
ydir) &
1707 - gsqrt(k,i ,j ,
i_xyz) * qflx_momy(k,i,j ,
ydir) ) * rfdy(j) * mapf(i,j,2,
i_xv) &
1708 + ( ( j13g(k+1,i,j ,
i_xvw) * ( qflx_momy(k+1,i,j ,
xdir) + qflx_momy(k+1,i-1,j,
xdir) ) &
1709 - j13g(k-1,i,j ,
i_xvw) * ( qflx_momy(k-1,i,j ,
xdir) + qflx_momy(k-1,i-1,j,
xdir) ) &
1710 ) * mapf(i,j,1,
i_xv) &
1711 + ( j23g(k+1,i,j+1,
i_xvw) * ( qflx_momy(k+1,i,j+1,
ydir) + qflx_momy(k+1,i ,j,
ydir) ) &
1712 - j23g(k-1,i,j+1,
i_xvw) * ( qflx_momy(k-1,i,j+1,
ydir) + qflx_momy(k-1,i ,j,
ydir) ) &
1713 ) * mapf(i,j,2,
i_xv) &
1714 ) * 0.5_rp / ( fdz(k)+fdz(k-1) ) &
1715 + j33g * ( qflx_momy(k,i,j,
zdir) - qflx_momy(k-1,i,j,
zdir) ) * rcdz(k) ) &
1716 / gsqrt(k,i,j,
i_xvw)
1722 momy_t_tb(
ks,i,j) = &
1724 - gsqrt(
ks,i-1,j ,
i_uvz) * qflx_momy(
ks,i-1,j,
xdir) ) * rcdx(i) * mapf(i,j,1,
i_xv) &
1726 - gsqrt(
ks,i ,j ,
i_xyz) * qflx_momy(
ks,i,j ,
ydir) ) * rfdy(j) * mapf(i,j,2,
i_xv) &
1727 + ( ( j13g(
ks+1,i,j ,
i_xvw) * ( qflx_momy(
ks+1,i,j ,
xdir) + qflx_momy(
ks+1,i-1,j,
xdir) ) &
1729 ) * mapf(i,j,1,
i_xv) &
1732 ) * mapf(i,j,2,
i_xv) &
1733 ) * 0.5_rp * rfdz(
ks) &
1734 + j33g * ( qflx_momy(
ks,i,j,
zdir) ) * rcdz(
ks) ) &
1736 momy_t_tb(
ke,i,j) = &
1738 - gsqrt(
ke,i-1,j ,
i_uvz) * qflx_momy(
ke,i-1,j,
xdir) ) * rcdx(i) * mapf(i,j,1,
i_xv) &
1740 - gsqrt(
ke,i ,j ,
i_xyz) * qflx_momy(
ke,i,j ,
ydir) ) * rfdy(j) * mapf(i,j,2,
i_xv) &
1743 ) * mapf(i,j,1,
i_xv) &
1746 ) * mapf(i,j,2,
i_xv) &
1747 ) * 0.5_rp * rfdz(
ke-1) &
1748 - j33g * ( qflx_momy(
ke-1,i,j,
zdir) ) * rcdz(
ke) ) &
1759 GSQRT, J13G, J23G, J33G, MAPF, &
1760 IIS, IIE, JJS, JJE )
1776 real(RP),
intent(out) :: phi_t_TB(
ka,
ia,
ja)
1777 real(RP),
intent(in) :: QFLX_phi(
ka,
ia,
ja,3)
1778 real(RP),
intent(in) :: GSQRT(
ka,
ia,
ja,7)
1779 real(RP),
intent(in) :: MAPF(
ia,
ja,2,4)
1780 real(RP),
intent(in) :: J13G(
ka,
ia,
ja,7)
1781 real(RP),
intent(in) :: J23G(
ka,
ia,
ja,7)
1782 real(RP),
intent(in) :: J33G
1783 integer ,
intent(in) :: IIS
1784 integer ,
intent(in) :: IIE
1785 integer ,
intent(in) :: JJS
1786 integer ,
intent(in) :: JJE
1794 - ( ( gsqrt(k,i ,j,
i_uyz) * qflx_phi(k,i ,j,
xdir) &
1795 - gsqrt(k,i-1,j,
i_uvz) * qflx_phi(k,i-1,j,
xdir) ) * rcdx(i) * mapf(i,j,1,
i_xy) &
1796 + ( gsqrt(k,i,j ,
i_xvz) * qflx_phi(k,i,j ,
ydir) &
1797 - gsqrt(k,i,j-1,
i_xvz) * qflx_phi(k,i,j-1,
ydir) ) * rcdy(j) * mapf(i,j,2,
i_xy) &
1798 + ( ( j13g(k+1,i,j,
i_xyw) * ( qflx_phi(k+1,i,j,
xdir) + qflx_phi(k+1,i-1,j,
xdir) ) &
1799 - j13g(k-1,i,j,
i_xyw) * ( qflx_phi(k-1,i,j,
xdir) + qflx_phi(k-1,i-1,j,
xdir) ) &
1800 ) * mapf(i,j,1,
i_xy) &
1801 + ( j23g(k+1,i,j,
i_xyw) * ( qflx_phi(k+1,i,j,
ydir) + qflx_phi(k+1,i,j-1,
ydir) ) &
1802 - j23g(k-1,i,j,
i_xyw) * ( qflx_phi(k-1,i,j,
ydir) + qflx_phi(k-1,i,j-1,
ydir) ) &
1803 ) * mapf(i,j,2,
i_xy) &
1804 ) * 0.5_rp / ( fdz(k) + fdz(k-1) ) &
1805 + j33g * ( qflx_phi(k,i,j,
zdir) - qflx_phi(k-1,i,j,
zdir) ) * rcdz(k) ) &
1806 / gsqrt(k,i,j,
i_xyz)
1812 phi_t_tb(
ks,i,j) = &
1814 - gsqrt(
ks,i-1,j,
i_uvz) * qflx_phi(
ks,i-1,j,
xdir) ) * rcdx(i) * mapf(i,j,1,
i_xy) &
1816 - gsqrt(
ks,i,j-1,
i_xvz) * qflx_phi(
ks,i,j-1,
ydir) ) * rcdy(j) * mapf(i,j,2,
i_xy) &
1819 ) * mapf(i,j,1,
i_xy) &
1822 ) * mapf(i,j,2,
i_xy) &
1823 ) * 0.5_rp * rfdz(
ks) &
1824 + j33g * ( qflx_phi(
ks,i,j,
zdir) ) * rcdz(
ks) ) &
1826 phi_t_tb(
ke,i,j) = &
1828 - gsqrt(
ke,i-1,j,
i_uvz) * qflx_phi(
ke,i-1,j,
xdir) ) * rcdx(i) * mapf(i,j,1,
i_xy) &
1830 - gsqrt(
ke,i,j-1,
i_xvz) * qflx_phi(
ke,i,j-1,
ydir) ) * rcdy(j) * mapf(i,j,2,
i_xy) &
1833 ) * mapf(i,j,1,
i_xy) &
1836 ) * mapf(i,j,2,
i_xy) &
1837 ) * 0.5_rp * rfdz(
ke-1) &
1838 - j33g * ( qflx_phi(
ke-1,i,j,
zdir) ) * rcdz(
ke) ) &
integer, public is
start point of inner domain: x, local
subroutine, public atmos_phy_tb_calc_tend_momz(MOMZ_t_TB, QFLX_MOMZ, GSQRT, J13G, J23G, J33G, MAPF, IIS, IIE, JJS, JJE)
integer, public je
end point of inner domain: y, local
real(rp), dimension(:), allocatable, public grid_rcdy
reciprocal of center-dy
subroutine, public atmos_phy_tb_diffusion_solver(phi, a, b, c, d, KE_TB)
subroutine, public atmos_phy_tb_calc_flux_phi(qflx_phi, DENS, PHI, Kh, FACT, GSQRT, J13G, J23G, J33G, MAPF, a, b, c, dt, implicit, IIS, IIE, JJS, JJE)
real(rp), dimension(:), allocatable, public grid_fdy
y-length of grid(j+1) to grid(j) [m]
integer, public iblock
block size for cache blocking: x
subroutine, public atmos_phy_tb_calc_tend_momy(MOMY_t_TB, QFLX_MOMY, GSQRT, J13G, J23G, J33G, MAPF, IIS, IIE, JJS, JJE)
real(rp), dimension(:), allocatable, public grid_rcdx
reciprocal of center-dx
integer, parameter, public zdir
integer, parameter, public ydir
integer, public ke
end point of inner domain: z, local
integer, parameter, public xdir
real(rp), dimension(:), allocatable, public grid_rfdy
reciprocal of face-dy
subroutine, public check(current_line, v)
Undefined value checker.
real(rp), dimension(:), allocatable, public grid_rcdz
reciprocal of center-dz
real(rp), public const_undef
subroutine, public atmos_phy_tb_calc_tend_phi(phi_t_TB, QFLX_phi, GSQRT, J13G, J23G, J33G, MAPF, IIS, IIE, JJS, JJE)
integer, public ia
of x whole cells (local, with HALO)
subroutine, public atmos_phy_tb_calc_tend_momx(MOMX_t_TB, QFLX_MOMX, GSQRT, J13G, J23G, J33G, MAPF, IIS, IIE, JJS, JJE)
real(rp), dimension(:), allocatable, public grid_fdz
z-length of grid(k+1) to grid(k) [m]
integer, public ka
of z whole cells (local, with HALO)
integer, public jblock
block size for cache blocking: y
integer, public js
start point of inner domain: y, local
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)
integer, parameter, public const_undef2
undefined value (INT2)
integer, public ks
start point of inner domain: z, local
module ATMOSPHERE / Physics Turbulence
integer, public ie
end point of inner domain: x, local
real(rp), dimension(:), allocatable, public grid_cdz
z-length of control volume [m]
real(rp), dimension(:), allocatable, public grid_fdx
x-length of grid(i+1) to grid(i) [m]
real(rp), dimension(:), allocatable, public grid_rfdx
reciprocal of face-dx
real(rp), dimension(:), allocatable, public grid_rfdz
reciprocal of face-dz
integer, public ja
of y whole cells (local, with HALO)