Go to the documentation of this file.
61 S33_C, S11_C, S22_C, &
62 S31_C, S12_C, S23_C, &
63 S12_Z, S23_X, S31_Y, &
65 DENS, MOMZ, MOMX, MOMY, &
66 GSQRT, J13G, J23G, J33G, MAPF )
79 real(
rp),
intent(out) :: s33_c (
ka,
ia,
ja)
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)
86 real(
rp),
intent(out) :: s23_x (
ka,
ia,
ja)
87 real(
rp),
intent(out) :: s31_y (
ka,
ia,
ja)
88 real(
rp),
intent(out) :: s2 (
ka,
ia,
ja)
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)
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)
115 integer :: iis, iie, jjs, jje
130 velz_c(:,:,:) = undef
131 velz_xy(:,:,:) = undef
132 velx_c(:,:,:) = undef
133 velx_yz(:,:,:) = undef
134 vely_c(:,:,:) = undef
135 vely_zx(:,:,:) = undef
137 work_v(:,:,:) = undef
138 work_z(:,:,:) = undef
139 work_x(:,:,:) = undef
140 work_y(:,:,:) = undef
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) )
153 velz_xy(
k,i,j) = 2.0_rp * momz(
k,i,j) / ( dens(
k+1,i,j)+dens(
k,i,j) )
158 i = iundef; j = iundef;
k = iundef
163 velz_xy(
ke,i,j) = 0.0_rp
167 i = iundef; j = iundef;
k = iundef
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) )
178 velz_c(
k,i,j) = 0.5_rp * ( momz(
k,i,j) + momz(
k-1,i,j) ) / dens(
k,i,j)
183 i = iundef; j = iundef;
k = iundef
189 call check( __line__, momz(
ks,i,j) )
190 call check( __line__, dens(
ks,i,j) )
192 velz_c(
ks,i,j) = 0.5_rp * momz(
ks,i,j) / dens(
ks,i,j)
197 i = iundef; j = iundef;
k = iundef
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) )
209 velx_yz(
k,i,j) = 2.0_rp * momx(
k,i,j) / ( dens(
k,i+1,j)+dens(
k,i,j) )
214 i = iundef; j = iundef;
k = iundef
219 velx_yz(
ke+1,i,j) = 0.0_rp
223 i = iundef; j = iundef;
k = iundef
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) )
234 velx_c(
k,i,j) = 0.5_rp * ( momx(
k,i,j) + momx(
k,i-1,j) ) / dens(
k,i,j)
239 i = iundef; j = iundef;
k = iundef
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) )
252 vely_zx(
k,i,j) = 2.0_rp * momy(
k,i,j) / ( dens(
k,i,j+1)+dens(
k,i,j) )
257 i = iundef; j = iundef;
k = iundef
262 vely_zx(
ke+1,i,j) = 0.0_rp
266 i = iundef; j = iundef;
k = iundef
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) )
277 vely_c(
k,i,j) = 0.5_rp * ( momy(
k,i,j) + momy(
k,i,j-1) ) / dens(
k,i,j)
282 i = iundef; j = iundef;
k = iundef
291 work_z(:,:,:) = undef; work_x(:,:,:) = undef; work_y(:,:,:) = undef; work_v(:,:,:) = undef
302 call check( __line__, velz_c(
k,i+1,j) )
303 call check( __line__, velz_c(
k,i,j) )
305 work_x(
k,i,j) = 0.5_rp * ( velz_c(
k,i+1,j) + velz_c(
k,i,j) )
310 i = iundef; j = iundef;
k = iundef
315 work_x(
ke+1,i,j) = 0.0_rp
319 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
340 work_y(
ke+1,i,j) = 0.0_rp
344 i = iundef; j = iundef;
k = iundef
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) )
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)
364 i = iundef; j = iundef;
k = iundef
370 call check( __line__, velz_xy(
ks,i,j) )
372 call check( __line__, rcdz(
ks) )
374 s33_c(
ks,i,j) = velz_xy(
ks,i,j) * rcdz(
ks) &
379 i = iundef; j = iundef;
k = iundef
389 call check( __line__, velz_c(
k,i+1,j) )
390 call check( __line__, velz_c(
k,i-1,j) )
393 call check( __line__, velz_xy(
k,i,j) )
394 call check( __line__, velz_xy(
k-1,i,j) )
397 call check( __line__, fdx(i) )
398 call check( __line__, fdx(i-1) )
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) &
409 i = iundef; j = iundef;
k = iundef
415 call check( __line__, velz_c(
ks,i+1,j) )
416 call check( __line__, velz_c(
ks,i-1,j) )
419 call check( __line__, velz_xy(
ks,i,j) )
421 call check( __line__, velz_c(
ke,i+1,j) )
422 call check( __line__, velz_c(
ke,i-1,j) )
425 call check( __line__, velz_xy(
ke,i,j) )
427 call check( __line__, fdx(i) )
428 call check( __line__, fdx(i-1) )
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) &
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) &
441 i = iundef; j = iundef;
k = iundef
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) )
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) &
462 i = iundef; j = iundef;
k = iundef
472 call check( __line__, velz_c(
k,i,j+1) )
473 call check( __line__, velz_c(
k,i,j-1) )
476 call check( __line__, velz_xy(
k,i,j) )
477 call check( __line__, velz_xy(
k-1,i,j) )
480 call check( __line__, fdy(j) )
481 call check( __line__, fdy(j-1) )
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) &
491 i = iundef; j = iundef;
k = iundef
497 call check( __line__, velz_c(
ks,i,j+1) )
498 call check( __line__, velz_c(
ks,i,j-1) )
501 call check( __line__, velz_xy(
ks,i,j) )
503 call check( __line__, velz_c(
ke,i,j+1) )
504 call check( __line__, velz_c(
ke,i,j-1) )
507 call check( __line__, velz_xy(
ke,i,j) )
509 call check( __line__, fdy(j) )
510 call check( __line__, fdy(j-1) )
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) &
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) &
523 i = iundef; j = iundef;
k = iundef
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) )
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) &
544 i = iundef; j = iundef;
k = iundef
548 work_z(:,:,:) = undef; work_x(:,:,:) = undef; work_y(:,:,:) = undef; work_v(:,:,:) = undef
557 call check( __line__, velx_c(
k+1,i,j) )
558 call check( __line__, velx_c(
k,i,j) )
560 work_z(
k,i,j) = 0.5_rp * ( velx_c(
k+1,i,j) + velx_c(
k,i,j) )
565 i = iundef; j = iundef;
k = iundef
575 call check( __line__, velx_c(
k,i,j+1) )
576 call check( __line__, velx_c(
k,i,j) )
578 work_y(
k,i,j) = 0.5_rp * ( velx_c(
k,i,j+1) + velx_c(
k,i,j) )
583 i = iundef; j = iundef;
k = iundef
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) )
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) )
609 i = iundef; j = iundef;
k = iundef
619 call check( __line__, velx_yz(
k,i,j) )
620 call check( __line__, velx_yz(
k,i-1,j) )
623 call check( __line__, work_z(
k,i,j) )
624 call check( __line__, work_z(
k-1,i,j) )
628 call check( __line__, rcdx(i) )
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) &
638 i = iundef; j = iundef;
k = iundef
644 call check( __line__, velx_yz(
ks,i,j) )
645 call check( __line__, velx_yz(
ks,i-1,j) )
648 call check( __line__, velx_c(
ks+1,i,j) )
649 call check( __line__, velx_c(
ks,i,j) )
653 call check( __line__, velx_yz(
ke,i,j) )
654 call check( __line__, velx_yz(
ke,i-1,j) )
657 call check( __line__, velx_c(
ke,i,j) )
658 call check( __line__, velx_c(
ke-1,i,j) )
662 call check( __line__, rcdx(i) )
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) &
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) &
675 i = iundef; j = iundef;
k = iundef
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) )
691 s31_c(
k,i,j) = ( s31_c(
k,i,j) &
692 + 0.5_rp * ( velx_c(
k+1,i,j) - velx_c(
k-1,i,j) ) * j33g / ( fdz(
k) + fdz(
k-1) ) &
698 i = iundef; j = iundef;
k = iundef
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) )
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) &
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) &
722 i = iundef; j = iundef;
k = iundef
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) )
735 s31_y(
k,i,j) = ( s31_y(
k,i,j) &
736 + 0.5_rp * ( velx_yz(
k+1,i,j) - velx_yz(
k,i,j) ) * j33g * rfdz(
k) &
742 i = iundef; j = iundef;
k = iundef
752 call check( __line__, velx_c(
k,i,j+1) )
753 call check( __line__, velx_c(
k,i,j-1) )
756 call check( __line__, work_z(
k,i,j) )
757 call check( __line__, work_z(
k-1,i,j) )
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) &
769 i = iundef; j = iundef;
k = iundef
775 call check( __line__, velx_c(
ks,i,j+1) )
776 call check( __line__, velx_c(
ks,i,j-1) )
779 call check( __line__, velx_c(
ks+1,i,j) )
780 call check( __line__, velx_c(
ks,i,j) )
784 call check( __line__, velx_c(
ke,i,j+1) )
785 call check( __line__, velx_c(
ke,i,j-1) )
788 call check( __line__, velx_c(
ke,i,j) )
789 call check( __line__, velx_c(
ke-1,i,j) )
793 call check( __line__, fdy(j) )
794 call check( __line__, fdy(j-1) )
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) &
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) &
807 i = iundef; j = iundef;
k = iundef
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) )
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) &
830 i = iundef; j = iundef;
k = iundef
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) )
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) )
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) &
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) &
862 i = iundef; j = iundef;
k = iundef
866 work_z(:,:,:) = undef; work_x(:,:,:) = undef; work_y(:,:,:) = undef; work_v(:,:,:) = undef
875 call check( __line__, vely_c(
k+1,i,j) )
876 call check( __line__, vely_c(
k,i,j) )
878 work_z(
k,i,j) = 0.5_rp * ( vely_c(
k+1,i,j) + vely_c(
k,i,j) )
883 i = iundef; j = iundef;
k = iundef
891 call check( __line__, vely_c(
k,i+1,j) )
892 call check( __line__, vely_c(
k,i,j) )
894 work_x(
k,i,j) = 0.5_rp * ( velz_c(
k,i+1,j) + velz_c(
k,i,j) )
899 i = iundef; j = iundef;
k = iundef
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) )
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) )
923 i = iundef; j = iundef;
k = iundef
933 call check( __line__, vely_zx(
k,i,j) )
934 call check( __line__, vely_zx(
k,i,j-1) )
937 call check( __line__, work_z(
k,i,j) )
938 call check( __line__, work_z(
k-1,i,j) )
941 call check( __line__, rcdy(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) &
951 i = iundef; j = iundef;
k = iundef
957 call check( __line__, vely_zx(
ks,i,j) )
958 call check( __line__, vely_zx(
ks,i,j-1) )
961 call check( __line__, vely_c(
ks+1,i,j) )
962 call check( __line__, vely_c(
ks,i,j) )
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) )
970 call check( __line__, vely_c(
ke,i,j) )
971 call check( __line__, vely_c(
ke-1,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) &
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) &
986 i = iundef; j = iundef;
k = iundef
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) )
1001 call check( __line__, work_z(
k,i,j) )
1002 call check( __line__, work_z(
k-1,i,j) )
1006 call check( __line__, fdx(i) )
1007 call check( __line__, fdx(i-1) )
1009 s12_c(
k,i,j) = ( s12_c(
k,i,j) &
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) &
1018 i = iundef; j = iundef;
k = iundef
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) )
1029 call check( __line__, vely_c(
ks+1,i,j) )
1030 call check( __line__, vely_c(
ks,i,j) )
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) )
1039 call check( __line__, vely_c(
ke,i,j) )
1040 call check( __line__, vely_c(
ke-1,i,j) )
1044 call check( __line__, fdx(i) )
1045 call check( __line__, fdx(i-1) )
1047 s12_c(
ks,i,j) = ( s12_c(
ks,i,j) &
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) &
1053 s12_c(
ke,i,j) = ( s12_c(
ke,i,j) &
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) &
1062 i = iundef; j = iundef;
k = iundef
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) )
1077 s12_z(
k,i,j) = ( s12_z(
k,i,j) &
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) &
1086 i = iundef; j = iundef;
k = iundef
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) )
1104 s12_z(
ks,i,j) = ( s12_z(
ks,i,j) &
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) &
1110 s12_z(
ke,i,j) = ( s12_z(
ke,i,j) &
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) &
1119 i = iundef; j = iundef;
k = iundef
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) )
1135 s23_c(
k,i,j) = ( s23_c(
k,i,j) &
1136 + 0.5_rp * ( vely_c(
k+1,i,j) - vely_c(
k-1,i,j) ) * j33g / ( fdz(
k) + fdz(
k-1) ) &
1142 i = iundef; j = iundef;
k = iundef
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) )
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) &
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) &
1166 i = iundef; j = iundef;
k = iundef
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) )
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) &
1187 i = iundef; j = iundef;
k = iundef
1193 work_z(:,:,:) = undef; work_x(:,:,:) = undef; work_y(:,:,:) = undef
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) )
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 )
1216 i = iundef; j = iundef;
k = iundef
1228 DENS, PHI, Kh, FACT, &
1229 GSQRT, J13G, J23G, J33G, MAPF, &
1233 IIS, IIE, JJS, JJE )
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
1268 if ( horizontal )
then
1270 qflx_phi(:,:,:,
zdir) = 0.0_rp
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) )
1286 qflx_phi(
k,i,j,
zdir) = - 0.25_rp &
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 &
1295 i = iundef; j = iundef;
k = iundef
1300 qflx_phi(
ks-1,i,j,
zdir) = 0.0_rp
1301 qflx_phi(
ke ,i,j,
zdir) = 0.0_rp
1305 i = iundef; j = iundef;
k = iundef
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) )
1325 qflx_phi(
k,i,j,
xdir) = - 0.25_rp &
1326 * ( dens(
k,i,j) + dens(
k,i+1,j) ) &
1327 * ( kh(
k,i,j) + kh(
k,i+1,j) ) * fact &
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) &
1339 i = iundef; j = iundef;
k = iundef
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) )
1353 qflx_phi(
ks,i,j,
xdir) = - 0.25_rp &
1354 * ( dens(
ks,i,j) + dens(
ks,i+1,j) ) &
1355 * ( kh(
ks,i,j) + kh(
ks,i+1,j) ) * fact &
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) &
1363 qflx_phi(
ke,i,j,
xdir) = - 0.25_rp &
1364 * ( dens(
ke,i,j) + dens(
ke,i+1,j) ) &
1365 * ( kh(
ke,i,j) + kh(
ke,i+1,j) ) * fact &
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) &
1376 i = iundef; j = iundef;
k = iundef
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) )
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 &
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) &
1409 i = iundef; j = iundef;
k = iundef
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) )
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 &
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) &
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 &
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) &
1446 i = iundef; j = iundef;
k = iundef
1449 if ( (.not. horizontal) .and.
implicit )
then
1452 gsqrt, j13g, j23g, j33g, mapf, &
1453 iis, iie, jjs, jje )
1466 a(:,i,j), b(:,i,j), c(:,i,j), d, &
1470 qflx_phi(
k,i,j,
zdir) = qflx_phi(
k,i,j,
zdir) &
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 &
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
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
1513 phi(ke_tb) = ( d(ke_tb) - c(ke_tb)*f(ke_tb-1) ) / ( b(ke_tb) + c(ke_tb)*e(ke_tb-1) )
1515 do k = ke_tb-1,
ks, -1
1516 phi(
k) = e(
k) * phi(
k+1) + f(
k)
1532 GSQRT, J13G, J23G, J33G, MAPF, &
1533 IIS, IIE, JJS, JJE )
1540 real(
rp),
intent(out) :: momz_t_tb(
ka,
ia,
ja)
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
1553 real(
rp) :: fluxz(
ka)
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) &
1569 + j33g * qflx_momz(
k,i,j,
zdir)
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) &
1592 GSQRT, J13G, J23G, J33G, MAPF, &
1593 IIS, IIE, JJS, JJE )
1599 real(
rp),
intent(out) :: momx_t_tb(
ka,
ia,
ja)
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
1612 real(
rp) :: fluxz(
ka)
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) &
1628 + j33g * qflx_momx(
k,i,j,
zdir)
1630 fluxz(
ks-1) = 0.0_rp
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) &
1651 GSQRT, J13G, J23G, J33G, MAPF, &
1652 IIS, IIE, JJS, JJE )
1659 real(
rp),
intent(out) :: momy_t_tb(
ka,
ia,
ja)
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
1672 real(
rp) :: fluxz(
ka)
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) &
1688 + j33g * qflx_momy(
k,i,j,
zdir)
1690 fluxz(
ks-1) = 0.0_rp
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) &
1711 GSQRT, J13G, J23G, J33G, MAPF, &
1712 IIS, IIE, JJS, JJE )
1720 real(
rp),
intent(out) :: phi_t_tb(
ka,
ia,
ja)
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
1733 real(
rp) :: fluxz(0:
ka)
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) &
1750 + j33g * qflx_phi(
k,i,j,
zdir)
1752 fluxz(
ks-1) = 0.0_rp
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) &
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) &
integer, public ke
end point of inner domain: z, local
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rfdx
reciprocal of face-dx
integer, parameter, public const_undef2
undefined value (INT2)
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, public jblock
block size for cache blocking: y
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rcdx
reciprocal of center-dx
subroutine, public atmos_phy_tb_calc_tend_momx(MOMX_t_TB, QFLX_MOMX, GSQRT, J13G, J23G, J33G, MAPF, IIS, IIE, JJS, JJE)
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)
subroutine, public atmos_phy_tb_diffusion_solver(phi, a, b, c, d, KE_TB)
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fdy
y-length of grid(j+1) to grid(j) [m]
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rcdz
reciprocal of center-dz
integer, parameter, public rp
integer, public ie
end point of inner domain: x, local
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rcdy
reciprocal of center-dy
module atmosphere / grid / cartesC index
subroutine, public check(current_line, v)
Undefined value checker.
integer, parameter, public zdir
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)
integer, public is
start point of inner domain: x, local
integer, parameter, public dp
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fdx
x-length of grid(i+1) to grid(i) [m]
subroutine, public atmos_phy_tb_calc_tend_phi(phi_t_TB, QFLX_phi, GSQRT, J13G, J23G, J33G, MAPF, IIS, IIE, JJS, JJE)
integer, parameter, public xdir
integer, public ks
start point of inner domain: z, local
module ATMOSPHERE / Physics Turbulence
integer, public js
start point of inner domain: y, local
integer, parameter, public ydir
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rfdz
reciprocal of face-dz
module atmosphere / grid / cartesC
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fdz
z-length of grid(i+1) to grid(i) [m]
real(rp), public const_undef
subroutine, public atmos_phy_tb_calc_tend_momz(MOMZ_t_TB, QFLX_MOMZ, GSQRT, J13G, J23G, J33G, MAPF, IIS, IIE, JJS, JJE)
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rfdy
reciprocal of face-dy
integer, public je
end point of inner domain: y, local