180 real(RP),
intent(out) :: qflx_sgs_momz(
ka,
ia,
ja,3)
181 real(RP),
intent(out) :: qflx_sgs_momx(
ka,
ia,
ja,3)
182 real(RP),
intent(out) :: qflx_sgs_momy(
ka,
ia,
ja,3)
183 real(RP),
intent(out) :: qflx_sgs_rhot(
ka,
ia,
ja,3)
184 real(RP),
intent(out) :: qflx_sgs_rhoq(
ka,
ia,
ja,3,
qa)
186 real(RP),
intent(inout) :: tke(
ka,
ia,
ja)
187 real(RP),
intent(out) :: tke_t(
ka,
ia,
ja)
188 real(RP),
intent(out) :: km (
ka,
ia,
ja)
189 real(RP),
intent(out) :: ri (
ka,
ia,
ja)
190 real(RP),
intent(out) :: pr (
ka,
ia,
ja)
191 real(RP),
intent(out) :: n2 (
ka,
ia,
ja)
193 real(RP),
intent(in) :: momz(
ka,
ia,
ja)
194 real(RP),
intent(in) :: momx(
ka,
ia,
ja)
195 real(RP),
intent(in) :: momy(
ka,
ia,
ja)
196 real(RP),
intent(in) :: rhot(
ka,
ia,
ja)
197 real(RP),
intent(in) :: dens(
ka,
ia,
ja)
198 real(RP),
intent(in) :: qtrc(
ka,
ia,
ja,
qa)
200 real(RP),
intent(in) :: sflx_mw(
ia,
ja)
201 real(RP),
intent(in) :: sflx_mu(
ia,
ja)
202 real(RP),
intent(in) :: sflx_mv(
ia,
ja)
203 real(RP),
intent(in) :: sflx_sh(
ia,
ja)
204 real(RP),
intent(in) :: sflx_qv(
ia,
ja)
206 real(RP),
intent(in) :: gsqrt (
ka,
ia,
ja,7)
207 real(RP),
intent(in) :: j13g (
ka,
ia,
ja,7)
208 real(RP),
intent(in) :: j23g (
ka,
ia,
ja,7)
209 real(RP),
intent(in) :: j33g
210 real(RP),
intent(in) :: mapf(
ia,
ja,2,4)
211 real(DP),
intent(in) :: dt
214 real(RP) :: pott (
ka,
ia,
ja)
217 real(RP) :: s33_c (
ka,
ia,
ja)
218 real(RP) :: s11_c (
ka,
ia,
ja)
219 real(RP) :: s22_c (
ka,
ia,
ja)
220 real(RP) :: s31_c (
ka,
ia,
ja)
221 real(RP) :: s12_c (
ka,
ia,
ja)
222 real(RP) :: s23_c (
ka,
ia,
ja)
223 real(RP) :: s12_z (
ka,
ia,
ja)
224 real(RP) :: s23_x (
ka,
ia,
ja)
225 real(RP) :: s31_y (
ka,
ia,
ja)
231 real(RP) :: qflx_tke(
ka,
ia,
ja,3)
233 real(RP) :: tend(
ka,
ia,
ja)
243 integer :: k, i, j, iq
246 if( io_l )
write(io_fid_log,*)
'*** Physics step: Turbulence(D1980)' 249 qflx_sgs_momz(:,:,:,:) = undef
250 qflx_sgs_momx(:,:,:,:) = undef
251 qflx_sgs_momy(:,:,:,:) = undef
252 qflx_sgs_rhot(:,:,:,:) = undef
253 qflx_sgs_rhoq(:,:,:,:,:) = undef
271 call check( __line__, rhot(k,i,j) )
272 call check( __line__, dens(k,i,j) )
274 pott(k,i,j) = rhot(k,i,j) / dens(k,i,j)
279 i = iundef; j = iundef; k = iundef
284 call calc_strain_tensor( &
285 s33_c, s11_c, s22_c, &
286 s31_c, s12_c, s23_c, &
287 s12_z, s23_x, s31_y, &
289 dens, momz, momx, momy, &
290 gsqrt, j13g, j23g, j33g, mapf )
302 call check( __line__, pott(k+1,i,j) )
303 call check( __line__, pott(k,i,j) )
304 call check( __line__, pott(k-1,i,j) )
305 call check( __line__, fdz(k) )
306 call check( __line__, fdz(k-1) )
307 call check( __line__, s2(k,i,j) )
309 n2(k,i,j) = grav * ( pott(k+1,i,j) - pott(k-1,i,j) ) * j33g &
310 / ( ( fdz(k) + fdz(k-1) ) * gsqrt(k,i,j,
i_xyz) * pott(k,i,j) )
311 ri(k,i,j) = n2(k,i,j) / s2(k,i,j)
316 i = iundef; j = iundef; k = iundef
321 call check( __line__, pott(
ks+1,i,j) )
322 call check( __line__, pott(
ks,i,j) )
323 call check( __line__, rfdz(
ks) )
324 call check( __line__, s2(
ks,i,j) )
326 n2(
ks,i,j) = grav * ( pott(
ks+1,i,j) - pott(
ks,i,j) ) * j33g &
328 ri(
ks,i,j) = grav * ( pott(
ks+1,i,j) - pott(
ks,i,j) ) * j33g * rfdz(
ks) &
329 / ( gsqrt(
ks,i,j,
i_xyz) * pott(
ks,i,j) * s2(
ks,i,j) )
333 i = iundef; j = iundef; k = iundef
338 call check( __line__, pott(
ke,i,j) )
339 call check( __line__, pott(
ke-1,i,j) )
340 call check( __line__, rfdz(
ke-1) )
341 call check( __line__, s2(
ke,i,j) )
343 n2(
ke,i,j) = grav * ( pott(
ke,i,j) - pott(
ke-1,i,j) ) * j33g &
344 / ( fdz(
ke-1) * gsqrt(
ke,i,j,
i_xyz) * pott(
ke,i,j) )
345 ri(
ke,i,j) = grav * ( pott(
ke,i,j) - pott(
ke-1,i,j) ) * j33g * rfdz(
ke-1) &
346 / ( gsqrt(
ke,i,j,
i_xyz) * pott(
ke,i,j) * s2(
ke,i,j) )
350 i = iundef; j = iundef; k = iundef
357 if ( n2(k,i,j) > 1e-10_rp )
then 358 l(k,i,j) = max( min( 0.76_rp * sqrt(tke(k,i,j)/n2(k,i,j)), delta(k,i,j) ), 1e-10_rp )
360 l(k,i,j) = delta(k,i,j)
366 i = iundef; j = iundef; k = iundef
372 km(k,i,j) = min( 0.10_rp * l(k,i,j) * sqrt(tke(k,i,j)), atmos_phy_tb_d1980_km_max )
373 pr(k,i,j) = 1.0_rp / ( 1.0_rp + 2.0_rp * l(k,i,j) / delta(k,i,j) )
374 kh(k,i,j) = km(k,i,j) / pr(k,i,j)
379 i = iundef; j = iundef; k = iundef
388 call check( __line__, dens(k,i,j) )
389 call check( __line__, km(k,i,j) )
390 call check( __line__, s33_c(k,i,j) )
392 qflx_sgs_momz(k,i,j,
zdir) = dens(k,i,j) * ( - 2.0_rp * km(k,i,j) * s33_c(k,i,j) )
397 i = iundef; j = iundef; k = iundef
401 qflx_sgs_momz(
ks,i,j,
zdir) = 0.0_rp
402 qflx_sgs_momz(
ke,i,j,
zdir) = 0.0_rp
406 i = iundef; j = iundef; k = iundef
413 call check( __line__, dens(k,i,j) )
414 call check( __line__, dens(k,i+1,j) )
415 call check( __line__, dens(k+1,i,j) )
416 call check( __line__, dens(k+1,i+1,j) )
417 call check( __line__, km(k,i,j) )
418 call check( __line__, km(k,i+1,j) )
419 call check( __line__, km(k+1,i,j) )
420 call check( __line__, km(k+1,i+1,j) )
421 call check( __line__, s31_y(k,i,j) )
423 qflx_sgs_momz(k,i,j,
xdir) = - 0.125_rp &
424 * ( dens(k,i,j)+dens(k+1,i,j)+dens(k,i+1,j)+dens(k+1,i+1,j) ) &
425 * ( km(k,i,j)+km(k+1,i,j)+km(k,i+1,j)+km(k+1,i+1,j)) &
431 i = iundef; j = iundef; k = iundef
438 call check( __line__, dens(k,i,j) )
439 call check( __line__, dens(k,i,j+1) )
440 call check( __line__, dens(k+1,i,j) )
441 call check( __line__, dens(k+1,i,j+1) )
442 call check( __line__, km(k,i,j) )
443 call check( __line__, km(k,i,j+1) )
444 call check( __line__, km(k+1,i,j) )
445 call check( __line__, km(k+1,i,j+1) )
446 call check( __line__, s23_x(k,i,j) )
448 qflx_sgs_momz(k,i,j,
ydir) = - 0.125_rp &
449 * ( dens(k,i,j)+dens(k+1,i,j)+dens(k,i,j+1)+dens(k+1,i,j+1) ) &
450 * ( km(k,i,j)+km(k+1,i,j)+km(k,i,j+1)+km(k+1,i,j+1) ) &
456 i = iundef; j = iundef; k = iundef
459 if ( atmos_phy_tb_d1980_implicit )
then 460 call calc_tend_momz( tend, &
462 gsqrt, j13g, j23g, j33g, mapf, &
468 ap = - fouroverthree * dt &
469 * dens(
ks+1,i,j)*km(
ks+1,i,j) &
473 b(
ks,i,j) = - a(
ks,i,j) + 0.5_rp * ( dens(
ks,i,j)+dens(
ks+1,i,j) )
475 c(k,i,j) = ap * rfdz(k+1) / gsqrt(k+1,i,j,
i_xyw)
476 ap = - fouroverthree * dt &
477 * dens(k+1,i,j)*km(k+1,i,j) &
478 * rcdz(k+1) / gsqrt(k+1,i,j,
i_xyz)
479 a(k,i,j) = ap * rfdz(k) / gsqrt(k,i,j,
i_xyw)
480 b(k,i,j) = - a(k,i,j) - c(k,i,j) + 0.5_rp * ( dens(k,i,j)+dens(k+1,i,j) )
484 b(
ke-1,i,j) = - c(
ke-1,i,j) + 0.5_rp * ( dens(
ke-1,i,j)+dens(
ke,i,j) )
490 call diffusion_solver( &
492 a(:,i,j), b(:,i,j), c(:,i,j), d, &
496 qflx_sgs_momz(k,i,j,
zdir) = qflx_sgs_momz(k,i,j,
zdir) &
497 - fouroverthree * dens(k,i,j) * km(k,i,j) * dt &
498 * ( tend(k,i,j) - tend(k-1,i,j) ) * rcdz(k) / gsqrt(k,i,j,
i_xyz)
512 call check( __line__, dens(k,i,j) )
513 call check( __line__, dens(k,i+1,j) )
514 call check( __line__, dens(k+1,i,j) )
515 call check( __line__, dens(k+1,i+1,j) )
516 call check( __line__, km(k,i,j) )
517 call check( __line__, km(k,i+1,j) )
518 call check( __line__, km(k+1,i,j) )
519 call check( __line__, km(k+1,i+1,j) )
520 call check( __line__, s31_y(k,i,j) )
522 qflx_sgs_momx(k,i,j,
zdir) = - 0.125_rp &
523 * ( dens(k,i,j)+dens(k+1,i,j)+dens(k,i+1,j)+dens(k+1,i+1,j) ) &
524 * ( km(k,i,j)+km(k+1,i,j)+km(k,i+1,j)+km(k+1,i+1,j) ) &
530 i = iundef; j = iundef; k = iundef
534 qflx_sgs_momx(
ks-1,i,j,
zdir) = 0.0_rp
535 qflx_sgs_momx(
ke ,i,j,
zdir) = 0.0_rp
539 i = iundef; j = iundef; k = iundef
546 call check( __line__, dens(k,i,j) )
547 call check( __line__, km(k,i,j) )
548 call check( __line__, s11_c(k,i,j) )
550 qflx_sgs_momx(k,i,j,
xdir) = dens(k,i,j) * ( - 2.0_rp * km(k,i,j) * s11_c(k,i,j) )
555 i = iundef; j = iundef; k = iundef
562 call check( __line__, dens(k,i,j) )
563 call check( __line__, dens(k,i+1,j) )
564 call check( __line__, dens(k,i,j+1) )
565 call check( __line__, dens(k,i+1,j+1) )
566 call check( __line__, km(k,i,j) )
567 call check( __line__, km(k,i+1,j) )
568 call check( __line__, km(k,i,j+1) )
569 call check( __line__, km(k,i+1,j+1) )
570 call check( __line__, s12_z(k,i,j) )
572 qflx_sgs_momx(k,i,j,
ydir) = - 0.125_rp &
573 * ( dens(k,i,j)+dens(k,i+1,j)+dens(k,i,j+1)+dens(k,i+1,j+1) ) &
574 * ( km(k,i,j)+km(k,i+1,j)+km(k,i,j+1)+km(k,i+1,j+1) ) &
580 i = iundef; j = iundef; k = iundef
583 if ( atmos_phy_tb_d1980_implicit )
then 584 call calc_tend_momx( tend, &
586 gsqrt, j13g, j23g, j33g, mapf, &
592 ap = - dt * 0.25_rp * ( dens(
ks ,i ,j)*km(
ks ,i ,j) &
593 + dens(
ks+1,i ,j)*km(
ks+1,i ,j) &
594 + dens(
ks ,i+1,j)*km(
ks ,i+1,j) &
595 + dens(
ks+1,i+1,j)*km(
ks+1,i+1,j) ) &
599 b(
ks,i,j) = - a(
ks,i,j) + 0.5_rp * ( dens(
ks,i,j)+dens(
ks,i+1,j) )
601 c(k,i,j) = ap * rcdz(k) / gsqrt(k,i,j,
i_uyz)
602 ap = - dt * 0.25_rp * ( dens(k ,i ,j)*km(k ,i ,j) &
603 + dens(k+1,i ,j)*km(k+1,i ,j) &
604 + dens(k ,i+1,j)*km(k ,i+1,j) &
605 + dens(k+1,i+1,j)*km(k+1,i+1,j) ) &
606 * rfdz(k) / gsqrt(k,i,j,
i_uyw)
607 a(k,i,j) = ap * rcdz(k) / gsqrt(k,i,j,
i_uyz)
608 b(k,i,j) = - a(k,i,j) - c(k,i,j) + 0.5_rp * ( dens(k,i,j)+dens(k,i+1,j) )
612 b(
ke,i,j) = - c(
ke,i,j) + 0.5_rp * ( dens(
ke,i,j)+dens(
ke,i+1,j) )
618 call diffusion_solver( &
620 a(:,i,j), b(:,i,j), c(:,i,j), d, &
624 qflx_sgs_momx(k,i,j,
zdir) = qflx_sgs_momx(k,i,j,
zdir) &
625 - 0.25_rp * ( dens(k ,i ,j)*km(k ,i ,j) &
626 + dens(k+1,i ,j)*km(k+1,i ,j) &
627 + dens(k ,i+1,j)*km(k ,i+1,j) &
628 + dens(k+1,i+1,j)*km(k+1,i+1,j) ) &
629 * dt * ( tend(k+1,i,j) - tend(k,i,j) ) * rfdz(k) / gsqrt(k,i,j,
i_uyw)
643 call check( __line__, dens(k,i,j) )
644 call check( __line__, dens(k,i,j+1) )
645 call check( __line__, dens(k+1,i,j) )
646 call check( __line__, dens(k+1,i,j+1) )
647 call check( __line__, km(k,i,j) )
648 call check( __line__, km(k,i,j+1) )
649 call check( __line__, km(k+1,i,j) )
650 call check( __line__, km(k+1,i,j+1) )
651 call check( __line__, s23_x(k,i,j) )
653 qflx_sgs_momy(k,i,j,
zdir) = - 0.125_rp &
654 * ( dens(k,i,j)+dens(k+1,i,j)+dens(k,i,j+1)+dens(k+1,i,j+1) ) &
655 * ( km(k,i,j)+km(k+1,i,j)+km(k,i,j+1)+km(k+1,i,j+1) ) &
661 i = iundef; j = iundef; k = iundef
665 qflx_sgs_momy(
ks-1,i,j,
zdir) = 0.0_rp
666 qflx_sgs_momy(
ke ,i,j,
zdir) = 0.0_rp
670 i = iundef; j = iundef; k = iundef
678 call check( __line__, dens(k,i,j) )
679 call check( __line__, dens(k,i+1,j) )
680 call check( __line__, dens(k,i,j+1) )
681 call check( __line__, dens(k,i+1,j+1) )
682 call check( __line__, km(k,i,j) )
683 call check( __line__, km(k,i+1,j) )
684 call check( __line__, km(k,i,j+1) )
685 call check( __line__, km(k,i+1,j+1) )
686 call check( __line__, s12_z(k,i,j) )
688 qflx_sgs_momy(k,i,j,
xdir) = - 0.125_rp &
689 * ( dens(k,i,j)+dens(k,i+1,j)+dens(k,i,j+1)+dens(k,i+1,j+1) ) &
690 * ( km(k,i,j)+km(k,i+1,j)+km(k,i,j+1)+km(k,i+1,j+1) ) &
696 i = iundef; j = iundef; k = iundef
704 call check( __line__, dens(k,i,j) )
705 call check( __line__, km(k,i,j) )
706 call check( __line__, s22_c(k,i,j) )
708 qflx_sgs_momy(k,i,j,
ydir) = dens(k,i,j) * ( - 2.0_rp * km(k,i,j) * s22_c(k,i,j) )
713 i = iundef; j = iundef; k = iundef
716 if ( atmos_phy_tb_d1980_implicit )
then 717 call calc_tend_momy( tend, &
719 gsqrt, j13g, j23g, j33g, mapf, &
725 ap = - dt * 0.25_rp * ( dens(
ks ,i,j )*km(
ks ,i,j ) &
726 + dens(
ks+1,i,j )*km(
ks+1,i,j ) &
727 + dens(
ks ,i,j+1)*km(
ks ,i,j+1) &
728 + dens(
ks+1,i,j+1)*km(
ks+1,i,j+1) ) &
732 b(
ks,i,j) = - a(
ks,i,j) + 0.5_rp * ( dens(
ks,i,j)+dens(
ks,i,j+1) )
734 c(k,i,j) = ap * rcdz(k) / gsqrt(k,i,j,
i_xvz)
735 ap = - dt * 0.25_rp * ( dens(k ,i,j )*km(k ,i,j ) &
736 + dens(k+1,i,j )*km(k+1,i,j ) &
737 + dens(k ,i,j+1)*km(k ,i,j+1) &
738 + dens(k+1,i,j+1)*km(k+1,i,j+1) ) &
739 * rfdz(k) / gsqrt(k,i,j,
i_xvw)
740 a(k,i,j) = ap * rcdz(k) / gsqrt(k,i,j,
i_xvz)
741 b(k,i,j) = - a(k,i,j) - c(k,i,j) + 0.5_rp * ( dens(k,i,j)+dens(k,i,j+1) )
745 b(
ke,i,j) = - c(
ke,i,j) + 0.5_rp * ( dens(
ke,i,j)+dens(
ke,i,j+1) )
751 call diffusion_solver( &
753 a(:,i,j), b(:,i,j), c(:,i,j), d, &
757 qflx_sgs_momy(k,i,j,
zdir) = qflx_sgs_momy(k,i,j,
zdir) &
758 - 0.25_rp * ( dens(k ,i,j )*km(k ,i,j ) &
759 + dens(k+1,i,j )*km(k+1,i,j ) &
760 + dens(k ,i,j+1)*km(k ,i,j+1) &
761 + dens(k+1,i,j+1)*km(k+1,i,j+1) ) &
762 * dt * ( tend(k+1,i,j) - tend(k,i,j) ) * rfdz(k) / gsqrt(k,i,j,
i_xvw)
772 if ( atmos_phy_tb_d1980_implicit )
then 777 ap = - dt * 0.25_rp * ( dens(
ks,i,j)+dens(
ks+1,i,j) ) &
778 * ( kh(
ks,i,j)+kh(
ks+1,i,j) ) &
782 b(
ks,i,j) = - a(
ks,i,j) + dens(
ks,i,j)
784 c(k,i,j) = ap * rcdz(k) / gsqrt(k,i,j,
i_xyz)
785 ap = - dt * 0.25_rp * ( dens(k,i,j)+dens(k+1,i,j) ) &
786 * ( kh(k,i,j)+kh(k+1,i,j) ) &
787 * rfdz(k) / gsqrt(k,i,j,
i_xyw)
788 a(k,i,j) = ap * rcdz(k) / gsqrt(k,i,j,
i_xyz)
789 b(k,i,j) = - a(k,i,j) - c(k,i,j) + dens(k,i,j)
793 b(
ke,i,j) = - c(
ke,i,j) + dens(
ke,i,j)
800 call calc_flux_phi( &
802 dens, pott, kh, 1.0_rp, &
803 gsqrt, j13g, j23g, j33g, mapf, &
805 atmos_phy_tb_d1980_implicit, &
820 call calc_flux_phi( &
821 qflx_sgs_rhoq(:,:,:,:,iq), &
822 dens, qtrc(:,:,:,iq), kh, 1.0_rp, &
823 gsqrt, j13g, j23g, j33g, mapf, &
825 atmos_phy_tb_d1980_implicit, &
831 iis = iundef; iie = iundef; jjs = iundef; jje = iundef
845 if ( atmos_phy_tb_d1980_implicit )
then 850 ap = - dt * 0.25_rp * ( dens(
ks,i,j)+dens(
ks+1,i,j) ) &
851 * 2.0_rp * ( km(
ks,i,j)+km(
ks+1,i,j) ) &
855 b(
ks,i,j) = - a(
ks,i,j) + dens(
ks,i,j)
857 c(k,i,j) = ap * rcdz(k) / gsqrt(k,i,j,
i_xyz)
858 ap = - dt * 0.25_rp * ( dens(k,i,j)+dens(k+1,i,j) ) &
859 * 2.0_rp * ( km(k,i,j)+km(k+1,i,j) ) &
860 * rfdz(k) / gsqrt(k,i,j,
i_xyw)
861 a(k,i,j) = ap * rcdz(k) / gsqrt(k,i,j,
i_xyz)
862 b(k,i,j) = - a(k,i,j) - c(k,i,j) + dens(k,i,j)
866 b(
ke,i,j) = - c(
ke,i,j) + dens(
ke,i,j)
873 call calc_flux_phi( &
875 dens, tke, km, 2.0_rp, &
876 gsqrt, j13g, j23g, j33g, mapf, &
878 atmos_phy_tb_d1980_implicit, &
884 qflx_tke(k,i,j,
zdir) = qflx_tke(k,i,j,
zdir) &
885 + 0.5_rp * ( momz(k,i,j) * ( tke(k+1,i,j)+tke(k,i,j) ) &
886 - abs(momz(k,i,j)) * ( tke(k+1,i,j)-tke(k,i,j) ) )
894 qflx_tke(k,i,j,
xdir) = qflx_tke(k,i,j,
xdir) &
895 + 0.5_rp * ( momx(k,i,j) * ( tke(k,i+1,j)+tke(k,i,j) ) &
896 - abs(momx(k,i,j)) * ( tke(k,i+1,j)-tke(k,i,j) ) )
904 qflx_tke(k,i,j,
ydir) = qflx_tke(k,i,j,
ydir) &
905 + 0.5_rp * ( momy(k,i,j) * ( tke(k,i,j+1)+tke(k,i,j) ) &
906 - abs(momy(k,i,j)) * ( tke(k,i,j+1)-tke(k,i,j) ) )
911 call calc_tend_phi( tke_t, &
913 gsqrt, j13g, j23g, j33g, mapf, &
918 tke_t(k,i,j) = min( &
920 + km(k,i,j) * s2(k,i,j) &
921 - kh(k,i,j) * n2(k,i,j) &
922 - ( c1 + c2*l(k,i,j)/delta(k,i,j) ) * tke(k,i,j)**(1.5_rp) / l(k,i,j), &
928 i = iundef; j = iundef; k = iundef
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
real(rp), dimension(:), allocatable, public grid_rcdz
reciprocal of center-dz
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
real(rp), public const_grav
standard acceleration of gravity [m/s2]
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, 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), public const_eps
small number
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)